escript  Revision_Unversioneddirectory
ripley/src/RipleyDomain.h
Go to the documentation of this file.
1 
2 /*****************************************************************************
3 *
4 * Copyright (c) 2003-2016 by The University of Queensland
5 * http://www.uq.edu.au
6 *
7 * Primary Business: Queensland, Australia
8 * Licensed under the Apache License, version 2.0
9 * http://www.apache.org/licenses/LICENSE-2.0
10 *
11 * Development until 2012 by Earth Systems Science Computational Center (ESSCC)
12 * Development 2012-2013 by School of Earth Sciences
13 * Development from 2014 by Centre for Geoscience Computing (GeoComp)
14 *
15 *****************************************************************************/
16 
17 #ifndef __RIPLEY_DOMAIN_H__
18 #define __RIPLEY_DOMAIN_H__
19 
20 #include <boost/python/tuple.hpp>
21 #include <boost/python/list.hpp>
22 
23 #include <ripley/Ripley.h>
24 #include <ripley/RipleyException.h>
25 #include <ripley/AbstractAssembler.h>
26 
27 #include <escript/AbstractContinuousDomain.h>
28 #include <escript/Data.h>
29 #include <escript/FunctionSpace.h>
30 #include <escript/SubWorld.h>
31 
32 #include <paso/SystemMatrix.h>
33 
34 namespace ripley {
35 
40 };
41 
43  SMT_PASO = 1024,
44  SMT_CUSP = 2048,
46 };
47 
53 {
55  std::vector<dim_t> first;
57  std::vector<dim_t> numValues;
60  std::vector<int> multiplier;
62  std::vector<int> reverse;
64  int byteOrder;
66  int dataType;
67 };
68 
73 struct DiracPoint
74 {
76  int tag;
77 };
78 
86 {
87 public:
93 
98  ~RipleyDomain();
99 
104  virtual int getMPISize() const { return m_mpiInfo->size; }
105 
110  virtual int getMPIRank() const { return m_mpiInfo->rank; }
111 
116  virtual void MPIBarrier() const {
117 #ifdef ESYS_MPI
118  MPI_Barrier(m_mpiInfo->comm);
119 #endif
120  }
121 
126  virtual bool onMasterProcessor() const { return getMPIRank()==0; }
127 
132  MPI_Comm getMPIComm() const { return m_mpiInfo->comm; }
133 
139  virtual bool isValidFunctionSpaceType(int fsType) const;
140 
145  virtual std::string functionSpaceTypeAsString(int fsType) const;
146 
151  virtual int getDim() const { return m_numDim; }
152 
156  virtual bool operator==(const escript::AbstractDomain& other) const;
157 
161  virtual bool operator!=(const escript::AbstractDomain& other) const {
162  return !(operator==(other));
163  }
164 
171  virtual std::pair<int,dim_t> getDataShape(int fsType) const;
172 
179  int getTagFromSampleNo(int fsType, dim_t sampleNo) const;
180 
187  virtual void setTagMap(const std::string& name, int tag) {
188  m_tagMap[name] = tag;
189  }
190 
196  virtual int getTag(const std::string& name) const {
197  if (m_tagMap.find(name) != m_tagMap.end()) {
198  return m_tagMap.find(name)->second;
199  } else {
200  throw RipleyException("getTag: invalid tag name");
201  }
202  }
203 
209  virtual bool isValidTagName(const std::string& name) const {
210  return (m_tagMap.find(name)!=m_tagMap.end());
211  }
212 
217  virtual std::string showTagNames() const;
218 
224  virtual void setNewX(const escript::Data& arg);
225 
231  virtual void interpolateOnDomain(escript::Data& target, const escript::Data& source) const;
232 
238  virtual bool probeInterpolationOnDomain(int fsType_source, int fsType_target) const;
239 
247  virtual signed char preferredInterpolationOnDomain(int fsType_source,
248  int fsType_target) const;
249 
256  bool
257  commonFunctionSpace(const std::vector<int>& fs, int& resultcode) const;
258 
264  virtual void interpolateAcross(escript::Data& target,
265  const escript::Data& source) const;
266 
271  virtual bool probeInterpolationAcross(int, const escript::AbstractDomain&, int) const;
272 
277  virtual escript::Data getX() const;
278 
283  virtual escript::Data getNormal() const;
284 
288  virtual escript::Data getSize() const;
289 
295  virtual void setToX(escript::Data& arg) const;
296 
303  virtual void setToGradient(escript::Data& out, const escript::Data& in) const;
304 
310  virtual void setTags(int fsType, int newTag, const escript::Data& mask) const;
311 
317  virtual bool isCellOriented(int fsType) const;
318 
325  virtual StatusType getStatus() const { return m_status; }
326 
331  virtual int getNumberOfTagsInUse(int fsType) const;
332 
337  virtual const int* borrowListOfTagsInUse(int fsType) const;
338 
343  virtual bool canTag(int fsType) const;
344 
349  virtual int getApproximationOrder(int fsType) const { return 1; }
350 
355  virtual bool supportsContactElements() const { return false; }
356 
361  virtual int getContinuousFunctionCode() const { return Nodes; }
362 
367  virtual int getReducedContinuousFunctionCode() const { return ReducedNodes; }
368 
373  virtual int getFunctionCode() const { return Elements; }
374 
379  virtual int getReducedFunctionCode() const { return ReducedElements; }
380 
385  virtual int getFunctionOnBoundaryCode() const { return FaceElements; }
386 
393 
398  virtual int getFunctionOnContactZeroCode() const {
399  throw RipleyException("Ripley does not support contact elements");
400  }
401 
407  throw RipleyException("Ripley does not support contact elements");
408  }
409 
414  virtual int getFunctionOnContactOneCode() const {
415  throw RipleyException("Ripley does not support contact elements");
416  }
417 
423  throw RipleyException("Ripley does not support contact elements");
424  }
425 
430  virtual int getSolutionCode() const { return DegreesOfFreedom; }
431 
436  virtual int getReducedSolutionCode() const { return ReducedDegreesOfFreedom; }
437 
442  virtual int getDiracDeltaFunctionsCode() const { return Points; }
443 
452  virtual int getSystemMatrixTypeId(const boost::python::object& options) const;
453 
463  virtual int getTransportTypeId(int solver, int preconditioner, int package,
464  bool symmetry) const;
465 
471  virtual void setToIntegrals(DoubleVector& integrals, const escript::Data& arg) const;
472 
473 
479  virtual void addToSystem(escript::AbstractSystemMatrix& mat,
480  escript::Data& rhs, const DataMap& data,
481  Assembler_ptr assembler) const;
482 
487  virtual void addToSystemFromPython(escript::AbstractSystemMatrix& mat,
488  escript::Data& rhs, const boost::python::list& data,
489  Assembler_ptr assembler) const;
490 
496  virtual void addToRHS(escript::Data& rhs, const DataMap& data,
497  Assembler_ptr assembler) const;
498 
503  virtual void addToRHSFromPython(escript::Data& rhs,
504  const boost::python::list& data,
505  Assembler_ptr assembler) const;
506 
511  virtual void addPDEToTransportProblem(escript::AbstractTransportProblem& tp,
512  escript::Data& source, const DataMap& data,
513  Assembler_ptr assembler) const;
517  virtual void addPDEToTransportProblem(
519  const escript::Data& M,
520  const escript::Data& A, const escript::Data& B, const escript::Data& C,const escript::Data& D,
521  const escript::Data& X,const escript::Data& Y,
522  const escript::Data& d, const escript::Data& y,
523  const escript::Data& d_contact,const escript::Data& y_contact,
524  const escript::Data& d_dirac,const escript::Data& y_dirac) const;
525 
530  void addPDEToTransportProblemFromPython(
532  escript::Data& source, const boost::python::list& data,
533  Assembler_ptr assembler) const;
534 
539  virtual escript::ASM_ptr newSystemMatrix(int row_blocksize,
540  const escript::FunctionSpace& row_functionspace,
541  int column_blocksize,
542  const escript::FunctionSpace& column_functionspace, int type) const;
543 
548  virtual escript::ATP_ptr newTransportProblem(int blocksize,
549  const escript::FunctionSpace& functionspace,
550  int type) const;
551 
557  virtual void Print_Mesh_Info(bool full=false) const;
558 
559 
560  /************************************************************************/
561 
567  virtual void write(const std::string& filename) const = 0;
568 
573  virtual std::string getDescription() const = 0;
574 
580  void dump(const std::string& filename) const = 0;
581 
587  const dim_t* borrowSampleReferenceIDs(int fsType) const = 0;
588 
595  virtual void setToNormal(escript::Data& out) const = 0;
596 
602  virtual void setToSize(escript::Data& out) const = 0;
603 
608  virtual void readNcGrid(escript::Data& out, std::string filename,
609  std::string varname, const ReaderParameters& params) const = 0;
610 
615  virtual void readBinaryGrid(escript::Data& out, std::string filename,
616  const ReaderParameters& params) const = 0;
617 
618 #ifdef USE_BOOSTIO
619 
623  virtual void readBinaryGridFromZipped(escript::Data& out,
624  std::string filename, const ReaderParameters& params) const = 0;
625 #endif
626 
631  virtual void writeBinaryGrid(const escript::Data& in, std::string filename,
632  int byteOrder, int dataType) const = 0;
633 
638  virtual bool ownSample(int fsType, index_t id) const = 0;
639 
644  virtual dim_t getNumDataPointsGlobal() const = 0;
645 
650  virtual const dim_t* getNumNodesPerDim() const = 0;
651 
656  virtual const dim_t* getNumElementsPerDim() const = 0;
657 
663  virtual const dim_t* getNumFacesPerBoundary() const = 0;
664 
669  virtual IndexVector getNodeDistribution() const = 0;
670 
675  virtual const int* getNumSubdivisionsPerDim() const = 0;
676 
681  virtual double getLocalCoordinate(index_t index, int dim) const = 0;
682 
687  virtual boost::python::tuple getGridParameters() const = 0;
688 
693  virtual const double *getLength() const = 0;
694 
699  virtual const double *getElementLength() const = 0;
700 
706  virtual RankVector getOwnerVector(int fsType) const = 0;
707 
713  virtual bool supportsFilter(const boost::python::tuple& t) const;
714 
718  virtual Assembler_ptr createAssembler(std::string type,
719  const DataMap& options) const {
720  throw RipleyException("Domain does not support custom assemblers");
721  }
722 
726  Assembler_ptr createAssemblerFromPython(std::string type,
727  const boost::python::list& options) const;
728 
729 
730 protected:
731  int m_numDim;
735  mutable std::vector<int> m_nodeTags, m_nodeTagsInUse;
736  mutable std::vector<int> m_elementTags, m_elementTagsInUse;
737  mutable std::vector<int> m_faceTags, m_faceTagsInUse;
738  std::vector<DiracPoint> m_diracPoints;
739  IndexVector m_diracPointNodeIDs; //for borrowSampleID
741 
743  void copyData(escript::Data& out, const escript::Data& in) const;
744 
746  void averageData(escript::Data& out, const escript::Data& in) const;
747 
749  void multiplyData(escript::Data& out, const escript::Data& in) const;
750 
751  // this is const because setTags is const
752  void updateTagsInUse(int fsType) const;
753 
755  paso::Pattern_ptr createPasoPattern(const std::vector<IndexVector>& indices,
756  dim_t N) const;
757 
758  void addToSystemMatrix(escript::AbstractSystemMatrix* mat,
759  const IndexVector& nodes, dim_t numEq,
760  const DoubleVector& array) const;
761 
762  void addPoints(const std::vector<double>& coords,
763  const std::vector<int>& tags);
764 
765  /***********************************************************************/
766 
768  virtual dim_t getNumNodes() const = 0;
769 
771  virtual dim_t getNumElements() const = 0;
772 
774  virtual dim_t getNumDOF() const = 0;
775 
777  virtual dim_t getNumFaceElements() const = 0;
778 
780  virtual IndexVector getDiagonalIndices(bool upperOnly) const = 0;
781 
783  virtual void assembleCoordinates(escript::Data& arg) const = 0;
784 
786  virtual void assembleGradient(escript::Data& out, const escript::Data& in) const = 0;
787 
789  virtual void assembleIntegrate(DoubleVector& integrals, const escript::Data& arg) const = 0;
790 
792  virtual paso::SystemMatrixPattern_ptr getPasoMatrixPattern(
793  bool reducedRowOrder,
794  bool reducedColOrder) const = 0;
795 
797  virtual void interpolateNodesOnElements(escript::Data& out,
798  const escript::Data& in,
799  bool reduced) const = 0;
800 
802  virtual void interpolateNodesOnFaces(escript::Data& out,
803  const escript::Data& in,
804  bool reduced) const = 0;
805 
807  virtual void nodesToDOF(escript::Data& out, const escript::Data& in) const = 0;
808 
810  virtual void dofToNodes(escript::Data& out, const escript::Data& in) const = 0;
811 
812  virtual dim_t getDofOfNode(dim_t node) const = 0;
813 
814 private:
816  void addToSystemMatrix(paso::SystemMatrix_ptr in, const IndexVector& nodes,
817  dim_t numEq, const DoubleVector& array) const;
818 
820  void assemblePDE(escript::AbstractSystemMatrix* mat, escript::Data& rhs,
821  const DataMap& coefs, Assembler_ptr assembler) const;
822 
825  void assemblePDEBoundary(escript::AbstractSystemMatrix* mat,
826  escript::Data& rhs, const DataMap& coefs,
827  Assembler_ptr assembler) const;
828 
829  void assemblePDEDirac(escript::AbstractSystemMatrix* mat,
830  escript::Data& rhs, const DataMap& coefs,
831  Assembler_ptr assembler) const;
832 
834  virtual dim_t findNode(const double *coords) const = 0;
835 };
836 
837 } // end of namespace ripley
838 
839 #endif // __RIPLEY_DOMAIN_H__
840 
AbstractContinuousDomain, base class for continuous domains.
Definition: AbstractContinuousDomain.h:46
Definition: FunctionSpace.h:34
boost::shared_ptr< Pattern > Pattern_ptr
Definition: Pattern.h:36
int dataType
data type in the file (used by binary reader only)
Definition: ripley/src/RipleyDomain.h:66
SystemMatrixType
Definition: ripley/src/RipleyDomain.h:42
dim_t node
Definition: ripley/src/RipleyDomain.h:75
bool probeInterpolationAcross(int fsType_source, const escript::AbstractDomain &domain, int fsType_target, int dim)
Definition: CrossDomainCoupler.cpp:36
virtual bool onMasterProcessor() const
returns true if on MPI processor 0, else false
Definition: ripley/src/RipleyDomain.h:126
virtual int getFunctionOnBoundaryCode() const
returns a function on boundary FunctionSpace code
Definition: ripley/src/RipleyDomain.h:385
Definition: Ripley.h:49
escript::Data readNcGrid(std::string filename, std::string varname, escript::FunctionSpace fs, const object &pyShape, double fill, const object &pyFirst, const object &pyNum, const object &pyMultiplier, const object &pyReverse)
Definition: ripleycpp.cpp:117
virtual int getFunctionOnContactOneCode() const
returns a FunctionOnContactOne code
Definition: ripley/src/RipleyDomain.h:414
int getLength(const escript::Data *data)
Definition: DataC.cpp:92
virtual bool operator!=(const escript::AbstractDomain &other) const
inequality operator
Definition: ripley/src/RipleyDomain.h:161
Definition: ripley/src/RipleyDomain.h:38
int byteOrder
byte order in the file (used by binary reader only)
Definition: ripley/src/RipleyDomain.h:64
esysUtils::JMPI m_mpiInfo
Definition: ripley/src/RipleyDomain.h:733
std::vector< int > m_faceTagsInUse
Definition: ripley/src/RipleyDomain.h:737
A struct to contain a dirac point&#39;s information.
Definition: ripley/src/RipleyDomain.h:73
Definition: ripley/src/RipleyDomain.h:37
virtual int getReducedContinuousFunctionCode() const
returns a continuous on reduced order nodes FunctionSpace code
Definition: ripley/src/RipleyDomain.h:367
Definition: Ripley.h:52
assembler_t assembler_type
Definition: ripley/src/RipleyDomain.h:740
Definition: ripley/src/RipleyDomain.h:44
RipleyException exception class.
Definition: RipleyException.h:29
std::vector< int > m_nodeTagsInUse
Definition: ripley/src/RipleyDomain.h:735
Definition: Ripley.h:50
Definition: ripley/src/RipleyDomain.h:43
static dim_t M
Definition: SparseMatrix_saveHB.cpp:37
virtual bool supportsContactElements() const
returns true if this domain supports contact elements, false otherwise
Definition: ripley/src/RipleyDomain.h:355
std::vector< int > reverse
if non-zero, values are written from last index to first index
Definition: ripley/src/RipleyDomain.h:62
escript::Data readBinaryGrid(std::string filename, escript::FunctionSpace fs, const object &pyShape, double fill, int byteOrder, int dataType, const object &pyFirst, const object &pyNum, const object &pyMultiplier, const object &pyReverse)
Definition: ripleycpp.cpp:65
assembler_t
Definition: ripley/src/RipleyDomain.h:36
virtual int getFunctionOnContactZeroCode() const
return a FunctionOnContactZero code
Definition: ripley/src/RipleyDomain.h:398
boost::shared_ptr< SystemMatrixPattern > SystemMatrixPattern_ptr
Definition: SystemMatrixPattern.h:38
boost::shared_ptr< AbstractTransportProblem > ATP_ptr
Definition: AbstractTransportProblem.h:162
int MPI_Comm
Definition: Esys_MPI.h:38
virtual int getDim() const
returns the number of spatial dimensions of the domain
Definition: ripley/src/RipleyDomain.h:151
boost::shared_ptr< SystemMatrix > SystemMatrix_ptr
Definition: SystemMatrix.h:38
std::map< std::string, int > TagMap
Definition: Ripley.h:41
IndexVector m_diracPointNodeIDs
Definition: ripley/src/RipleyDomain.h:739
virtual int getMPIRank() const
returns the MPI rank of this processor
Definition: ripley/src/RipleyDomain.h:110
Definition: Ripley.h:45
virtual int getApproximationOrder(int fsType) const
returns the approximation order used for a function space
Definition: ripley/src/RipleyDomain.h:349
std::vector< index_t > IndexVector
Definition: Ripley.h:38
StatusType m_status
Definition: ripley/src/RipleyDomain.h:732
std::map< std::string, escript::Data > DataMap
Definition: ripley/src/domainhelpers.h:24
boost::shared_ptr< SubWorld > SubWorld_ptr
Definition: SubWorld.h:142
virtual int getReducedFunctionCode() const
returns a function with reduced integration order FunctionSpace code
Definition: ripley/src/RipleyDomain.h:379
Structure that wraps parameters for the grid reading routines.
Definition: ripley/src/RipleyDomain.h:52
std::vector< int > multiplier
Definition: ripley/src/RipleyDomain.h:60
virtual int getReducedFunctionOnContactOneCode() const
returns a FunctionOnContactOne code with reduced integration order
Definition: ripley/src/RipleyDomain.h:422
int StatusType
Definition: AbstractDomain.h:48
virtual int getDiracDeltaFunctionsCode() const
returns a DiracDeltaFunctions FunctionSpace code
Definition: ripley/src/RipleyDomain.h:442
std::vector< dim_t > numValues
the number of values to read from file
Definition: ripley/src/RipleyDomain.h:57
Data represents a collection of datapoints.
Definition: Data.h:68
int m_numDim
Definition: ripley/src/RipleyDomain.h:731
Definition: ripley/src/RipleyDomain.h:45
virtual Assembler_ptr createAssembler(std::string type, const DataMap &options) const
Definition: ripley/src/RipleyDomain.h:718
virtual int getReducedFunctionOnContactZeroCode() const
returns a FunctionOnContactZero code with reduced integration order
Definition: ripley/src/RipleyDomain.h:406
virtual void setTagMap(const std::string &name, int tag)
sets a map from a clear tag name to a tag key
Definition: ripley/src/RipleyDomain.h:187
static dim_t N
Definition: SparseMatrix_saveHB.cpp:37
RipleyDomain extends the AbstractContinuousDomain interface for the Ripley library and is the base cl...
Definition: ripley/src/RipleyDomain.h:85
virtual bool isValidTagName(const std::string &name) const
returns true if name is a defined tag name
Definition: ripley/src/RipleyDomain.h:209
Definition: Ripley.h:48
virtual int getReducedFunctionOnBoundaryCode() const
returns a function on boundary with reduced integration order FunctionSpace code
Definition: ripley/src/RipleyDomain.h:392
virtual int getMPISize() const
returns the number of processors used for this domain
Definition: ripley/src/RipleyDomain.h:104
virtual int getTag(const std::string &name) const
returns the tag key for tag name
Definition: ripley/src/RipleyDomain.h:196
int index_t
Definition: types.h:24
virtual int getReducedSolutionCode() const
returns a ReducedSolution FunctionSpace code
Definition: ripley/src/RipleyDomain.h:436
virtual int getFunctionCode() const
returns a function FunctionSpace code
Definition: ripley/src/RipleyDomain.h:373
int tag
Definition: ripley/src/RipleyDomain.h:76
MPI_Comm getMPIComm() const
returns the MPI communicator
Definition: ripley/src/RipleyDomain.h:132
std::vector< Esys_MPI_rank > RankVector
Definition: Ripley.h:40
Give a short description of what AbstractTransportProblem does.
Definition: AbstractTransportProblem.h:45
Definition: ripley/src/RipleyDomain.h:39
index_t dim_t
Definition: types.h:27
Definition: Ripley.h:46
Definition: Ripley.h:51
std::vector< dim_t > first
the (global) offset into the data object to start writing into
Definition: ripley/src/RipleyDomain.h:55
Base class for escript system matrices.
Definition: AbstractSystemMatrix.h:37
Definition: ripley/src/AbstractAssembler.h:25
std::vector< DiracPoint > m_diracPoints
Definition: ripley/src/RipleyDomain.h:738
virtual StatusType getStatus() const
returns a status indicator of the domain. The status identifier should be unique over the lifetime of...
Definition: ripley/src/RipleyDomain.h:325
std::vector< int > m_elementTagsInUse
Definition: ripley/src/RipleyDomain.h:736
boost::shared_ptr< AbstractSystemMatrix > ASM_ptr
Definition: AbstractSystemMatrix.h:170
Base class for all escript domains.
Definition: AbstractDomain.h:45
TagMap m_tagMap
Definition: ripley/src/RipleyDomain.h:734
virtual int getSolutionCode() const
returns a Solution FunctionSpace code
Definition: ripley/src/RipleyDomain.h:430
Definition: Ripley.h:44
#define RIPLEY_DLL_API
Definition: ripley/src/system_dep.h:22
boost::shared_ptr< JMPI_ > JMPI
Definition: Esys_MPI.h:79
Definition: Ripley.h:47
std::vector< double > DoubleVector
Definition: Ripley.h:39
virtual int getContinuousFunctionCode() const
returns a continuous FunctionSpace code
Definition: ripley/src/RipleyDomain.h:361
virtual void MPIBarrier() const
if compiled for MPI then executes an MPI_Barrier, else does nothing
Definition: ripley/src/RipleyDomain.h:116