escript  Revision_
speckley/src/SpeckleyDomain.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 __Speckley_DOMAIN_H__
18 #define __Speckley_DOMAIN_H__
19 
20 #include <speckley/Speckley.h>
21 #include <speckley/SpeckleyException.h>
22 #include <speckley/AbstractAssembler.h>
23 #include <speckley/domainhelpers.h>
24 
25 #include <escript/AbstractContinuousDomain.h>
26 #include <escript/Data.h>
27 #include <escript/FunctionSpace.h>
28 #include <escript/SubWorld.h>
29 
30 #include <boost/python/tuple.hpp>
31 #include <boost/python/list.hpp>
32 
33 namespace speckley {
34 
37 };
38 
39 /* There is no particular significance to this type,
40 It is here as a typedef because a bug in clang++ prevents
41 that compiler from recognising it as a valid part of
42 a constant expression.
43 */
44 typedef std::map<std::string, int> simap_t;
45 
46 
52 {
54  std::vector<dim_t> first;
56  std::vector<dim_t> numValues;
59  std::vector<int> multiplier;
61  std::vector<int> reverse;
63  int byteOrder;
65  int dataType;
66 };
67 
72 struct DiracPoint
73 {
75  int tag;
76 };
77 
85 {
86 public:
92 
97  ~SpeckleyDomain();
98 
103  virtual escript::JMPI getMPI() const { return m_mpiInfo; }
104 
109  virtual int getMPISize() const { return m_mpiInfo->size; }
110 
115  virtual int getMPIRank() const { return m_mpiInfo->rank; }
116 
121  virtual void MPIBarrier() const {
122 #ifdef ESYS_MPI
123  MPI_Barrier(m_mpiInfo->comm);
124 #endif
125  }
126 
131  virtual bool onMasterProcessor() const { return getMPIRank()==0; }
132 
137  MPI_Comm getMPIComm() const { return m_mpiInfo->comm; }
138 
144  virtual bool isValidFunctionSpaceType(int fsType) const;
145 
150  virtual std::string functionSpaceTypeAsString(int fsType) const;
151 
156  virtual int getDim() const { return m_numDim; }
157 
161  virtual bool operator==(const escript::AbstractDomain& other) const;
162 
166  virtual bool operator!=(const escript::AbstractDomain& other) const {
167  return !(operator==(other));
168  }
169 
176  virtual std::pair<int,dim_t> getDataShape(int fsType) const;
177 
184  int getTagFromSampleNo(int fsType, dim_t sampleNo) const;
185 
192  virtual void setTagMap(const std::string& name, int tag) {
193  m_tagMap[name] = tag;
194  }
195 
201  virtual int getTag(const std::string& name) const {
202  if (m_tagMap.find(name) != m_tagMap.end()) {
203  return m_tagMap.find(name)->second;
204  } else {
205  throw SpeckleyException("getTag: invalid tag name");
206  }
207  }
208 
214  virtual bool isValidTagName(const std::string& name) const {
215  return (m_tagMap.find(name)!=m_tagMap.end());
216  }
217 
222  virtual std::string showTagNames() const;
223 
229  virtual void setNewX(const escript::Data& arg);
230 
236  virtual void interpolateOnDomain(escript::Data& target,
237  const escript::Data& source) const;
238 
244  virtual bool probeInterpolationOnDomain(int fsType_source,
245  int fsType_target) const;
246 
254  virtual signed char preferredInterpolationOnDomain(int fsType_source,
255  int fsType_target) const;
256 
263  bool
264  commonFunctionSpace(const std::vector<int>& fs, int& resultcode) const;
265 
271  virtual void interpolateAcross(escript::Data& target,
272  const escript::Data& source) const = 0;
273 
278  virtual bool probeInterpolationAcross(int, const escript::AbstractDomain&,
279  int) const = 0;
280 
285  virtual escript::Data getX() const;
286 
291  virtual escript::Data getNormal() const;
292 
296  virtual escript::Data getSize() const;
297 
303  virtual void setToX(escript::Data& arg) const;
304 
311  virtual void setToGradient(escript::Data& out,
312  const escript::Data& in) const;
313 
319  virtual void setTags(int fsType, int newTag, const escript::Data& mask) const;
320 
326  virtual bool isCellOriented(int fsType) const;
327 
334  virtual StatusType getStatus() const { return m_status; }
335 
340  virtual int getNumberOfTagsInUse(int fsType) const;
341 
346  virtual const int* borrowListOfTagsInUse(int fsType) const;
347 
352  virtual bool canTag(int fsType) const;
353 
358  virtual int getApproximationOrder(int fsType) const { return 1; }
359 
364  virtual bool supportsContactElements() const { return false; }
365 
370  virtual int getContinuousFunctionCode() const { return Nodes; }
371 
376  virtual int getReducedContinuousFunctionCode() const {
377  throw SpeckleyException("Speckley does not support reduced functionspaces");
378  }
379 
384  virtual int getFunctionCode() const { return Elements; }
385 
390  virtual int getReducedFunctionCode() const {
391  return ReducedElements;
392  }
393 
398  virtual int getFunctionOnBoundaryCode() const {
399  throw SpeckleyException("Speckley does not support face elements");
400  }
401 
407  virtual int getReducedFunctionOnBoundaryCode() const {
408  throw SpeckleyException("Speckley does not support face elements");
409  }
410 
415  virtual int getFunctionOnContactZeroCode() const {
416  throw SpeckleyException("Speckley does not support contact elements");
417  }
418 
424  throw SpeckleyException("Speckley does not support contact elements");
425  }
426 
431  virtual int getFunctionOnContactOneCode() const {
432  throw SpeckleyException("Speckley does not support contact elements");
433  }
434 
440  throw SpeckleyException("Speckley does not support contact elements");
441  }
442 
447  virtual int getSolutionCode() const {
448  return DegreesOfFreedom;
449  }
450 
455  virtual int getReducedSolutionCode() const {
456  throw SpeckleyException("Speckley does not support reduced function spaces");
457  }
458 
463  virtual int getDiracDeltaFunctionsCode() const { return Points; }
464 
473  virtual int getSystemMatrixTypeId(const boost::python::object& options) const;
474 
484  virtual int getTransportTypeId(int solver, int preconditioner, int package,
485  bool symmetry) const;
486 
492  virtual void setToIntegrals(DoubleVector& integrals, const escript::Data& arg) const;
493 
494 
500  virtual void addToSystem(escript::AbstractSystemMatrix& mat,
501  escript::Data& rhs, const DataMap& data,
502  Assembler_ptr assembler) const;
503 
508  virtual void addToSystemFromPython(escript::AbstractSystemMatrix& mat,
509  escript::Data& rhs, const boost::python::list& data,
510  Assembler_ptr assembler) const;
511 
517  virtual void addToRHS(escript::Data& rhs, const DataMap& data,
518  Assembler_ptr assembler) const;
519 
524  virtual void addToRHSFromPython(escript::Data& rhs,
525  const boost::python::list& data,
526  Assembler_ptr assembler) const;
527 
533  virtual void addPDEToTransportProblem(escript::AbstractTransportProblem& tp,
534  escript::Data& source, const DataMap& data,
535  Assembler_ptr assembler) const;
540  void addPDEToTransportProblemFromPython(
542  escript::Data& source, const boost::python::list& data,
543  Assembler_ptr assembler) const;
544 
549  virtual escript::ASM_ptr newSystemMatrix(int row_blocksize,
550  const escript::FunctionSpace& row_functionspace,
551  int column_blocksize,
552  const escript::FunctionSpace& column_functionspace, int type) const;
553 
558  virtual escript::ATP_ptr newTransportProblem(int blocksize,
559  const escript::FunctionSpace& functionspace,
560  int type) const;
561 
567  virtual void Print_Mesh_Info(bool full=false) const;
568 
569 
570  /************************************************************************/
571 
577  virtual void write(const std::string& filename) const = 0;
578 
583  virtual std::string getDescription() const = 0;
584 
590  void dump(const std::string& filename) const = 0;
591 
597  const index_t* borrowSampleReferenceIDs(int fsType) const = 0;
598 
605  virtual void setToNormal(escript::Data& out) const = 0;
606 
612  virtual void setToSize(escript::Data& out) const = 0;
613 
618  virtual void readNcGrid(escript::Data& out, std::string filename,
619  std::string varname, const ReaderParameters& params) const = 0;
620 
625  virtual void readBinaryGrid(escript::Data& out, std::string filename,
626  const ReaderParameters& params) const = 0;
627 
632  virtual void readBinaryGridFromZipped(escript::Data& out,
633  std::string filename, const ReaderParameters& params) const = 0;
634 
639  virtual void writeBinaryGrid(const escript::Data& in, std::string filename,
640  int byteOrder, int dataType) const = 0;
641 
646  virtual bool ownSample(int fsType, index_t id) const = 0;
647 
652  virtual dim_t getNumDataPointsGlobal() const = 0;
653 
658  virtual const dim_t* getNumNodesPerDim() const = 0;
659 
664  virtual const dim_t* getNumElementsPerDim() const = 0;
665 
671  virtual const dim_t* getNumFacesPerBoundary() const = 0;
672 
677  virtual IndexVector getNodeDistribution() const = 0;
678 
683  virtual const int* getNumSubdivisionsPerDim() const = 0;
684 
689  virtual double getLocalCoordinate(dim_t index, int dim) const = 0;
690 
695  virtual boost::python::tuple getGridParameters() const = 0;
696 
702  virtual bool supportsFilter(const boost::python::tuple& t) const;
703 
707  virtual Assembler_ptr createAssembler(const std::string type,
708  const DataMap& options) const {
709  throw SpeckleyException("Domain does not support custom assemblers");
710  }
711 
712  Assembler_ptr createAssemblerFromPython(const std::string type,
713  const boost::python::list& options) const;
714 
719  virtual const double *getLength() const = 0;
720 
725  inline int getOrder() const { return m_order;}
726 
727 protected:
728  int m_numDim;
732  mutable std::vector<int> m_nodeTags, m_nodeTagsInUse;
733  mutable std::vector<int> m_elementTags, m_elementTagsInUse;
734  std::vector<DiracPoint> m_diracPoints;
735  IndexVector m_diracPointNodeIDs; //for borrowSampleID
738  int m_order;
739 
741  void copyData(escript::Data& out, const escript::Data& in) const;
742 
743  // this is const because setTags is const
744  void updateTagsInUse(int fsType) const;
745 
746  void addToSystemMatrix(escript::AbstractSystemMatrix* mat,
747  const IndexVector& nodes, dim_t numEq,
748  const DoubleVector& array) const;
749 
750  void addPoints(const std::vector<double>& coords,
751  const std::vector<int>& tags);
752 
754  void multiplyData(escript::Data& out, const escript::Data& in) const;
755 
756  /***********************************************************************/
757 
759  virtual dim_t getNumNodes() const = 0;
760 
762  virtual dim_t getNumElements() const = 0;
763 
765  virtual dim_t getNumDOF() const = 0;
766 
768  virtual void assembleCoordinates(escript::Data& arg) const = 0;
769 
771  virtual void assembleGradient(escript::Data& out,
772  const escript::Data& in) const = 0;
773 
775  virtual void assembleIntegrate(DoubleVector& integrals,
776  const escript::Data& arg) const = 0;
777 
779  virtual void interpolateNodesOnElements(escript::Data& out,
780  const escript::Data& in,
781  bool reduced) const = 0;
782 
784  virtual void interpolateElementsOnNodes(escript::Data& out,
785  const escript::Data& in) const = 0;
786 
787  virtual dim_t getDofOfNode(dim_t node) const = 0;
788 
790  virtual void reduceElements(escript::Data& out, const escript::Data& in) const = 0;
791 
792 #ifdef ESYS_MPI
793  virtual void balanceNeighbours(escript::Data& data, bool average) const = 0;
795 #endif
796 
797 private:
799  void assemblePDE(escript::AbstractSystemMatrix* mat, escript::Data& rhs,
800  const DataMap& coefs, Assembler_ptr assembler) const;
801 
804  void assemblePDEBoundary(escript::AbstractSystemMatrix* mat,
805  escript::Data& rhs, const DataMap& coefs,
806  Assembler_ptr assembler) const;
807 
808  void assemblePDEDirac(escript::AbstractSystemMatrix* mat,
809  escript::Data& rhs, const DataMap& coefs,
810  Assembler_ptr assembler) const;
811 
813  virtual dim_t findNode(const double *coords) const = 0;
814 };
815 
816 } // end of namespace speckley
817 
818 #endif // __Speckley_DOMAIN_H__
819 
MPI_Comm getMPIComm() const
returns the MPI communicator
Definition: speckley/src/SpeckleyDomain.h:137
AbstractContinuousDomain, base class for continuous domains.
Definition: AbstractContinuousDomain.h:45
std::vector< int > m_elementTagsInUse
Definition: speckley/src/SpeckleyDomain.h:733
virtual int getFunctionOnBoundaryCode() const
returns a function on boundary FunctionSpace code
Definition: speckley/src/SpeckleyDomain.h:398
Definition: FunctionSpace.h:34
Definition: speckley/src/SpeckleyDomain.h:36
Definition: AbstractAssembler.cpp:18
Definition: Speckley.h:52
virtual bool operator!=(const escript::AbstractDomain &other) const
inequality operator
Definition: speckley/src/SpeckleyDomain.h:166
IndexVector m_diracPointNodeIDs
Definition: speckley/src/SpeckleyDomain.h:735
SpeckleyDomain extends the AbstractContinuousDomain interface for the Speckley library and is the bas...
Definition: speckley/src/SpeckleyDomain.h:84
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: speckleycpp.cpp:85
bool probeInterpolationAcross(int fsType_source, const escript::AbstractDomain &domain, int fsType_target, int dim)
Definition: CrossDomainCoupler.cpp:31
virtual escript::JMPI getMPI() const
returns a reference to the MPI information wrapper for this domain
Definition: speckley/src/SpeckleyDomain.h:103
std::vector< real_t > DoubleVector
Definition: Speckley.h:43
virtual int getMPIRank() const
returns the MPI rank of this processor
Definition: speckley/src/SpeckleyDomain.h:115
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:157
virtual int getFunctionOnContactOneCode() const
returns a FunctionOnContactOne code
Definition: speckley/src/SpeckleyDomain.h:431
StatusType m_status
Definition: speckley/src/SpeckleyDomain.h:729
virtual bool supportsContactElements() const
returns true if this domain supports contact elements, false otherwise
Definition: speckley/src/SpeckleyDomain.h:364
Definition: Speckley.h:47
virtual int getReducedFunctionOnBoundaryCode() const
returns a function on boundary with reduced integration order FunctionSpace code
Definition: speckley/src/SpeckleyDomain.h:407
virtual int getDiracDeltaFunctionsCode() const
returns a DiracDeltaFunctions FunctionSpace code
Definition: speckley/src/SpeckleyDomain.h:463
int byteOrder
byte order in the file (used by binary reader only)
Definition: speckley/src/SpeckleyDomain.h:63
boost::shared_ptr< JMPI_ > JMPI
Definition: EsysMPI.h:71
std::map< std::string, escript::Data > DataMap
Definition: speckley/src/domainhelpers.h:24
virtual int getDim() const
returns the number of spatial dimensions of the domain
Definition: speckley/src/SpeckleyDomain.h:156
Definition: Speckley.h:55
TagMap m_tagMap
Definition: speckley/src/SpeckleyDomain.h:731
assembler_t assembler_type
Definition: speckley/src/SpeckleyDomain.h:736
virtual int getReducedSolutionCode() const
returns a ReducedSolution FunctionSpace code
Definition: speckley/src/SpeckleyDomain.h:455
std::vector< dim_t > first
the (global) offset into the data object to start writing into
Definition: speckley/src/SpeckleyDomain.h:54
#define Speckley_DLL_API
Definition: speckley/src/system_dep.h:22
boost::shared_ptr< AbstractTransportProblem > ATP_ptr
Definition: AbstractTransportProblem.h:160
Structure that wraps parameters for the grid reading routines.
Definition: speckley/src/SpeckleyDomain.h:51
std::vector< int > m_nodeTagsInUse
Definition: speckley/src/SpeckleyDomain.h:732
virtual int getReducedContinuousFunctionCode() const
returns a continuous on reduced order nodes FunctionSpace code
Definition: speckley/src/SpeckleyDomain.h:376
virtual int getFunctionCode() const
returns a function FunctionSpace code
Definition: speckley/src/SpeckleyDomain.h:384
virtual int getApproximationOrder(int fsType) const
returns the approximation order used for a function space
Definition: speckley/src/SpeckleyDomain.h:358
virtual int getMPISize() const
returns the number of processors used for this domain
Definition: speckley/src/SpeckleyDomain.h:109
int m_order
element order (will be m_order + 1 quad points in each axis)
Definition: speckley/src/SpeckleyDomain.h:738
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:192
int index_t
type for array/matrix indices used both globally and on each rank
Definition: DataTypes.h:59
virtual void MPIBarrier() const
if compiled for MPI then executes an MPI_Barrier, else does nothing
Definition: speckley/src/SpeckleyDomain.h:121
std::vector< int > reverse
if non-zero, values are written from last index to first index
Definition: speckley/src/SpeckleyDomain.h:61
A struct to contain a dirac point&#39;s information.
Definition: speckley/src/SpeckleyDomain.h:72
boost::shared_ptr< SubWorld > SubWorld_ptr
Definition: SubWorld.h:146
virtual bool onMasterProcessor() const
returns true if on MPI processor 0, else false
Definition: speckley/src/SpeckleyDomain.h:131
Definition: Speckley.h:49
int getOrder() const
returns the order of the domain
Definition: speckley/src/SpeckleyDomain.h:725
virtual int getReducedFunctionCode() const
returns a function with reduced integration order FunctionSpace code
Definition: speckley/src/SpeckleyDomain.h:390
int StatusType
Definition: AbstractDomain.h:48
virtual bool isValidTagName(const std::string &name) const
returns true if name is a defined tag name
Definition: speckley/src/SpeckleyDomain.h:214
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: speckleycpp.cpp:114
Data represents a collection of datapoints.
Definition: Data.h:63
virtual int getTag(const std::string &name) const
returns the tag key for tag name
Definition: speckley/src/SpeckleyDomain.h:201
std::vector< index_t > IndexVector
Definition: Speckley.h:42
SpeckleyException exception class.
Definition: SpeckleyException.h:29
virtual int getReducedFunctionOnContactOneCode() const
returns a FunctionOnContactOne code with reduced integration order
Definition: speckley/src/SpeckleyDomain.h:439
int dataType
data type in the file (used by binary reader only)
Definition: speckley/src/SpeckleyDomain.h:65
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:334
assembler_t
Definition: speckley/src/SpeckleyDomain.h:35
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: speckleycpp.cpp:60
escript::JMPI m_mpiInfo
Definition: speckley/src/SpeckleyDomain.h:730
int MPI_Comm
Definition: EsysMPI.h:41
std::map< std::string, int > simap_t
Definition: speckley/src/SpeckleyDomain.h:44
std::map< std::string, int > TagMap
Definition: Speckley.h:44
dim_t node
Definition: speckley/src/SpeckleyDomain.h:74
Give a short description of what AbstractTransportProblem does.
Definition: AbstractTransportProblem.h:43
virtual int getFunctionOnContactZeroCode() const
return a FunctionOnContactZero code
Definition: speckley/src/SpeckleyDomain.h:415
Base class for escript system matrices.
Definition: AbstractSystemMatrix.h:42
int m_numDim
Definition: speckley/src/SpeckleyDomain.h:728
std::vector< int > multiplier
Definition: speckley/src/SpeckleyDomain.h:59
virtual int getReducedFunctionOnContactZeroCode() const
returns a FunctionOnContactZero code with reduced integration order
Definition: speckley/src/SpeckleyDomain.h:423
boost::shared_ptr< AbstractSystemMatrix > ASM_ptr
Definition: AbstractSystemMatrix.h:32
std::vector< dim_t > numValues
the number of values to read from file
Definition: speckley/src/SpeckleyDomain.h:56
virtual int getContinuousFunctionCode() const
returns a continuous FunctionSpace code
Definition: speckley/src/SpeckleyDomain.h:370
Definition: Speckley.h:51
Base class for all escript domains.
Definition: AbstractDomain.h:45
int tag
Definition: speckley/src/SpeckleyDomain.h:75
virtual int getSolutionCode() const
returns a Solution FunctionSpace code
Definition: speckley/src/SpeckleyDomain.h:447
std::vector< DiracPoint > m_diracPoints
Definition: speckley/src/SpeckleyDomain.h:734
index_t dim_t
Definition: DataTypes.h:64
virtual Assembler_ptr createAssembler(const std::string type, const DataMap &options) const
Definition: speckley/src/SpeckleyDomain.h:707