17 #ifndef __RIPLEY_DOMAIN_H__
18 #define __RIPLEY_DOMAIN_H__
20 #ifdef BADPYTHONMACROS
27 #undef BADPYTHONMACROS
30 #include <boost/python/tuple.hpp>
31 #include <boost/python/list.hpp>
33 #include <ripley/Ripley.h>
34 #include <ripley/RipleyException.h>
35 #include <ripley/AbstractAssembler.h>
37 #include <escript/AbstractContinuousDomain.h>
38 #include <escript/Data.h>
39 #include <escript/FunctionSpace.h>
41 #include <paso/SystemMatrix.h>
56 typedef std::map<std::string, int>
simap_t;
129 MPI_Barrier(m_mpiInfo->comm);
150 return m_mpiInfo->comm;
161 virtual bool isValidFunctionSpaceType(
int fsType)
const;
167 virtual std::string functionSpaceTypeAsString(
int fsType)
const;
173 virtual int getDim()
const {
return m_numDim; }
184 return !(operator==(other));
193 virtual std::pair<int,int> getDataShape(
int fsType)
const;
201 int getTagFromSampleNo(
int fsType,
int sampleNo)
const;
209 virtual void setTagMap(
const std::string& name,
int tag) {
210 m_tagMap[name] = tag;
218 virtual int getTag(
const std::string& name)
const {
219 if (m_tagMap.find(name) != m_tagMap.end()) {
220 return m_tagMap.find(name)->second;
232 return (m_tagMap.find(name)!=m_tagMap.end());
239 virtual std::string showTagNames()
const;
260 virtual bool probeInterpolationOnDomain(
int fsType_source,
int fsType_target)
const;
269 virtual signed char preferredInterpolationOnDomain(
int fsType_source,
int fsType_target)
const;
278 commonFunctionSpace(
const std::vector<int>& fs,
int& resultcode)
const;
330 virtual void setTags(
const int fsType,
const int newTag,
const escript::Data& mask)
const;
337 virtual bool isCellOriented(
int fsType)
const;
351 virtual int getNumberOfTagsInUse(
int fsType)
const;
357 virtual const int* borrowListOfTagsInUse(
int fsType)
const;
363 virtual bool canTag(
int fsType)
const;
474 virtual int getSystemMatrixTypeId(
const int solver,
const int preconditioner,
const int package,
const bool symmetry)
const;
485 virtual int getTransportTypeId(
const int solver,
const int preconditioner,
const int package,
const bool symmetry)
const;
512 escript::Data& rhs,std::map<std::string, escript::Data> data)
const;
535 std::map<std::string, escript::Data> data)
const;
542 boost::python::list data)
const;
564 const int column_blocksize,
573 const int type)
const;
580 virtual void Print_Mesh_Info(
const bool full=
false)
const;
596 virtual std::string getDescription()
const = 0;
603 void dump(
const std::string& filename)
const = 0;
610 const int* borrowSampleReferenceIDs(
int fsType)
const = 0;
637 virtual void readBinaryGridFromZipped(
escript::Data& out, std::string filename,
643 virtual void writeBinaryGrid(
const escript::Data& in, std::string filename,
644 int byteOrder,
int dataType)
const = 0;
650 virtual bool ownSample(
int fsType,
index_t id)
const = 0;
656 virtual int getNumDataPointsGlobal()
const = 0;
662 virtual const int* getNumNodesPerDim()
const = 0;
668 virtual const int* getNumElementsPerDim()
const = 0;
675 virtual const int* getNumFacesPerBoundary()
const = 0;
681 virtual IndexVector getNodeDistribution()
const = 0;
687 virtual const int* getNumSubdivisionsPerDim()
const = 0;
693 virtual double getLocalCoordinate(
int index,
int dim)
const = 0;
699 virtual boost::python::tuple getGridParameters()
const = 0;
706 virtual bool supportsFilter(
const boost::python::tuple& t)
const;
709 std::map<std::string, escript::Data> options) {
713 void setAssemblerFromPython(std::string type,
714 boost::python::list options);
738 void updateTagsInUse(
int fsType)
const;
750 void createCouplePatterns(
const std::vector<IndexVector>& colIndices,
751 const std::vector<IndexVector>& rowIndices,
760 void addPoints(
int numPoints,
const double* points_ptr,
761 const int* tags_ptr);
766 virtual dim_t getNumNodes()
const = 0;
769 virtual dim_t getNumElements()
const = 0;
772 virtual dim_t getNumDOF()
const = 0;
775 virtual dim_t getNumFaceElements()
const = 0;
782 virtual void assembleCoordinates(
escript::Data& arg)
const = 0;
792 bool reducedColOrder)
const = 0;
797 bool reduced)
const = 0;
802 bool reduced)
const = 0;
810 virtual int getDofOfNode(
int node)
const = 0;
815 std::map<std::string, escript::Data> coefs)
const;
820 std::map<std::string, escript::Data> coefs)
const;
823 std::map<std::string, escript::Data> coefs)
const;
826 virtual int findNode(
const double *coords)
const = 0;
831 #endif // __RIPLEY_DOMAIN_H__
AbstractContinuousDomain, base class for continuous domains.
Definition: AbstractContinuousDomain.h:46
Give a short description of what FunctionSpace does.
Definition: FunctionSpace.h:46
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:77
virtual void setAssembler(std::string type, std::map< std::string, escript::Data > options)
Definition: ripley/src/RipleyDomain.h:708
std::vector< struct DiracPoint > m_diracPoints
Definition: ripley/src/RipleyDomain.h:724
IndexVector m_faceTagsInUse
Definition: ripley/src/RipleyDomain.h:722
virtual bool onMasterProcessor() const
returns true if on MPI processor 0, else false
Definition: ripley/src/RipleyDomain.h:137
virtual int getFunctionOnBoundaryCode() const
returns a function on boundary FunctionSpace code
Definition: ripley/src/RipleyDomain.h:405
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:109
virtual int getFunctionOnContactOneCode() const
returns a FunctionOnContactOne code
Definition: ripley/src/RipleyDomain.h:434
unsigned int getMPIComm() const
returns the MPI communicator
Definition: ripley/src/RipleyDomain.h:148
virtual bool operator!=(const escript::AbstractDomain &other) const
inequality operator
Definition: ripley/src/RipleyDomain.h:183
Definition: ripley/src/RipleyDomain.h:47
Definition: AbstractAssembler.h:35
int node
Definition: ripley/src/RipleyDomain.h:86
Struct that holds MPI communicator, rank, size and a tag counter.
Definition: Esys_MPI.h:48
int byteOrder
byte order in the file (used by binary reader only)
Definition: ripley/src/RipleyDomain.h:75
A struct to contain a dirac point's information.
Definition: ripley/src/RipleyDomain.h:84
Definition: ripley/src/RipleyDomain.h:46
virtual int getReducedContinuousFunctionCode() const
returns a continuous on reduced order nodes FunctionSpace code
Definition: ripley/src/RipleyDomain.h:387
assembler_t assembler_type
Definition: ripley/src/RipleyDomain.h:726
RipleyException exception class.
Definition: RipleyException.h:29
Esys_MPIInfo * m_mpiInfo
Definition: ripley/src/RipleyDomain.h:718
static dim_t M
Definition: SparseMatrix_saveHB.cpp:36
virtual bool supportsContactElements() const
returns true if this domain supports contact elements, false otherwise
Definition: ripley/src/RipleyDomain.h:375
std::vector< int > numValues
the number of values to read from file
Definition: ripley/src/RipleyDomain.h:68
std::vector< int > reverse
if non-zero, values are written from last index to first index
Definition: ripley/src/RipleyDomain.h:73
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:57
assembler_t
Definition: ripley/src/RipleyDomain.h:45
virtual int getFunctionOnContactZeroCode() const
return a FunctionOnContactZero code
Definition: ripley/src/RipleyDomain.h:418
virtual int getApproximationOrder(const int fsType) const
returns the approximation order used for a function space
Definition: ripley/src/RipleyDomain.h:369
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:29
virtual int getDim() const
returns the number of spatial dimensions of the domain
Definition: ripley/src/RipleyDomain.h:173
boost::shared_ptr< SystemMatrix > SystemMatrix_ptr
Definition: SystemMatrix.h:38
IndexVector m_diracPointNodeIDs
Definition: ripley/src/RipleyDomain.h:725
virtual int getMPIRank() const
returns the MPI rank of this processor
Definition: ripley/src/RipleyDomain.h:121
std::vector< int > first
the (global) offset into the data object to start writing into
Definition: ripley/src/RipleyDomain.h:66
std::map< std::string, int > simap_t
Definition: ripley/src/RipleyDomain.h:56
std::vector< index_t > IndexVector
Definition: Ripley.h:38
StatusType m_status
Definition: ripley/src/RipleyDomain.h:717
IndexVector m_nodeTagsInUse
Definition: ripley/src/RipleyDomain.h:720
virtual int getReducedFunctionCode() const
returns a function with reduced integration order FunctionSpace code
Definition: ripley/src/RipleyDomain.h:399
Structure that wraps parameters for the grid reading routines.
Definition: ripley/src/RipleyDomain.h:63
std::vector< int > multiplier
Definition: ripley/src/RipleyDomain.h:71
virtual int getReducedFunctionOnContactOneCode() const
returns a FunctionOnContactOne code with reduced integration order
Definition: ripley/src/RipleyDomain.h:442
int StatusType
Definition: AbstractDomain.h:78
virtual int getDiracDeltaFunctionsCode() const
returns a DiracDeltaFunctions FunctionSpace code
Definition: ripley/src/RipleyDomain.h:462
Data represents a collection of datapoints.
Definition: Data.h:71
virtual int getReducedFunctionOnContactZeroCode() const
returns a FunctionOnContactZero code with reduced integration order
Definition: ripley/src/RipleyDomain.h:426
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:209
static dim_t N
Definition: SparseMatrix_saveHB.cpp:36
RipleyDomain extends the AbstractContinuousDomain interface for the Ripley library and is the base cl...
Definition: ripley/src/RipleyDomain.h:96
virtual bool isValidTagName(const std::string &name) const
returns true if name is a defined tag name
Definition: ripley/src/RipleyDomain.h:231
virtual int getReducedFunctionOnBoundaryCode() const
returns a function on boundary with reduced integration order FunctionSpace code
Definition: ripley/src/RipleyDomain.h:412
virtual int getMPISize() const
returns the number of processors used for this domain
Definition: ripley/src/RipleyDomain.h:115
virtual int getTag(const std::string &name) const
returns the tag key for tag name
Definition: ripley/src/RipleyDomain.h:218
int index_t
Definition: types.h:25
AbstractAssembler * assembler
Definition: ripley/src/RipleyDomain.h:723
virtual int getReducedSolutionCode() const
returns a ReducedSolution FunctionSpace code
Definition: ripley/src/RipleyDomain.h:456
virtual int getFunctionCode() const
returns a function FunctionSpace code
Definition: ripley/src/RipleyDomain.h:393
int tag
Definition: ripley/src/RipleyDomain.h:87
Give a short description of what AbstractTransportProblem does.
Definition: AbstractTransportProblem.h:45
Definition: ripley/src/RipleyDomain.h:48
Give a short description of what AbstractSystemMatrix does.
Definition: AbstractSystemMatrix.h:44
IndexVector m_elementTagsInUse
Definition: ripley/src/RipleyDomain.h:721
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:345
boost::shared_ptr< AbstractSystemMatrix > ASM_ptr
Definition: AbstractSystemMatrix.h:201
Definition: AbstractDomain.h:62
TagMap m_tagMap
Definition: ripley/src/RipleyDomain.h:719
virtual int getSolutionCode() const
returns a Solution FunctionSpace code
Definition: ripley/src/RipleyDomain.h:450
std::map< std::string, index_t > TagMap
Definition: Ripley.h:41
dim_t m_numDim
Definition: ripley/src/RipleyDomain.h:716
#define RIPLEY_DLL_API
Definition: ripley/src/system_dep.h:35
int dim_t
Definition: types.h:24
std::vector< double > DoubleVector
Definition: Ripley.h:39
virtual int getContinuousFunctionCode() const
returns a continuous FunctionSpace code
Definition: ripley/src/RipleyDomain.h:381
virtual void MPIBarrier() const
if compiled for MPI then executes an MPI_Barrier, else does nothing
Definition: ripley/src/RipleyDomain.h:127