escript  Revision_
ripley/src/RipleyDomain.h
Go to the documentation of this file.
1 
2 /*****************************************************************************
3 *
4 * Copyright (c) 2003-2020 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-2017 by Centre for Geoscience Computing (GeoComp)
14 * Development from 2019 by School of Earth and Environmental Sciences
15 **
16 *****************************************************************************/
17 
18 #ifndef __RIPLEY_DOMAIN_H__
19 #define __RIPLEY_DOMAIN_H__
20 
21 #include <ripley/Ripley.h>
22 #include <ripley/AbstractAssembler.h>
23 #include <ripley/RipleyException.h>
24 
25 #include <escript/AbstractContinuousDomain.h>
26 #include <escript/Data.h>
27 #include <escript/FunctionSpace.h>
28 #include <escript/SubWorld.h>
29 
30 #ifdef ESYS_HAVE_PASO
31 #include <paso/Coupler.h>
32 #include <paso/SystemMatrix.h>
33 #endif
34 
35 #ifdef ESYS_HAVE_TRILINOS
36 #include <trilinoswrap/types.h>
37 #endif
38 
39 #include <boost/python/list.hpp>
40 #include <boost/python/tuple.hpp>
41 
42 namespace ripley {
43 
48 };
49 
51  SMT_PASO = 1<<8,
52  SMT_CUSP = 1<<9,
53  SMT_TRILINOS = 1<<10,
54  SMT_SYMMETRIC = 1<<15,
55  SMT_COMPLEX = 1<<16,
56  SMT_UNROLL = 1<<17
57 };
58 
63 };
64 
70 {
72  std::vector<dim_t> first;
74  std::vector<dim_t> numValues;
77  std::vector<int> multiplier;
79  std::vector<int> reverse;
81  int byteOrder;
83  int dataType;
84 };
85 
90 struct DiracPoint
91 {
93  int tag;
94 };
95 
103 {
104 public:
110 
115  ~RipleyDomain();
116 
117  static void setDecompositionPolicy(DecompositionPolicy value);
118  static DecompositionPolicy getDecompositionPolicy();
119 
124  virtual escript::JMPI getMPI() const { return m_mpiInfo; }
125 
130  virtual int getMPISize() const { return m_mpiInfo->size; }
131 
136  virtual int getMPIRank() const { return m_mpiInfo->rank; }
137 
142  virtual void MPIBarrier() const {
143 #ifdef ESYS_MPI
144  MPI_Barrier(m_mpiInfo->comm);
145 #endif
146  }
147 
152  virtual bool onMasterProcessor() const { return getMPIRank()==0; }
153 
158  MPI_Comm getMPIComm() const { return m_mpiInfo->comm; }
159 
165  virtual bool isValidFunctionSpaceType(int fsType) const;
166 
171  virtual std::string functionSpaceTypeAsString(int fsType) const;
172 
177  virtual int getDim() const { return m_numDim; }
178 
182  virtual bool operator==(const escript::AbstractDomain& other) const;
183 
187  virtual bool operator!=(const escript::AbstractDomain& other) const {
188  return !(operator==(other));
189  }
190 
197  virtual std::pair<int,dim_t> getDataShape(int fsType) const;
198 
205  int getTagFromSampleNo(int fsType, dim_t sampleNo) const;
206 
213  virtual void setTagMap(const std::string& name, int tag) {
214  m_tagMap[name] = tag;
215  }
216 
222  virtual int getTag(const std::string& name) const {
223  if (m_tagMap.find(name) != m_tagMap.end()) {
224  return m_tagMap.find(name)->second;
225  } else {
226  throw escript::ValueError("getTag: invalid tag name");
227  }
228  }
229 
235  virtual bool isValidTagName(const std::string& name) const {
236  return (m_tagMap.find(name)!=m_tagMap.end());
237  }
238 
243  virtual std::string showTagNames() const;
244 
250  virtual void setNewX(const escript::Data& arg);
251 
257  virtual void interpolateOnDomain(escript::Data& target, const escript::Data& source) const;
258 
264  virtual bool probeInterpolationOnDomain(int fsType_source, int fsType_target) const;
265 
273  virtual signed char preferredInterpolationOnDomain(int fsType_source,
274  int fsType_target) const;
275 
282  bool
283  commonFunctionSpace(const std::vector<int>& fs, int& resultcode) const;
284 
290  virtual void interpolateAcross(escript::Data& target,
291  const escript::Data& source) const;
292 
297  virtual bool probeInterpolationAcross(int, const escript::AbstractDomain&, int) const;
298 
303  virtual escript::Data getX() const;
304 
305 #ifdef ESYS_HAVE_BOOST_NUMPY
310  virtual boost::python::numpy::ndarray getNumpyX() const;
311 
315  virtual boost::python::numpy::ndarray getConnectivityInfo() const;
316 #endif
317 
322  virtual escript::Data getNormal() const;
323 
327  virtual escript::Data getSize() const;
328 
334  virtual void setToX(escript::Data& arg) const;
335 
342  virtual void setToGradient(escript::Data& out, const escript::Data& in) const;
343 
349  virtual void setTags(int fsType, int newTag, const escript::Data& mask) const;
350 
356  virtual bool isCellOriented(int fsType) const;
357 
364  virtual StatusType getStatus() const { return m_status; }
365 
370  virtual int getNumberOfTagsInUse(int fsType) const;
371 
376  virtual const int* borrowListOfTagsInUse(int fsType) const;
377 
382  virtual bool canTag(int fsType) const;
383 
388  virtual int getApproximationOrder(int fsType) const { return 1; }
389 
394  virtual bool supportsContactElements() const { return false; }
395 
400  virtual int getContinuousFunctionCode() const { return Nodes; }
401 
406  virtual int getReducedContinuousFunctionCode() const { return ReducedNodes; }
407 
412  virtual int getFunctionCode() const { return Elements; }
413 
418  virtual int getReducedFunctionCode() const { return ReducedElements; }
419 
424  virtual int getFunctionOnBoundaryCode() const { return FaceElements; }
425 
432 
437  virtual int getFunctionOnContactZeroCode() const {
438  throw escript::NotImplementedError("Ripley does not support contact elements");
439  }
440 
446  throw escript::NotImplementedError("Ripley does not support contact elements");
447  }
448 
453  virtual int getFunctionOnContactOneCode() const {
454  throw escript::NotImplementedError("Ripley does not support contact elements");
455  }
456 
462  throw escript::NotImplementedError("Ripley does not support contact elements");
463  }
464 
469  virtual int getSolutionCode() const { return DegreesOfFreedom; }
470 
475  virtual int getReducedSolutionCode() const { return ReducedDegreesOfFreedom; }
476 
481  virtual int getDiracDeltaFunctionsCode() const { return Points; }
482 
491  virtual int getSystemMatrixTypeId(const boost::python::object& options) const;
492 
502  virtual int getTransportTypeId(int solver, int preconditioner, int package,
503  bool symmetry) const;
504 
510  virtual void setToIntegrals(std::vector<real_t>& integrals,
511  const escript::Data& arg) const;
512  virtual void setToIntegrals(std::vector<cplx_t>& integrals,
513  const escript::Data& arg) const;
514 
515 
521  virtual void addToSystem(escript::AbstractSystemMatrix& mat,
522  escript::Data& rhs, const DataMap& data,
523  Assembler_ptr assembler) const;
524 
529  virtual void addToSystemFromPython(escript::AbstractSystemMatrix& mat,
530  escript::Data& rhs, const boost::python::list& data,
531  Assembler_ptr assembler) const;
532 
538  virtual void addToRHS(escript::Data& rhs, const DataMap& data,
539  Assembler_ptr assembler) const;
540 
545  virtual void addToRHSFromPython(escript::Data& rhs,
546  const boost::python::list& data,
547  Assembler_ptr assembler) const;
548 
553  virtual void addPDEToTransportProblem(escript::AbstractTransportProblem& tp,
554  escript::Data& source, const DataMap& data,
555  Assembler_ptr assembler) const;
559  virtual void addPDEToTransportProblem(
561  const escript::Data& M,
562  const escript::Data& A, const escript::Data& B, const escript::Data& C,const escript::Data& D,
563  const escript::Data& X,const escript::Data& Y,
564  const escript::Data& d, const escript::Data& y,
565  const escript::Data& d_contact,const escript::Data& y_contact,
566  const escript::Data& d_dirac,const escript::Data& y_dirac) const;
567 
572  void addPDEToTransportProblemFromPython(
574  escript::Data& source, const boost::python::list& data,
575  Assembler_ptr assembler) const;
576 
581  virtual escript::ASM_ptr newSystemMatrix(int row_blocksize,
582  const escript::FunctionSpace& row_functionspace,
583  int column_blocksize,
584  const escript::FunctionSpace& column_functionspace, int type) const;
585 
590  virtual escript::ATP_ptr newTransportProblem(int blocksize,
591  const escript::FunctionSpace& functionspace,
592  int type) const;
593 
599  virtual void Print_Mesh_Info(bool full=false) const;
600 
601 
602  /************************************************************************/
603 
609  virtual void write(const std::string& filename) const = 0;
610 
615  virtual std::string getDescription() const = 0;
616 
622  void dump(const std::string& filename) const = 0;
623 
629  const dim_t* borrowSampleReferenceIDs(int fsType) const = 0;
630 
637  virtual void setToNormal(escript::Data& out) const = 0;
638 
644  virtual void setToSize(escript::Data& out) const = 0;
645 
650  virtual void readNcGrid(escript::Data& out, std::string filename,
651  std::string varname, const ReaderParameters& params) const = 0;
652 
657  virtual void readBinaryGrid(escript::Data& out, std::string filename,
658  const ReaderParameters& params) const = 0;
659 
665  std::string filename, const ReaderParameters& params) const = 0;
666 
671  virtual void writeBinaryGrid(const escript::Data& in, std::string filename,
672  int byteOrder, int dataType) const = 0;
673 
678  virtual bool ownSample(int fsType, index_t id) const = 0;
679 
684  virtual dim_t getNumDataPointsGlobal() const = 0;
685 
690  virtual const dim_t* getNumNodesPerDim() const = 0;
691 
696  virtual const dim_t* getNumElementsPerDim() const = 0;
697 
703  virtual const dim_t* getNumFacesPerBoundary() const = 0;
704 
709  virtual IndexVector getNodeDistribution() const = 0;
710 
715  virtual const int* getNumSubdivisionsPerDim() const = 0;
716 
721  virtual double getLocalCoordinate(index_t index, int dim) const = 0;
722 
727  virtual boost::python::tuple getGridParameters() const = 0;
728 
733  virtual const double *getLength() const = 0;
734 
739  virtual const double *getElementLength() const = 0;
740 
746  virtual RankVector getOwnerVector(int fsType) const = 0;
747 
753  virtual bool supportsFilter(const boost::python::tuple& t) const;
754 
758  virtual Assembler_ptr createAssembler(std::string type,
759  const DataMap& options) const {
760  throw escript::NotImplementedError("Domain does not support custom assemblers");
761  }
762 
766  Assembler_ptr createAssemblerFromPython(std::string type,
767  const boost::python::list& options) const;
768 
769 
770 protected:
771  int m_numDim;
775  mutable std::vector<int> m_nodeTags, m_nodeTagsInUse;
776  mutable std::vector<int> m_elementTags, m_elementTagsInUse;
777  mutable std::vector<int> m_faceTags, m_faceTagsInUse;
778  std::vector<DiracPoint> m_diracPoints;
779  IndexVector m_diracPointNodeIDs; //for borrowSampleID
781 
783  template<typename Scalar>
784  void copyData(escript::Data& out, const escript::Data& in) const;
785 
787  template<typename Scalar>
788  void averageData(escript::Data& out, const escript::Data& in) const;
789 
791  template<typename Scalar>
792  void multiplyData(escript::Data& out, const escript::Data& in) const;
793 
794  // this is const because setTags is const
795  void updateTagsInUse(int fsType) const;
796 
797 #ifdef ESYS_HAVE_PASO
799  void createPasoConnector(const RankVector& neighbour,
800  const IndexVector& offsetInSharedSend,
801  const IndexVector& offsetInSharedRecv,
802  const IndexVector& sendShared,
803  const IndexVector& recvShared);
804 
807  paso::Connector_ptr getPasoConnector() const { return m_connector; }
808 
810  paso::Pattern_ptr createPasoPattern(const std::vector<IndexVector>& indices,
811  dim_t N) const;
812 #endif
813 
814 #ifdef ESYS_HAVE_TRILINOS
817  esys_trilinos::const_TrilinosGraph_ptr createTrilinosGraph(
818  const IndexVector& myRows, const IndexVector& myColumns) const;
819 #endif
820 
821  template<typename Scalar>
823  const IndexVector& nodes, dim_t numEq,
824  const std::vector<Scalar>& array) const;
825 
826  void addPoints(const std::vector<double>& coords,
827  const std::vector<int>& tags);
828 
829  /***********************************************************************/
830 
832  virtual dim_t getNumNodes() const = 0;
833 
835  virtual dim_t getNumElements() const = 0;
836 
838  virtual dim_t getNumDOF() const = 0;
839 
841  virtual dim_t getNumFaceElements() const = 0;
842 
844  virtual IndexVector getDiagonalIndices(bool upperOnly) const = 0;
845 
847  virtual void assembleCoordinates(escript::Data& arg) const = 0;
848 
850  virtual void assembleGradient(escript::Data& out, const escript::Data& in) const = 0;
851 
853  virtual void assembleIntegrate(std::vector<real_t>& integrals, const escript::Data& arg) const = 0;
854  virtual void assembleIntegrate(std::vector<cplx_t>& integrals, const escript::Data& arg) const = 0;
855 
856 #ifdef ESYS_HAVE_TRILINOS
858  virtual esys_trilinos::const_TrilinosGraph_ptr getTrilinosGraph() const = 0;
859 #endif
860 
862  virtual std::vector<IndexVector> getConnections(bool includeShared) const = 0;
863 
864 #ifdef ESYS_HAVE_PASO
866  virtual paso::SystemMatrixPattern_ptr getPasoMatrixPattern(
867  bool reducedRowOrder,
868  bool reducedColOrder) const = 0;
869 #endif
870 
873  const escript::Data& in,
874  bool reduced) const = 0;
875 
878  const escript::Data& in,
879  bool reduced) const = 0;
880 
882  virtual void nodesToDOF(escript::Data& out, const escript::Data& in) const = 0;
883 
885  template<typename Scalar>
886  void dofToNodes(escript::Data& out, const escript::Data& in) const;
887 
888  virtual dim_t getDofOfNode(dim_t node) const = 0;
889 
890 private:
892 
893 #ifdef ESYS_HAVE_PASO
894  // Paso connector used by the system matrix and to interpolate DOF to
895  // nodes
896  paso::Connector_ptr m_connector;
897 
899  template <typename T>
900  void addToPasoMatrix(paso::SystemMatrix<T>* in, const IndexVector& nodes,
901  dim_t numEq, const std::vector<T>& array) const;
902 #endif
903 
905  void assemblePDE(escript::AbstractSystemMatrix* mat, escript::Data& rhs,
906  const DataMap& coefs, Assembler_ptr assembler) const;
907 
910  void assemblePDEBoundary(escript::AbstractSystemMatrix* mat,
911  escript::Data& rhs, const DataMap& coefs,
912  Assembler_ptr assembler) const;
913 
914  void assemblePDEDirac(escript::AbstractSystemMatrix* mat,
915  escript::Data& rhs, const DataMap& coefs,
916  Assembler_ptr assembler) const;
917 
918  template<typename Scalar>
919  void setToIntegralsWorker(std::vector<Scalar>& integrals,
920  const escript::Data& arg) const;
921 
923  virtual dim_t findNode(const double *coords) const = 0;
924 };
925 
926 } // end of namespace ripley
927 
928 #endif // __RIPLEY_DOMAIN_H__
929 
int MPI_Comm
Definition: EsysMPI.h:44
AbstractContinuousDomain, base class for continuous domains.
Definition: AbstractContinuousDomain.h:47
Base class for all escript domains.
Definition: AbstractDomain.h:51
int StatusType
Definition: AbstractDomain.h:53
Base class for escript system matrices.
Definition: AbstractSystemMatrix.h:44
Give a short description of what AbstractTransportProblem does.
Definition: AbstractTransportProblem.h:45
Data represents a collection of datapoints.
Definition: Data.h:64
Definition: FunctionSpace.h:36
An exception class for features which are not (yet) implemented.
Definition: EsysException.h:90
An exception class that signals an invalid argument value.
Definition: EsysException.h:100
this class holds a (distributed) stiffness matrix
Definition: SystemMatrix.h:50
RipleyDomain extends the AbstractContinuousDomain interface for the Ripley library and is the base cl...
Definition: ripley/src/RipleyDomain.h:103
virtual int getReducedFunctionOnContactOneCode() const
returns a FunctionOnContactOne code with reduced integration order
Definition: ripley/src/RipleyDomain.h:461
virtual bool onMasterProcessor() const
returns true if on MPI processor 0, else false
Definition: ripley/src/RipleyDomain.h:152
virtual dim_t getNumNodes() const =0
returns the number of nodes per MPI rank
virtual void readBinaryGrid(escript::Data &out, std::string filename, const ReaderParameters &params) const =0
reads grid data from a raw binary file into a Data object
virtual boost::python::tuple getGridParameters() const =0
returns the tuple (origin, spacing, number_of_elements)
std::vector< DiracPoint > m_diracPoints
Definition: ripley/src/RipleyDomain.h:778
MPI_Comm getMPIComm() const
returns the MPI communicator
Definition: ripley/src/RipleyDomain.h:158
void addToSystemMatrix(escript::AbstractSystemMatrix *mat, const IndexVector &nodes, dim_t numEq, const std::vector< Scalar > &array) const
virtual void interpolateNodesOnFaces(escript::Data &out, const escript::Data &in, bool reduced) const =0
interpolates data on nodes in 'in' onto (reduced) face elements in 'out'
virtual bool operator!=(const escript::AbstractDomain &other) const
inequality operator
Definition: ripley/src/RipleyDomain.h:187
virtual const int * getNumSubdivisionsPerDim() const =0
returns the number of spatial subdivisions in each dimension
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:213
virtual void assembleGradient(escript::Data &out, const escript::Data &in) const =0
computes the gradient of 'in' and puts the result in 'out'
virtual RankVector getOwnerVector(int fsType) const =0
returns a vector of rank numbers where vec[i]=n means that rank n 'owns' element/face element i.
virtual void setToNormal(escript::Data &out) const =0
copies the surface normals at data points into out. The actual function space to be considered is def...
virtual std::string getDescription() const =0
returns a description for this domain
virtual int getTag(const std::string &name) const
returns the tag key for tag name
Definition: ripley/src/RipleyDomain.h:222
virtual int getContinuousFunctionCode() const
returns a continuous FunctionSpace code
Definition: ripley/src/RipleyDomain.h:400
virtual int getFunctionCode() const
returns a function FunctionSpace code
Definition: ripley/src/RipleyDomain.h:412
virtual std::vector< IndexVector > getConnections(bool includeShared) const =0
returns occupied matrix column indices for all matrix rows
std::vector< int > m_nodeTags
Definition: ripley/src/RipleyDomain.h:775
escript::JMPI m_mpiInfo
Definition: ripley/src/RipleyDomain.h:773
virtual void assembleCoordinates(escript::Data &arg) const =0
populates the data object 'arg' with the node coordinates
virtual void interpolateNodesOnElements(escript::Data &out, const escript::Data &in, bool reduced) const =0
interpolates data on nodes in 'in' onto (reduced) elements in 'out'
virtual double getLocalCoordinate(index_t index, int dim) const =0
returns the index'th coordinate value in given dimension for this rank
virtual int getReducedFunctionOnBoundaryCode() const
returns a function on boundary with reduced integration order FunctionSpace code
Definition: ripley/src/RipleyDomain.h:431
virtual bool isValidTagName(const std::string &name) const
returns true if name is a defined tag name
Definition: ripley/src/RipleyDomain.h:235
virtual void readBinaryGridFromZipped(escript::Data &out, std::string filename, const ReaderParameters &params) const =0
reads grid data from a compressed raw binary file into a Data object
virtual int getSolutionCode() const
returns a Solution FunctionSpace code
Definition: ripley/src/RipleyDomain.h:469
virtual const double * getElementLength() const =0
returns the lengths of an element
virtual bool supportsContactElements() const
returns true if this domain supports contact elements, false otherwise
Definition: ripley/src/RipleyDomain.h:394
void dump(const std::string &filename) const =0
dumps the mesh to a file with the given name
virtual dim_t findNode(const double *coords) const =0
finds the node that the given point coordinates belong to
virtual void setToSize(escript::Data &out) const =0
copies the size of samples into out. The actual function space to be considered is defined by out....
virtual void writeBinaryGrid(const escript::Data &in, std::string filename, int byteOrder, int dataType) const =0
writes a Data object to a file in raw binary format
virtual void nodesToDOF(escript::Data &out, const escript::Data &in) const =0
converts data on nodes in 'in' to degrees of freedom in 'out'
IndexVector m_diracPointNodeIDs
Definition: ripley/src/RipleyDomain.h:779
int m_numDim
Definition: ripley/src/RipleyDomain.h:771
virtual Assembler_ptr createAssembler(std::string type, const DataMap &options) const
Definition: ripley/src/RipleyDomain.h:758
StatusType m_status
Definition: ripley/src/RipleyDomain.h:772
virtual bool ownSample(int fsType, index_t id) const =0
returns true if this rank owns the sample id on given function space
virtual int getMPISize() const
returns the number of processors used for this domain
Definition: ripley/src/RipleyDomain.h:130
assembler_t assembler_type
Definition: ripley/src/RipleyDomain.h:780
virtual const dim_t * getNumNodesPerDim() const =0
returns the number of nodes per MPI rank in each dimension
virtual int getApproximationOrder(int fsType) const
returns the approximation order used for a function space
Definition: ripley/src/RipleyDomain.h:388
virtual const dim_t * getNumElementsPerDim() const =0
returns the number of elements per MPI rank in each dimension
virtual int getMPIRank() const
returns the MPI rank of this processor
Definition: ripley/src/RipleyDomain.h:136
std::vector< int > m_elementTags
Definition: ripley/src/RipleyDomain.h:776
virtual dim_t getNumFaceElements() const =0
returns the number of face elements on current MPI rank
virtual void MPIBarrier() const
if compiled for MPI then executes an MPI_Barrier, else does nothing
Definition: ripley/src/RipleyDomain.h:142
static DecompositionPolicy m_decompPolicy
Definition: ripley/src/RipleyDomain.h:891
virtual dim_t getNumDataPointsGlobal() const =0
returns the number of data points summed across all MPI processes
std::vector< int > m_faceTags
Definition: ripley/src/RipleyDomain.h:777
virtual int getReducedContinuousFunctionCode() const
returns a continuous on reduced order nodes FunctionSpace code
Definition: ripley/src/RipleyDomain.h:406
virtual int getFunctionOnContactOneCode() const
returns a FunctionOnContactOne code
Definition: ripley/src/RipleyDomain.h:453
virtual dim_t getNumDOF() const =0
returns the number of degrees of freedom per MPI rank
virtual int getReducedFunctionOnContactZeroCode() const
returns a FunctionOnContactZero code with reduced integration order
Definition: ripley/src/RipleyDomain.h:445
virtual void assembleIntegrate(std::vector< cplx_t > &integrals, const escript::Data &arg) const =0
virtual IndexVector getNodeDistribution() const =0
returns the node distribution vector
virtual IndexVector getDiagonalIndices(bool upperOnly) const =0
returns the indices of the occupied matrix diagonals
virtual int getDiracDeltaFunctionsCode() const
returns a DiracDeltaFunctions FunctionSpace code
Definition: ripley/src/RipleyDomain.h:481
virtual int getFunctionOnBoundaryCode() const
returns a function on boundary FunctionSpace code
Definition: ripley/src/RipleyDomain.h:424
virtual void assembleIntegrate(std::vector< real_t > &integrals, const escript::Data &arg) const =0
copies the integrals of the function defined by 'arg' into 'integrals'
virtual int getFunctionOnContactZeroCode() const
return a FunctionOnContactZero code
Definition: ripley/src/RipleyDomain.h:437
virtual int getReducedFunctionCode() const
returns a function with reduced integration order FunctionSpace code
Definition: ripley/src/RipleyDomain.h:418
virtual const double * getLength() const =0
returns the lengths of the domain
virtual int getReducedSolutionCode() const
returns a ReducedSolution FunctionSpace code
Definition: ripley/src/RipleyDomain.h:475
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:364
virtual void write(const std::string &filename) const =0
writes the current mesh to a file with the given name
virtual dim_t getDofOfNode(dim_t node) const =0
TagMap m_tagMap
Definition: ripley/src/RipleyDomain.h:774
virtual const dim_t * getNumFacesPerBoundary() const =0
returns the number of face elements in the order (left,right,bottom,top,[front,back]) on current MPI ...
const dim_t * borrowSampleReferenceIDs(int fsType) const =0
returns the array of reference numbers for a function space type
virtual void readNcGrid(escript::Data &out, std::string filename, std::string varname, const ReaderParameters &params) const =0
reads grid data from a netCDF file into a Data object
virtual dim_t getNumElements() const =0
returns the number of elements per MPI rank
virtual escript::JMPI getMPI() const
returns a reference to the MPI information wrapper for this domain
Definition: ripley/src/RipleyDomain.h:124
virtual int getDim() const
returns the number of spatial dimensions of the domain
Definition: ripley/src/RipleyDomain.h:177
index_t dim_t
Definition: DataTypes.h:66
int index_t
type for array/matrix indices used both globally and on each rank
Definition: DataTypes.h:61
boost::shared_ptr< SubWorld > SubWorld_ptr
Definition: SubWorld.h:147
boost::shared_ptr< AbstractTransportProblem > ATP_ptr
Definition: AbstractTransportProblem.h:161
boost::shared_ptr< AbstractSystemMatrix > ASM_ptr
Definition: AbstractSystemMatrix.h:33
boost::shared_ptr< JMPI_ > JMPI
Definition: EsysMPI.h:74
static dim_t N
Definition: SparseMatrix_saveHB.cpp:37
boost::shared_ptr< SystemMatrixPattern > SystemMatrixPattern_ptr
Definition: SystemMatrixPattern.h:40
boost::shared_ptr< Pattern > Pattern_ptr
Definition: Pattern.h:39
boost::shared_ptr< Connector > Connector_ptr
Definition: Coupler.h:38
static dim_t M
Definition: SparseMatrix_saveHB.cpp:37
Definition: ripley/src/AbstractAssembler.h:26
std::vector< index_t > IndexVector
Definition: Ripley.h:44
void setDecompositionPolicy(int policy)
Definition: ripleycpp.cpp:140
@ ReducedElements
Definition: Ripley.h:55
@ Points
Definition: Ripley.h:58
@ ReducedNodes
Definition: Ripley.h:53
@ FaceElements
Definition: Ripley.h:56
@ ReducedFaceElements
Definition: Ripley.h:57
@ DegreesOfFreedom
Definition: Ripley.h:50
@ Nodes
Definition: Ripley.h:52
@ Elements
Definition: Ripley.h:54
@ ReducedDegreesOfFreedom
Definition: Ripley.h:51
SystemMatrixType
Definition: ripley/src/RipleyDomain.h:50
@ SMT_CUSP
Definition: ripley/src/RipleyDomain.h:52
@ SMT_PASO
Definition: ripley/src/RipleyDomain.h:51
@ SMT_SYMMETRIC
Definition: ripley/src/RipleyDomain.h:54
@ SMT_UNROLL
Definition: ripley/src/RipleyDomain.h:56
@ SMT_TRILINOS
Definition: ripley/src/RipleyDomain.h:53
@ SMT_COMPLEX
Definition: ripley/src/RipleyDomain.h:55
DecompositionPolicy
Definition: ripley/src/RipleyDomain.h:59
@ DECOMP_ADD_ELEMENTS
Definition: ripley/src/RipleyDomain.h:60
@ DECOMP_EXPAND
Definition: ripley/src/RipleyDomain.h:61
@ DECOMP_STRICT
Definition: ripley/src/RipleyDomain.h:62
assembler_t
Definition: ripley/src/RipleyDomain.h:44
@ LAME_ASSEMBLER
Definition: ripley/src/RipleyDomain.h:47
@ WAVE_ASSEMBLER
Definition: ripley/src/RipleyDomain.h:46
@ DEFAULT_ASSEMBLER
Definition: ripley/src/RipleyDomain.h:45
std::map< std::string, int > TagMap
Definition: Ripley.h:47
std::vector< int > RankVector
Definition: Ripley.h:46
std::map< std::string, escript::Data > DataMap
Definition: ripley/src/domainhelpers.h:25
bool probeInterpolationAcross(int fsType_source, const escript::AbstractDomain &domain, int fsType_target, int dim)
Definition: CrossDomainCoupler.cpp:32
#define RIPLEY_DLL_API
Definition: ripley/src/system_dep.h:21
A struct to contain a dirac point's information.
Definition: ripley/src/RipleyDomain.h:91
int tag
Definition: ripley/src/RipleyDomain.h:93
dim_t node
Definition: ripley/src/RipleyDomain.h:92
Structure that wraps parameters for the grid reading routines.
Definition: ripley/src/RipleyDomain.h:70
int byteOrder
byte order in the file (used by binary reader only)
Definition: ripley/src/RipleyDomain.h:81
std::vector< dim_t > first
the (global) offset into the data object to start writing into
Definition: ripley/src/RipleyDomain.h:72
std::vector< int > multiplier
Definition: ripley/src/RipleyDomain.h:77
int dataType
data type in the file (used by binary reader only)
Definition: ripley/src/RipleyDomain.h:83
std::vector< int > reverse
if non-zero, values are written from last index to first index
Definition: ripley/src/RipleyDomain.h:79
std::vector< dim_t > numValues
the number of values to read from file
Definition: ripley/src/RipleyDomain.h:74