escript  Revision_
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 <ripley/Ripley.h>
21 #include <ripley/AbstractAssembler.h>
22 #include <ripley/RipleyException.h>
23 
24 #include <escript/AbstractContinuousDomain.h>
25 #include <escript/Data.h>
26 #include <escript/FunctionSpace.h>
27 #include <escript/SubWorld.h>
28 
29 #ifdef ESYS_HAVE_PASO
30 #include <paso/Coupler.h>
31 #include <paso/SystemMatrix.h>
32 #endif
33 
34 #ifdef ESYS_HAVE_TRILINOS
35 #include <trilinoswrap/types.h>
36 #endif
37 
38 #include <boost/python/list.hpp>
39 #include <boost/python/tuple.hpp>
40 
41 namespace ripley {
42 
47 };
48 
50  SMT_PASO = 1<<8,
51  SMT_CUSP = 1<<9,
52  SMT_TRILINOS = 1<<10,
53  SMT_SYMMETRIC = 1<<15,
54  SMT_COMPLEX = 1<<16,
55  SMT_UNROLL = 1<<17
56 };
57 
62 };
63 
69 {
71  std::vector<dim_t> first;
73  std::vector<dim_t> numValues;
76  std::vector<int> multiplier;
78  std::vector<int> reverse;
80  int byteOrder;
82  int dataType;
83 };
84 
89 struct DiracPoint
90 {
92  int tag;
93 };
94 
102 {
103 public:
109 
114  ~RipleyDomain();
115 
116  static void setDecompositionPolicy(DecompositionPolicy value);
117  static DecompositionPolicy getDecompositionPolicy();
118 
123  virtual escript::JMPI getMPI() const { return m_mpiInfo; }
124 
129  virtual int getMPISize() const { return m_mpiInfo->size; }
130 
135  virtual int getMPIRank() const { return m_mpiInfo->rank; }
136 
141  virtual void MPIBarrier() const {
142 #ifdef ESYS_MPI
143  MPI_Barrier(m_mpiInfo->comm);
144 #endif
145  }
146 
151  virtual bool onMasterProcessor() const { return getMPIRank()==0; }
152 
157  MPI_Comm getMPIComm() const { return m_mpiInfo->comm; }
158 
164  virtual bool isValidFunctionSpaceType(int fsType) const;
165 
170  virtual std::string functionSpaceTypeAsString(int fsType) const;
171 
176  virtual int getDim() const { return m_numDim; }
177 
181  virtual bool operator==(const escript::AbstractDomain& other) const;
182 
186  virtual bool operator!=(const escript::AbstractDomain& other) const {
187  return !(operator==(other));
188  }
189 
196  virtual std::pair<int,dim_t> getDataShape(int fsType) const;
197 
204  int getTagFromSampleNo(int fsType, dim_t sampleNo) const;
205 
212  virtual void setTagMap(const std::string& name, int tag) {
213  m_tagMap[name] = tag;
214  }
215 
221  virtual int getTag(const std::string& name) const {
222  if (m_tagMap.find(name) != m_tagMap.end()) {
223  return m_tagMap.find(name)->second;
224  } else {
225  throw escript::ValueError("getTag: invalid tag name");
226  }
227  }
228 
234  virtual bool isValidTagName(const std::string& name) const {
235  return (m_tagMap.find(name)!=m_tagMap.end());
236  }
237 
242  virtual std::string showTagNames() const;
243 
249  virtual void setNewX(const escript::Data& arg);
250 
256  virtual void interpolateOnDomain(escript::Data& target, const escript::Data& source) const;
257 
263  virtual bool probeInterpolationOnDomain(int fsType_source, int fsType_target) const;
264 
272  virtual signed char preferredInterpolationOnDomain(int fsType_source,
273  int fsType_target) const;
274 
281  bool
282  commonFunctionSpace(const std::vector<int>& fs, int& resultcode) const;
283 
289  virtual void interpolateAcross(escript::Data& target,
290  const escript::Data& source) const;
291 
296  virtual bool probeInterpolationAcross(int, const escript::AbstractDomain&, int) const;
297 
302  virtual escript::Data getX() const;
303 
308  virtual escript::Data getNormal() const;
309 
313  virtual escript::Data getSize() const;
314 
320  virtual void setToX(escript::Data& arg) const;
321 
328  virtual void setToGradient(escript::Data& out, const escript::Data& in) const;
329 
335  virtual void setTags(int fsType, int newTag, const escript::Data& mask) const;
336 
342  virtual bool isCellOriented(int fsType) const;
343 
350  virtual StatusType getStatus() const { return m_status; }
351 
356  virtual int getNumberOfTagsInUse(int fsType) const;
357 
362  virtual const int* borrowListOfTagsInUse(int fsType) const;
363 
368  virtual bool canTag(int fsType) const;
369 
374  virtual int getApproximationOrder(int fsType) const { return 1; }
375 
380  virtual bool supportsContactElements() const { return false; }
381 
386  virtual int getContinuousFunctionCode() const { return Nodes; }
387 
392  virtual int getReducedContinuousFunctionCode() const { return ReducedNodes; }
393 
398  virtual int getFunctionCode() const { return Elements; }
399 
404  virtual int getReducedFunctionCode() const { return ReducedElements; }
405 
410  virtual int getFunctionOnBoundaryCode() const { return FaceElements; }
411 
418 
423  virtual int getFunctionOnContactZeroCode() const {
424  throw escript::NotImplementedError("Ripley does not support contact elements");
425  }
426 
432  throw escript::NotImplementedError("Ripley does not support contact elements");
433  }
434 
439  virtual int getFunctionOnContactOneCode() const {
440  throw escript::NotImplementedError("Ripley does not support contact elements");
441  }
442 
448  throw escript::NotImplementedError("Ripley does not support contact elements");
449  }
450 
455  virtual int getSolutionCode() const { return DegreesOfFreedom; }
456 
461  virtual int getReducedSolutionCode() const { return ReducedDegreesOfFreedom; }
462 
467  virtual int getDiracDeltaFunctionsCode() const { return Points; }
468 
477  virtual int getSystemMatrixTypeId(const boost::python::object& options) const;
478 
488  virtual int getTransportTypeId(int solver, int preconditioner, int package,
489  bool symmetry) const;
490 
496  virtual void setToIntegrals(DoubleVector& integrals, const escript::Data& arg) const;
497 
498 
504  virtual void addToSystem(escript::AbstractSystemMatrix& mat,
505  escript::Data& rhs, const DataMap& data,
506  Assembler_ptr assembler) const;
507 
512  virtual void addToSystemFromPython(escript::AbstractSystemMatrix& mat,
513  escript::Data& rhs, const boost::python::list& data,
514  Assembler_ptr assembler) const;
515 
521  virtual void addToRHS(escript::Data& rhs, const DataMap& data,
522  Assembler_ptr assembler) const;
523 
528  virtual void addToRHSFromPython(escript::Data& rhs,
529  const boost::python::list& data,
530  Assembler_ptr assembler) const;
531 
536  virtual void addPDEToTransportProblem(escript::AbstractTransportProblem& tp,
537  escript::Data& source, const DataMap& data,
538  Assembler_ptr assembler) const;
542  virtual void addPDEToTransportProblem(
544  const escript::Data& M,
545  const escript::Data& A, const escript::Data& B, const escript::Data& C,const escript::Data& D,
546  const escript::Data& X,const escript::Data& Y,
547  const escript::Data& d, const escript::Data& y,
548  const escript::Data& d_contact,const escript::Data& y_contact,
549  const escript::Data& d_dirac,const escript::Data& y_dirac) const;
550 
555  void addPDEToTransportProblemFromPython(
557  escript::Data& source, const boost::python::list& data,
558  Assembler_ptr assembler) const;
559 
564  virtual escript::ASM_ptr newSystemMatrix(int row_blocksize,
565  const escript::FunctionSpace& row_functionspace,
566  int column_blocksize,
567  const escript::FunctionSpace& column_functionspace, int type) const;
568 
573  virtual escript::ATP_ptr newTransportProblem(int blocksize,
574  const escript::FunctionSpace& functionspace,
575  int type) const;
576 
582  virtual void Print_Mesh_Info(bool full=false) const;
583 
584 
585  /************************************************************************/
586 
592  virtual void write(const std::string& filename) const = 0;
593 
598  virtual std::string getDescription() const = 0;
599 
605  void dump(const std::string& filename) const = 0;
606 
612  const dim_t* borrowSampleReferenceIDs(int fsType) const = 0;
613 
620  virtual void setToNormal(escript::Data& out) const = 0;
621 
627  virtual void setToSize(escript::Data& out) const = 0;
628 
633  virtual void readNcGrid(escript::Data& out, std::string filename,
634  std::string varname, const ReaderParameters& params) const = 0;
635 
640  virtual void readBinaryGrid(escript::Data& out, std::string filename,
641  const ReaderParameters& params) const = 0;
642 
647  virtual void readBinaryGridFromZipped(escript::Data& out,
648  std::string filename, const ReaderParameters& params) const = 0;
649 
654  virtual void writeBinaryGrid(const escript::Data& in, std::string filename,
655  int byteOrder, int dataType) const = 0;
656 
661  virtual bool ownSample(int fsType, index_t id) const = 0;
662 
667  virtual dim_t getNumDataPointsGlobal() const = 0;
668 
673  virtual const dim_t* getNumNodesPerDim() const = 0;
674 
679  virtual const dim_t* getNumElementsPerDim() const = 0;
680 
686  virtual const dim_t* getNumFacesPerBoundary() const = 0;
687 
692  virtual IndexVector getNodeDistribution() const = 0;
693 
698  virtual const int* getNumSubdivisionsPerDim() const = 0;
699 
704  virtual double getLocalCoordinate(index_t index, int dim) const = 0;
705 
710  virtual boost::python::tuple getGridParameters() const = 0;
711 
716  virtual const double *getLength() const = 0;
717 
722  virtual const double *getElementLength() const = 0;
723 
729  virtual RankVector getOwnerVector(int fsType) const = 0;
730 
736  virtual bool supportsFilter(const boost::python::tuple& t) const;
737 
741  virtual Assembler_ptr createAssembler(std::string type,
742  const DataMap& options) const {
743  throw escript::NotImplementedError("Domain does not support custom assemblers");
744  }
745 
749  Assembler_ptr createAssemblerFromPython(std::string type,
750  const boost::python::list& options) const;
751 
752 
753 protected:
754  int m_numDim;
758  mutable std::vector<int> m_nodeTags, m_nodeTagsInUse;
759  mutable std::vector<int> m_elementTags, m_elementTagsInUse;
760  mutable std::vector<int> m_faceTags, m_faceTagsInUse;
761  std::vector<DiracPoint> m_diracPoints;
762  IndexVector m_diracPointNodeIDs; //for borrowSampleID
764 
766  void copyData(escript::Data& out, const escript::Data& in) const;
767 
769  void averageData(escript::Data& out, const escript::Data& in) const;
770 
772  void multiplyData(escript::Data& out, const escript::Data& in) const;
773 
774  // this is const because setTags is const
775  void updateTagsInUse(int fsType) const;
776 
777 #ifdef ESYS_HAVE_PASO
778  void createPasoConnector(const RankVector& neighbour,
780  const IndexVector& offsetInSharedSend,
781  const IndexVector& offsetInSharedRecv,
782  const IndexVector& sendShared,
783  const IndexVector& recvShared);
784 
787  paso::Connector_ptr getPasoConnector() const { return m_connector; }
788 
790  paso::Pattern_ptr createPasoPattern(const std::vector<IndexVector>& indices,
791  dim_t N) const;
792 #endif
793 
794 #ifdef ESYS_HAVE_TRILINOS
795  esys_trilinos::const_TrilinosGraph_ptr createTrilinosGraph(
798  const IndexVector& myRows, const IndexVector& myColumns) const;
799 #endif
800 
801  template<typename Scalar>
802  void addToSystemMatrix(escript::AbstractSystemMatrix* mat,
803  const IndexVector& nodes, dim_t numEq,
804  const std::vector<Scalar>& array) const;
805 
806  void addPoints(const std::vector<double>& coords,
807  const std::vector<int>& tags);
808 
809  /***********************************************************************/
810 
812  virtual dim_t getNumNodes() const = 0;
813 
815  virtual dim_t getNumElements() const = 0;
816 
818  virtual dim_t getNumDOF() const = 0;
819 
821  virtual dim_t getNumFaceElements() const = 0;
822 
824  virtual IndexVector getDiagonalIndices(bool upperOnly) const = 0;
825 
827  virtual void assembleCoordinates(escript::Data& arg) const = 0;
828 
830  virtual void assembleGradient(escript::Data& out, const escript::Data& in) const = 0;
831 
833  virtual void assembleIntegrate(DoubleVector& integrals, const escript::Data& arg) const = 0;
834 
835 #ifdef ESYS_HAVE_TRILINOS
836  virtual esys_trilinos::const_TrilinosGraph_ptr getTrilinosGraph() const = 0;
838 #endif
839 
841  virtual std::vector<IndexVector> getConnections(bool includeShared) const = 0;
842 
843 #ifdef ESYS_HAVE_PASO
844  virtual paso::SystemMatrixPattern_ptr getPasoMatrixPattern(
846  bool reducedRowOrder,
847  bool reducedColOrder) const = 0;
848 #endif
849 
851  virtual void interpolateNodesOnElements(escript::Data& out,
852  const escript::Data& in,
853  bool reduced) const = 0;
854 
856  virtual void interpolateNodesOnFaces(escript::Data& out,
857  const escript::Data& in,
858  bool reduced) const = 0;
859 
861  virtual void nodesToDOF(escript::Data& out, const escript::Data& in) const = 0;
862 
864  virtual void dofToNodes(escript::Data& out, const escript::Data& in) const;
865 
866  virtual dim_t getDofOfNode(dim_t node) const = 0;
867 
868 private:
870 
871 #ifdef ESYS_HAVE_PASO
872  // Paso connector used by the system matrix and to interpolate DOF to
873  // nodes
874  paso::Connector_ptr m_connector;
875 
877  void addToPasoMatrix(paso::SystemMatrix* in, const IndexVector& nodes,
878  dim_t numEq, const DoubleVector& array) const;
879 #endif
880 
882  void assemblePDE(escript::AbstractSystemMatrix* mat, escript::Data& rhs,
883  const DataMap& coefs, Assembler_ptr assembler) const;
884 
887  void assemblePDEBoundary(escript::AbstractSystemMatrix* mat,
888  escript::Data& rhs, const DataMap& coefs,
889  Assembler_ptr assembler) const;
890 
891  void assemblePDEDirac(escript::AbstractSystemMatrix* mat,
892  escript::Data& rhs, const DataMap& coefs,
893  Assembler_ptr assembler) const;
894 
896  virtual dim_t findNode(const double *coords) const = 0;
897 };
898 
899 } // end of namespace ripley
900 
901 #endif // __RIPLEY_DOMAIN_H__
902 
AbstractContinuousDomain, base class for continuous domains.
Definition: AbstractContinuousDomain.h:45
Definition: FunctionSpace.h:34
boost::shared_ptr< Pattern > Pattern_ptr
Definition: Pattern.h:37
virtual bool onMasterProcessor() const
returns true if on MPI processor 0, else false
Definition: ripley/src/RipleyDomain.h:151
int dataType
data type in the file (used by binary reader only)
Definition: ripley/src/RipleyDomain.h:82
Definition: Ripley.h:52
SystemMatrixType
Definition: ripley/src/RipleyDomain.h:49
Definition: Ripley.h:55
dim_t node
Definition: ripley/src/RipleyDomain.h:91
bool probeInterpolationAcross(int fsType_source, const escript::AbstractDomain &domain, int fsType_target, int dim)
Definition: CrossDomainCoupler.cpp:31
Definition: ripley/src/RipleyDomain.h:59
Definition: Ripley.h:53
virtual bool supportsContactElements() const
returns true if this domain supports contact elements, false otherwise
Definition: ripley/src/RipleyDomain.h:380
std::vector< int > RankVector
Definition: Ripley.h:45
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:116
Definition: ripley/src/RipleyDomain.h:61
virtual int getReducedFunctionOnContactZeroCode() const
returns a FunctionOnContactZero code with reduced integration order
Definition: ripley/src/RipleyDomain.h:431
static DecompositionPolicy m_decompPolicy
Definition: ripley/src/RipleyDomain.h:869
virtual int getReducedFunctionCode() const
returns a function with reduced integration order FunctionSpace code
Definition: ripley/src/RipleyDomain.h:404
virtual Assembler_ptr createAssembler(std::string type, const DataMap &options) const
Definition: ripley/src/RipleyDomain.h:741
Definition: ripley/src/RipleyDomain.h:45
virtual int getMPIRank() const
returns the MPI rank of this processor
Definition: ripley/src/RipleyDomain.h:135
Definition: ripley/src/RipleyDomain.h:60
int byteOrder
byte order in the file (used by binary reader only)
Definition: ripley/src/RipleyDomain.h:80
virtual void MPIBarrier() const
if compiled for MPI then executes an MPI_Barrier, else does nothing
Definition: ripley/src/RipleyDomain.h:141
std::vector< int > m_faceTagsInUse
Definition: ripley/src/RipleyDomain.h:760
A struct to contain a dirac point&#39;s information.
Definition: ripley/src/RipleyDomain.h:89
Definition: ripley/src/RipleyDomain.h:44
boost::shared_ptr< JMPI_ > JMPI
Definition: EsysMPI.h:71
assembler_t assembler_type
Definition: ripley/src/RipleyDomain.h:763
Definition: ripley/src/RipleyDomain.h:51
std::vector< int > m_nodeTagsInUse
Definition: ripley/src/RipleyDomain.h:758
Definition: Ripley.h:51
Definition: ripley/src/RipleyDomain.h:50
static dim_t M
Definition: SparseMatrix_saveHB.cpp:37
std::vector< int > reverse
if non-zero, values are written from last index to first index
Definition: ripley/src/RipleyDomain.h:78
Definition: Ripley.h:56
virtual int getFunctionOnBoundaryCode() const
returns a function on boundary FunctionSpace code
Definition: ripley/src/RipleyDomain.h:410
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:350
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:62
assembler_t
Definition: ripley/src/RipleyDomain.h:43
boost::shared_ptr< SystemMatrixPattern > SystemMatrixPattern_ptr
Definition: SystemMatrixPattern.h:39
boost::shared_ptr< AbstractTransportProblem > ATP_ptr
Definition: AbstractTransportProblem.h:160
std::map< std::string, int > TagMap
Definition: Ripley.h:46
virtual int getApproximationOrder(int fsType) const
returns the approximation order used for a function space
Definition: ripley/src/RipleyDomain.h:374
IndexVector m_diracPointNodeIDs
Definition: ripley/src/RipleyDomain.h:762
int index_t
type for array/matrix indices used both globally and on each rank
Definition: DataTypes.h:59
std::vector< index_t > IndexVector
Definition: Ripley.h:43
StatusType m_status
Definition: ripley/src/RipleyDomain.h:755
Definition: ripley/src/RipleyDomain.h:52
std::map< std::string, escript::Data > DataMap
Definition: ripley/src/domainhelpers.h:24
virtual int getReducedFunctionOnContactOneCode() const
returns a FunctionOnContactOne code with reduced integration order
Definition: ripley/src/RipleyDomain.h:447
boost::shared_ptr< SubWorld > SubWorld_ptr
Definition: SubWorld.h:146
An exception class for features which are not (yet) implemented.
Definition: EsysException.h:78
virtual int getSolutionCode() const
returns a Solution FunctionSpace code
Definition: ripley/src/RipleyDomain.h:455
virtual int getDim() const
returns the number of spatial dimensions of the domain
Definition: ripley/src/RipleyDomain.h:176
Structure that wraps parameters for the grid reading routines.
Definition: ripley/src/RipleyDomain.h:68
Definition: Ripley.h:50
std::vector< int > multiplier
Definition: ripley/src/RipleyDomain.h:76
int StatusType
Definition: AbstractDomain.h:48
this class holds a (distributed) stiffness matrix
Definition: SystemMatrix.h:47
std::vector< dim_t > numValues
the number of values to read from file
Definition: ripley/src/RipleyDomain.h:73
virtual int getFunctionOnContactOneCode() const
returns a FunctionOnContactOne code
Definition: ripley/src/RipleyDomain.h:439
virtual int getReducedContinuousFunctionCode() const
returns a continuous on reduced order nodes FunctionSpace code
Definition: ripley/src/RipleyDomain.h:392
Data represents a collection of datapoints.
Definition: Data.h:63
DecompositionPolicy
Definition: ripley/src/RipleyDomain.h:58
int m_numDim
Definition: ripley/src/RipleyDomain.h:754
Definition: ripley/src/RipleyDomain.h:53
virtual bool operator!=(const escript::AbstractDomain &other) const
inequality operator
Definition: ripley/src/RipleyDomain.h:186
virtual int getReducedSolutionCode() const
returns a ReducedSolution FunctionSpace code
Definition: ripley/src/RipleyDomain.h:461
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:212
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:101
Definition: Ripley.h:57
Definition: Ripley.h:54
int MPI_Comm
Definition: EsysMPI.h:41
virtual int getReducedFunctionOnBoundaryCode() const
returns a function on boundary with reduced integration order FunctionSpace code
Definition: ripley/src/RipleyDomain.h:417
virtual int getFunctionOnContactZeroCode() const
return a FunctionOnContactZero code
Definition: ripley/src/RipleyDomain.h:423
int tag
Definition: ripley/src/RipleyDomain.h:92
virtual int getDiracDeltaFunctionsCode() const
returns a DiracDeltaFunctions FunctionSpace code
Definition: ripley/src/RipleyDomain.h:467
Give a short description of what AbstractTransportProblem does.
Definition: AbstractTransportProblem.h:43
An exception class that signals an invalid argument value.
Definition: EsysException.h:88
boost::shared_ptr< Connector > Connector_ptr
Definition: Coupler.h:37
Definition: ripley/src/RipleyDomain.h:46
virtual bool isValidTagName(const std::string &name) const
returns true if name is a defined tag name
Definition: ripley/src/RipleyDomain.h:234
std::vector< dim_t > first
the (global) offset into the data object to start writing into
Definition: ripley/src/RipleyDomain.h:71
Base class for escript system matrices.
Definition: AbstractSystemMatrix.h:42
virtual escript::JMPI getMPI() const
returns a reference to the MPI information wrapper for this domain
Definition: ripley/src/RipleyDomain.h:123
Definition: ripley/src/AbstractAssembler.h:25
std::vector< real_t > DoubleVector
Definition: Ripley.h:44
escript::Data readBinaryGridFromZipped(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:87
virtual int getMPISize() const
returns the number of processors used for this domain
Definition: ripley/src/RipleyDomain.h:129
std::vector< DiracPoint > m_diracPoints
Definition: ripley/src/RipleyDomain.h:761
boost::shared_ptr< AbstractSystemMatrix > ASM_ptr
Definition: AbstractSystemMatrix.h:32
virtual int getTag(const std::string &name) const
returns the tag key for tag name
Definition: ripley/src/RipleyDomain.h:221
std::vector< int > m_elementTagsInUse
Definition: ripley/src/RipleyDomain.h:759
Base class for all escript domains.
Definition: AbstractDomain.h:45
MPI_Comm getMPIComm() const
returns the MPI communicator
Definition: ripley/src/RipleyDomain.h:157
TagMap m_tagMap
Definition: ripley/src/RipleyDomain.h:757
virtual int getContinuousFunctionCode() const
returns a continuous FunctionSpace code
Definition: ripley/src/RipleyDomain.h:386
Definition: Ripley.h:49
#define RIPLEY_DLL_API
Definition: ripley/src/system_dep.h:20
Definition: ripley/src/RipleyDomain.h:54
escript::JMPI m_mpiInfo
Definition: ripley/src/RipleyDomain.h:756
index_t dim_t
Definition: DataTypes.h:64
Definition: ripley/src/RipleyDomain.h:55
virtual int getFunctionCode() const
returns a function FunctionSpace code
Definition: ripley/src/RipleyDomain.h:398