18 #ifndef __Speckley_DOMAIN_H__
19 #define __Speckley_DOMAIN_H__
21 #include <speckley/Speckley.h>
22 #include <speckley/SpeckleyException.h>
23 #include <speckley/AbstractAssembler.h>
24 #include <speckley/domainhelpers.h>
26 #include <escript/AbstractContinuousDomain.h>
27 #include <escript/Data.h>
28 #include <escript/FunctionSpace.h>
29 #include <escript/SubWorld.h>
31 #include <boost/python/tuple.hpp>
32 #include <boost/python/list.hpp>
45 typedef std::map<std::string, int>
simap_t;
124 MPI_Barrier(m_mpiInfo->comm);
145 virtual bool isValidFunctionSpaceType(
int fsType)
const;
151 virtual std::string functionSpaceTypeAsString(
int fsType)
const;
157 virtual int getDim()
const {
return m_numDim; }
168 return !(operator==(other));
177 virtual std::pair<int,dim_t> getDataShape(
int fsType)
const;
185 int getTagFromSampleNo(
int fsType,
dim_t sampleNo)
const;
193 virtual void setTagMap(
const std::string& name,
int tag) {
194 m_tagMap[name] = tag;
202 virtual int getTag(
const std::string& name)
const {
203 if (m_tagMap.find(name) != m_tagMap.end()) {
204 return m_tagMap.find(name)->second;
216 return (m_tagMap.find(name)!=m_tagMap.end());
223 virtual std::string showTagNames()
const;
245 virtual bool probeInterpolationOnDomain(
int fsType_source,
246 int fsType_target)
const;
255 virtual signed char preferredInterpolationOnDomain(
int fsType_source,
256 int fsType_target)
const;
265 commonFunctionSpace(
const std::vector<int>& fs,
int& resultcode)
const;
288 #ifdef ESYS_HAVE_BOOST_NUMPY
293 virtual boost::python::numpy::ndarray getNumpyX()
const;
298 virtual boost::python::numpy::ndarray getConnectivityInfo()
const;
333 virtual void setTags(
int fsType,
int newTag,
const escript::Data& mask)
const;
340 virtual bool isCellOriented(
int fsType)
const;
354 virtual int getNumberOfTagsInUse(
int fsType)
const;
360 virtual const int* borrowListOfTagsInUse(
int fsType)
const;
366 virtual bool canTag(
int fsType)
const;
487 virtual int getSystemMatrixTypeId(
const boost::python::object& options)
const;
498 virtual int getTransportTypeId(
int solver,
int preconditioner,
int package,
499 bool symmetry)
const;
506 virtual void setToIntegrals(std::vector<real_t>& integrals,
508 virtual void setToIntegrals(std::vector<cplx_t>& integrals,
519 Assembler_ptr assembler)
const;
527 Assembler_ptr assembler)
const;
535 Assembler_ptr assembler)
const;
542 const boost::python::list& data,
543 Assembler_ptr assembler)
const;
552 Assembler_ptr assembler)
const;
557 void addPDEToTransportProblemFromPython(
560 Assembler_ptr assembler)
const;
568 int column_blocksize,
584 virtual void Print_Mesh_Info(
bool full=
false)
const;
594 virtual void write(
const std::string& filename)
const = 0;
607 void dump(
const std::string& filename)
const = 0;
657 int byteOrder,
int dataType)
const = 0;
719 virtual bool supportsFilter(
const boost::python::tuple& t)
const;
725 const DataMap& options)
const {
729 Assembler_ptr createAssemblerFromPython(
const std::string type,
730 const boost::python::list& options)
const;
758 template<
typename Scalar>
762 void updateTagsInUse(
int fsType)
const;
768 void addPoints(
const std::vector<double>& coords,
769 const std::vector<int>& tags);
772 template<
typename Scalar>
802 bool reduced)
const = 0;
815 virtual void balanceNeighbours(
escript::Data& data,
bool average)
const = 0;
821 const DataMap& coefs, Assembler_ptr assembler)
const;
827 Assembler_ptr assembler)
const;
831 Assembler_ptr assembler)
const;
835 Assembler_ptr assembler)
const;
838 Assembler_ptr assembler)
const;
840 template<
typename Scalar>
841 void setToIntegralsWorker(std::vector<Scalar>& integrals,
int MPI_Comm
Definition: EsysMPI.h:44
AbstractContinuousDomain, base class for continuous domains.
Definition: AbstractContinuousDomain.h:47
virtual void addPDEToTransportProblem(AbstractTransportProblem &tp, escript::Data &source, const escript::Data &M, const escript::Data &A, const escript::Data &B, const escript::Data &C, const escript::Data &D, const escript::Data &X, const escript::Data &Y, const escript::Data &d, const escript::Data &y, const escript::Data &d_contact, const escript::Data &y_contact, const escript::Data &d_dirac, const escript::Data &y_dirac) const
adds a PDE onto a transport problem
Definition: AbstractContinuousDomain.cpp:166
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
SpeckleyDomain extends the AbstractContinuousDomain interface for the Speckley library and is the bas...
Definition: speckley/src/SpeckleyDomain.h:86
std::vector< int > m_nodeTags
Definition: speckley/src/SpeckleyDomain.h:749
virtual escript::JMPI getMPI() const
returns a reference to the MPI information wrapper for this domain
Definition: speckley/src/SpeckleyDomain.h:104
virtual int getReducedFunctionOnBoundaryCode() const
returns a function on boundary with reduced integration order FunctionSpace code
Definition: speckley/src/SpeckleyDomain.h:421
assembler_t assembler_type
Definition: speckley/src/SpeckleyDomain.h:753
virtual StatusType getStatus() const
returns a status indicator of the domain. The status identifier should be unique over the lifetime of...
Definition: speckley/src/SpeckleyDomain.h:348
virtual IndexVector getNodeDistribution() const =0
returns the node distribution vector
virtual int getReducedFunctionCode() const
returns a function with reduced integration order FunctionSpace code
Definition: speckley/src/SpeckleyDomain.h:404
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 void write(const std::string &filename) const =0
writes the current mesh to a file with the given name
virtual void interpolateNodesOnElements(escript::Data &out, const escript::Data &in, bool reduced) const =0
interpolates data on nodes in 'in' onto elements in '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 reduceElements(escript::Data &out, const escript::Data &in) const =0
interpolates from Element -> ReducedElement
std::vector< int > m_elementTags
Definition: speckley/src/SpeckleyDomain.h:750
virtual int getFunctionOnBoundaryCode() const
returns a function on boundary FunctionSpace code
Definition: speckley/src/SpeckleyDomain.h:412
virtual void MPIBarrier() const
if compiled for MPI then executes an MPI_Barrier, else does nothing
Definition: speckley/src/SpeckleyDomain.h:122
virtual dim_t getNumElements() const =0
returns the number of elements per MPI rank
virtual bool ownSample(int fsType, index_t id) const =0
returns true if this rank owns the sample id on given function space
virtual dim_t getNumDataPointsGlobal() const =0
returns the number of data points summed across all MPI processes
virtual const int * getNumSubdivisionsPerDim() const =0
returns the number of spatial subdivisions in each dimension
std::vector< DiracPoint > m_diracPoints
Definition: speckley/src/SpeckleyDomain.h:751
virtual int getMPIRank() const
returns the MPI rank of this processor
Definition: speckley/src/SpeckleyDomain.h:116
virtual int getReducedContinuousFunctionCode() const
returns a continuous on reduced order nodes FunctionSpace code
Definition: speckley/src/SpeckleyDomain.h:390
int getOrder() const
returns the order of the domain
Definition: speckley/src/SpeckleyDomain.h:742
virtual int getReducedSolutionCode() const
returns a ReducedSolution FunctionSpace code
Definition: speckley/src/SpeckleyDomain.h:469
int m_numDim
Definition: speckley/src/SpeckleyDomain.h:745
virtual double getLocalCoordinate(dim_t index, int dim) const =0
returns the index'th coordinate value in given dimension for this rank
virtual const dim_t * getNumElementsPerDim() const =0
returns the number of elements per MPI rank in each dimension
virtual bool isValidTagName(const std::string &name) const
returns true if name is a defined tag name
Definition: speckley/src/SpeckleyDomain.h:215
virtual dim_t getNumDOF() const =0
returns the number of degrees of freedom per MPI rank
virtual int getFunctionCode() const
returns a function FunctionSpace code
Definition: speckley/src/SpeckleyDomain.h:398
virtual bool supportsContactElements() const
returns true if this domain supports contact elements, false otherwise
Definition: speckley/src/SpeckleyDomain.h:378
virtual int getMPISize() const
returns the number of processors used for this domain
Definition: speckley/src/SpeckleyDomain.h:110
virtual int getReducedFunctionOnContactOneCode() const
returns a FunctionOnContactOne code with reduced integration order
Definition: speckley/src/SpeckleyDomain.h:453
virtual void readNcGrid(escript::Data &out, std::string filename, std::string varname, const ReaderParameters ¶ms) const =0
reads grid data from a netCDF file into a Data object
virtual boost::python::tuple getGridParameters() const =0
returns the tuple (origin, spacing, number_of_elements)
TagMap m_tagMap
Definition: speckley/src/SpeckleyDomain.h:748
virtual void assembleCoordinates(escript::Data &arg) const =0
populates the data object 'arg' with the node coordinates
StatusType m_status
Definition: speckley/src/SpeckleyDomain.h:746
const index_t * borrowSampleReferenceIDs(int fsType) const =0
returns the array of reference numbers for a function space type
virtual std::string getDescription() const =0
returns a description for this domain
virtual void readBinaryGridFromZipped(escript::Data &out, std::string filename, const ReaderParameters ¶ms) const =0
reads grid data from a compressed raw binary file into a Data object
IndexVector m_diracPointNodeIDs
Definition: speckley/src/SpeckleyDomain.h:752
virtual const dim_t * getNumNodesPerDim() const =0
returns the number of nodes per MPI rank in each dimension
virtual int getDim() const
returns the number of spatial dimensions of the domain
Definition: speckley/src/SpeckleyDomain.h:157
virtual int getReducedFunctionOnContactZeroCode() const
returns a FunctionOnContactZero code with reduced integration order
Definition: speckley/src/SpeckleyDomain.h:437
virtual int getApproximationOrder(int fsType) const
returns the approximation order used for a function space
Definition: speckley/src/SpeckleyDomain.h:372
virtual bool onMasterProcessor() const
returns true if on MPI processor 0, else false
Definition: speckley/src/SpeckleyDomain.h:132
escript::JMPI m_mpiInfo
Definition: speckley/src/SpeckleyDomain.h:747
virtual void interpolateElementsOnNodes(escript::Data &out, const escript::Data &in) const =0
interpolates data on elements in 'in' onto nodes in 'out'
MPI_Comm getMPIComm() const
returns the MPI communicator
Definition: speckley/src/SpeckleyDomain.h:138
virtual dim_t getNumNodes() const =0
returns the number of nodes per MPI rank
virtual bool probeInterpolationAcross(int, const escript::AbstractDomain &, int) const =0
determines whether interpolation from source to target is possible
virtual int getContinuousFunctionCode() const
returns a continuous FunctionSpace code
Definition: speckley/src/SpeckleyDomain.h:384
virtual int getDiracDeltaFunctionsCode() const
returns a DiracDeltaFunctions FunctionSpace code
Definition: speckley/src/SpeckleyDomain.h:477
virtual int getSolutionCode() const
returns a Solution FunctionSpace code
Definition: speckley/src/SpeckleyDomain.h:461
int m_order
element order (will be m_order + 1 quad points in each axis)
Definition: speckley/src/SpeckleyDomain.h:755
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 int getTag(const std::string &name) const
returns the tag key for tag name
Definition: speckley/src/SpeckleyDomain.h:202
virtual void interpolateAcross(escript::Data &target, const escript::Data &source) const =0
interpolates data given on source onto target where source and target are given on different domains
virtual bool operator!=(const escript::AbstractDomain &other) const
inequality operator
Definition: speckley/src/SpeckleyDomain.h:167
virtual dim_t findNode(const double *coords) const =0
finds the node that the given point coordinates belong to
virtual void setTagMap(const std::string &name, int tag)
sets a map from a clear tag name to a tag key
Definition: speckley/src/SpeckleyDomain.h:193
virtual Assembler_ptr createAssembler(const std::string type, const DataMap &options) const
Definition: speckley/src/SpeckleyDomain.h:724
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 int getFunctionOnContactZeroCode() const
return a FunctionOnContactZero code
Definition: speckley/src/SpeckleyDomain.h:429
virtual void assembleIntegrate(std::vector< cplx_t > &integrals, const escript::Data &arg) const =0
virtual int getFunctionOnContactOneCode() const
returns a FunctionOnContactOne code
Definition: speckley/src/SpeckleyDomain.h:445
virtual void assembleGradient(escript::Data &out, const escript::Data &in) const =0
computes the gradient of 'in' and puts the result in 'out'
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 ...
void dump(const std::string &filename) const =0
dumps the mesh to a file with the given name
virtual const double * getLength() const =0
returns the lengths of the domain
virtual dim_t getDofOfNode(dim_t node) const =0
virtual void readBinaryGrid(escript::Data &out, std::string filename, const ReaderParameters ¶ms) const =0
reads grid data from a raw binary file into a Data object
SpeckleyException exception class.
Definition: SpeckleyException.h:31
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
Definition: AbstractAssembler.cpp:19
std::map< std::string, int > simap_t
Definition: speckley/src/SpeckleyDomain.h:45
std::map< std::string, int > TagMap
Definition: Speckley.h:46
@ Points
Definition: Speckley.h:67
@ Elements
Definition: Speckley.h:63
@ Nodes
Definition: Speckley.h:61
@ DegreesOfFreedom
Definition: Speckley.h:59
@ ReducedElements
Definition: Speckley.h:64
std::vector< real_t > DoubleVector
Definition: Speckley.h:45
assembler_t
Definition: speckley/src/SpeckleyDomain.h:36
@ DEFAULT_ASSEMBLER
Definition: speckley/src/SpeckleyDomain.h:37
std::vector< index_t > IndexVector
Definition: Speckley.h:44
std::map< std::string, escript::Data > DataMap
Definition: speckley/src/domainhelpers.h:25
#define Speckley_DLL_API
Definition: speckley/src/system_dep.h:23
A struct to contain a dirac point's information.
Definition: speckley/src/SpeckleyDomain.h:74
dim_t node
Definition: speckley/src/SpeckleyDomain.h:75
int tag
Definition: speckley/src/SpeckleyDomain.h:76
Structure that wraps parameters for the grid reading routines.
Definition: speckley/src/SpeckleyDomain.h:53
int byteOrder
byte order in the file (used by binary reader only)
Definition: speckley/src/SpeckleyDomain.h:64
std::vector< int > multiplier
Definition: speckley/src/SpeckleyDomain.h:60
std::vector< dim_t > numValues
the number of values to read from file
Definition: speckley/src/SpeckleyDomain.h:57
std::vector< dim_t > first
the (global) offset into the data object to start writing into
Definition: speckley/src/SpeckleyDomain.h:55
int dataType
data type in the file (used by binary reader only)
Definition: speckley/src/SpeckleyDomain.h:66
std::vector< int > reverse
if non-zero, values are written from last index to first index
Definition: speckley/src/SpeckleyDomain.h:62