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 <boost/python/tuple.hpp>
21 #include <boost/python/list.hpp>
22 
23 #include <speckley/Speckley.h>
24 #include <speckley/SpeckleyException.h>
25 #include <speckley/AbstractAssembler.h>
26 #include <speckley/domainhelpers.h>
27 
28 #include <escript/AbstractContinuousDomain.h>
29 #include <escript/Data.h>
30 #include <escript/FunctionSpace.h>
31 #include <escript/SubWorld.h>
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 int getMPISize() const { return m_mpiInfo->size; }
104 
109  virtual int getMPIRank() const { return m_mpiInfo->rank; }
110 
115  virtual void MPIBarrier() const {
116 #ifdef ESYS_MPI
117  MPI_Barrier(m_mpiInfo->comm);
118 #endif
119  }
120 
125  virtual bool onMasterProcessor() const { return getMPIRank()==0; }
126 
131  MPI_Comm getMPIComm() const { return m_mpiInfo->comm; }
132 
138  virtual bool isValidFunctionSpaceType(int fsType) const;
139 
144  virtual std::string functionSpaceTypeAsString(int fsType) const;
145 
150  virtual int getDim() const { return m_numDim; }
151 
155  virtual bool operator==(const escript::AbstractDomain& other) const;
156 
160  virtual bool operator!=(const escript::AbstractDomain& other) const {
161  return !(operator==(other));
162  }
163 
170  virtual std::pair<int,dim_t> getDataShape(int fsType) const;
171 
178  int getTagFromSampleNo(int fsType, dim_t sampleNo) const;
179 
186  virtual void setTagMap(const std::string& name, int tag) {
187  m_tagMap[name] = tag;
188  }
189 
195  virtual int getTag(const std::string& name) const {
196  if (m_tagMap.find(name) != m_tagMap.end()) {
197  return m_tagMap.find(name)->second;
198  } else {
199  throw SpeckleyException("getTag: invalid tag name");
200  }
201  }
202 
208  virtual bool isValidTagName(const std::string& name) const {
209  return (m_tagMap.find(name)!=m_tagMap.end());
210  }
211 
216  virtual std::string showTagNames() const;
217 
223  virtual void setNewX(const escript::Data& arg);
224 
230  virtual void interpolateOnDomain(escript::Data& target,
231  const escript::Data& source) const;
232 
238  virtual bool probeInterpolationOnDomain(int fsType_source,
239  int fsType_target) const;
240 
248  virtual signed char preferredInterpolationOnDomain(int fsType_source,
249  int fsType_target) const;
250 
257  bool
258  commonFunctionSpace(const std::vector<int>& fs, int& resultcode) const;
259 
265  virtual void interpolateAcross(escript::Data& target,
266  const escript::Data& source) const = 0;
267 
272  virtual bool probeInterpolationAcross(int, const escript::AbstractDomain&,
273  int) const = 0;
274 
279  virtual escript::Data getX() const;
280 
285  virtual escript::Data getNormal() const;
286 
290  virtual escript::Data getSize() const;
291 
297  virtual void setToX(escript::Data& arg) const;
298 
305  virtual void setToGradient(escript::Data& out,
306  const escript::Data& in) const;
307 
313  virtual void setTags(int fsType, int newTag, const escript::Data& mask) const;
314 
320  virtual bool isCellOriented(int fsType) const;
321 
328  virtual StatusType getStatus() const { return m_status; }
329 
334  virtual int getNumberOfTagsInUse(int fsType) const;
335 
340  virtual const int* borrowListOfTagsInUse(int fsType) const;
341 
346  virtual bool canTag(int fsType) const;
347 
352  virtual int getApproximationOrder(int fsType) const { return 1; }
353 
358  virtual bool supportsContactElements() const { return false; }
359 
364  virtual int getContinuousFunctionCode() const { return Nodes; }
365 
370  virtual int getReducedContinuousFunctionCode() const {
371  throw SpeckleyException("Speckley does not support reduced functionspaces");
372  }
373 
378  virtual int getFunctionCode() const { return Elements; }
379 
384  virtual int getReducedFunctionCode() const {
385  return ReducedElements;
386  }
387 
392  virtual int getFunctionOnBoundaryCode() const {
393  throw SpeckleyException("Speckley does not support face elements");
394  }
395 
401  virtual int getReducedFunctionOnBoundaryCode() const {
402  throw SpeckleyException("Speckley does not support face elements");
403  }
404 
409  virtual int getFunctionOnContactZeroCode() const {
410  throw SpeckleyException("Speckley does not support contact elements");
411  }
412 
418  throw SpeckleyException("Speckley does not support contact elements");
419  }
420 
425  virtual int getFunctionOnContactOneCode() const {
426  throw SpeckleyException("Speckley does not support contact elements");
427  }
428 
434  throw SpeckleyException("Speckley does not support contact elements");
435  }
436 
441  virtual int getSolutionCode() const {
442  return DegreesOfFreedom;
443  }
444 
449  virtual int getReducedSolutionCode() const {
450  throw SpeckleyException("Speckley does not support reduced function spaces");
451  }
452 
457  virtual int getDiracDeltaFunctionsCode() const { return Points; }
458 
467  virtual int getSystemMatrixTypeId(const boost::python::object& options) const;
468 
478  virtual int getTransportTypeId(int solver, int preconditioner, int package,
479  bool symmetry) const;
480 
486  virtual void setToIntegrals(DoubleVector& integrals, const escript::Data& arg) const;
487 
488 
494  virtual void addToSystem(escript::AbstractSystemMatrix& mat,
495  escript::Data& rhs, const DataMap& data,
496  Assembler_ptr assembler) const;
497 
502  virtual void addToSystemFromPython(escript::AbstractSystemMatrix& mat,
503  escript::Data& rhs, const boost::python::list& data,
504  Assembler_ptr assembler) const;
505 
511  virtual void addToRHS(escript::Data& rhs, const DataMap& data,
512  Assembler_ptr assembler) const;
513 
518  virtual void addToRHSFromPython(escript::Data& rhs,
519  const boost::python::list& data,
520  Assembler_ptr assembler) const;
521 
527  virtual void addPDEToTransportProblem(escript::AbstractTransportProblem& tp,
528  escript::Data& source, const DataMap& data,
529  Assembler_ptr assembler) const;
534  void addPDEToTransportProblemFromPython(
536  escript::Data& source, const boost::python::list& data,
537  Assembler_ptr assembler) const;
538 
543  virtual escript::ASM_ptr newSystemMatrix(int row_blocksize,
544  const escript::FunctionSpace& row_functionspace,
545  int column_blocksize,
546  const escript::FunctionSpace& column_functionspace, int type) const;
547 
552  virtual escript::ATP_ptr newTransportProblem(int blocksize,
553  const escript::FunctionSpace& functionspace,
554  int type) const;
555 
561  virtual void Print_Mesh_Info(bool full=false) const;
562 
563 
564  /************************************************************************/
565 
571  virtual void write(const std::string& filename) const = 0;
572 
577  virtual std::string getDescription() const = 0;
578 
584  void dump(const std::string& filename) const = 0;
585 
591  const index_t* borrowSampleReferenceIDs(int fsType) const = 0;
592 
599  virtual void setToNormal(escript::Data& out) const = 0;
600 
606  virtual void setToSize(escript::Data& out) const = 0;
607 
612  virtual void readNcGrid(escript::Data& out, std::string filename,
613  std::string varname, const ReaderParameters& params) const = 0;
614 
619  virtual void readBinaryGrid(escript::Data& out, std::string filename,
620  const ReaderParameters& params) const = 0;
621 
622 #ifdef USE_BOOSTIO
623 
627  virtual void readBinaryGridFromZipped(escript::Data& out,
628  std::string filename, const ReaderParameters& params) const = 0;
629 #endif
630 
635  virtual void writeBinaryGrid(const escript::Data& in, std::string filename,
636  int byteOrder, int dataType) const = 0;
637 
642  virtual bool ownSample(int fsType, index_t id) const = 0;
643 
648  virtual dim_t getNumDataPointsGlobal() const = 0;
649 
654  virtual const dim_t* getNumNodesPerDim() const = 0;
655 
660  virtual const dim_t* getNumElementsPerDim() const = 0;
661 
667  virtual const dim_t* getNumFacesPerBoundary() const = 0;
668 
673  virtual IndexVector getNodeDistribution() const = 0;
674 
679  virtual const int* getNumSubdivisionsPerDim() const = 0;
680 
685  virtual double getLocalCoordinate(dim_t index, int dim) const = 0;
686 
691  virtual boost::python::tuple getGridParameters() const = 0;
692 
698  virtual bool supportsFilter(const boost::python::tuple& t) const;
699 
703  virtual Assembler_ptr createAssembler(const std::string type,
704  const DataMap& options) const {
705  throw SpeckleyException("Domain does not support custom assemblers");
706  }
707 
708  Assembler_ptr createAssemblerFromPython(const std::string type,
709  const boost::python::list& options) const;
710 
715  virtual const double *getLength() const = 0;
716 
721  inline int getOrder() const { return m_order;}
722 
723 protected:
724  int m_numDim;
725  StatusType m_status;
728  mutable std::vector<int> m_nodeTags, m_nodeTagsInUse;
729  mutable std::vector<int> m_elementTags, m_elementTagsInUse;
730  std::vector<DiracPoint> m_diracPoints;
731  IndexVector m_diracPointNodeIDs; //for borrowSampleID
734  int m_order;
735 
737  void copyData(escript::Data& out, const escript::Data& in) const;
738 
739  // this is const because setTags is const
740  void updateTagsInUse(int fsType) const;
741 
742  void addToSystemMatrix(escript::AbstractSystemMatrix* mat,
743  const IndexVector& nodes, dim_t numEq,
744  const DoubleVector& array) const;
745 
746  void addPoints(const std::vector<double>& coords,
747  const std::vector<int>& tags);
748 
750  void multiplyData(escript::Data& out, const escript::Data& in) const;
751 
752  /***********************************************************************/
753 
755  virtual dim_t getNumNodes() const = 0;
756 
758  virtual dim_t getNumElements() const = 0;
759 
761  virtual dim_t getNumDOF() const = 0;
762 
764  virtual void assembleCoordinates(escript::Data& arg) const = 0;
765 
767  virtual void assembleGradient(escript::Data& out,
768  const escript::Data& in) const = 0;
769 
771  virtual void assembleIntegrate(DoubleVector& integrals,
772  const escript::Data& arg) const = 0;
773 
775  virtual void interpolateNodesOnElements(escript::Data& out,
776  const escript::Data& in,
777  bool reduced) const = 0;
778 
780  virtual void interpolateElementsOnNodes(escript::Data& out,
781  const escript::Data& in) const = 0;
782 
783  virtual dim_t getDofOfNode(dim_t node) const = 0;
784 
786  virtual void reduceElements(escript::Data& out, const escript::Data& in) const = 0;
787 
788 #ifdef ESYS_MPI
789  virtual void balanceNeighbours(escript::Data& data, bool average) const = 0;
791 #endif
792 
793 private:
795  void assemblePDE(escript::AbstractSystemMatrix* mat, escript::Data& rhs,
796  const DataMap& coefs, Assembler_ptr assembler) const;
797 
800  void assemblePDEBoundary(escript::AbstractSystemMatrix* mat,
801  escript::Data& rhs, const DataMap& coefs,
802  Assembler_ptr assembler) const;
803 
804  void assemblePDEDirac(escript::AbstractSystemMatrix* mat,
805  escript::Data& rhs, const DataMap& coefs,
806  Assembler_ptr assembler) const;
807 
809  virtual dim_t findNode(const double *coords) const = 0;
810 };
811 
812 } // end of namespace speckley
813 
814 #endif // __Speckley_DOMAIN_H__
815 
AbstractContinuousDomain, base class for continuous domains.
Definition: AbstractContinuousDomain.h:46
std::vector< int > m_elementTagsInUse
Definition: speckley/src/SpeckleyDomain.h:729
Definition: FunctionSpace.h:34
Definition: speckley/src/SpeckleyDomain.h:36
Definition: AbstractAssembler.cpp:22
virtual int getMPIRank() const
returns the MPI rank of this processor
Definition: speckley/src/SpeckleyDomain.h:109
IndexVector m_diracPointNodeIDs
Definition: speckley/src/SpeckleyDomain.h:731
SpeckleyDomain extends the AbstractContinuousDomain interface for the Speckley library and is the bas...
Definition: speckley/src/SpeckleyDomain.h:84
bool probeInterpolationAcross(int fsType_source, const escript::AbstractDomain &domain, int fsType_target, int dim)
Definition: CrossDomainCoupler.cpp:36
virtual int getFunctionOnContactZeroCode() const
return a FunctionOnContactZero code
Definition: speckley/src/SpeckleyDomain.h:409
virtual Assembler_ptr createAssembler(const std::string type, const DataMap &options) const
Definition: speckley/src/SpeckleyDomain.h:703
virtual int getSolutionCode() const
returns a Solution FunctionSpace code
Definition: speckley/src/SpeckleyDomain.h:441
StatusType m_status
Definition: speckley/src/SpeckleyDomain.h:725
int getOrder() const
returns the order of the domain
Definition: speckley/src/SpeckleyDomain.h:721
int getLength(const escript::Data *data)
Definition: DataC.cpp:92
virtual int getReducedContinuousFunctionCode() const
returns a continuous on reduced order nodes FunctionSpace code
Definition: speckley/src/SpeckleyDomain.h:370
Definition: Speckley.h:44
virtual int getFunctionCode() const
returns a function FunctionSpace code
Definition: speckley/src/SpeckleyDomain.h:378
int byteOrder
byte order in the file (used by binary reader only)
Definition: speckley/src/SpeckleyDomain.h:63
std::map< std::string, escript::Data > DataMap
Definition: speckley/src/domainhelpers.h:24
TagMap m_tagMap
Definition: speckley/src/SpeckleyDomain.h:727
std::vector< double > DoubleVector
Definition: Speckley.h:39
assembler_t assembler_type
Definition: speckley/src/SpeckleyDomain.h:732
std::vector< dim_t > first
the (global) offset into the data object to start writing into
Definition: speckley/src/SpeckleyDomain.h:54
virtual bool operator!=(const escript::AbstractDomain &other) const
inequality operator
Definition: speckley/src/SpeckleyDomain.h:160
virtual int getDim() const
returns the number of spatial dimensions of the domain
Definition: speckley/src/SpeckleyDomain.h:150
#define Speckley_DLL_API
Definition: speckley/src/system_dep.h:22
boost::shared_ptr< AbstractTransportProblem > ATP_ptr
Definition: AbstractTransportProblem.h:162
Structure that wraps parameters for the grid reading routines.
Definition: speckley/src/SpeckleyDomain.h:51
int MPI_Comm
Definition: Esys_MPI.h:38
std::vector< int > m_nodeTagsInUse
Definition: speckley/src/SpeckleyDomain.h:728
int m_order
element order (will be m_order + 1 quad points in each axis)
Definition: speckley/src/SpeckleyDomain.h:734
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:186
virtual int getMPISize() const
returns the number of processors used for this domain
Definition: speckley/src/SpeckleyDomain.h:103
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
Definition: Speckley.h:49
boost::shared_ptr< SubWorld > SubWorld_ptr
Definition: SubWorld.h:142
virtual int getFunctionOnBoundaryCode() const
returns a function on boundary FunctionSpace code
Definition: speckley/src/SpeckleyDomain.h:392
Definition: Speckley.h:48
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:116
Data represents a collection of datapoints.
Definition: Data.h:68
virtual bool onMasterProcessor() const
returns true if on MPI processor 0, else false
Definition: speckley/src/SpeckleyDomain.h:125
Definition: Speckley.h:52
std::vector< index_t > IndexVector
Definition: Speckley.h:38
SpeckleyException exception class.
Definition: SpeckleyException.h:29
int dataType
data type in the file (used by binary reader only)
Definition: speckley/src/SpeckleyDomain.h:65
virtual int getReducedFunctionCode() const
returns a function with reduced integration order FunctionSpace code
Definition: speckley/src/SpeckleyDomain.h:384
Definition: Speckley.h:46
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:159
assembler_t
Definition: speckley/src/SpeckleyDomain.h:35
MPI_Comm getMPIComm() const
returns the MPI communicator
Definition: speckley/src/SpeckleyDomain.h:131
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:64
virtual void MPIBarrier() const
if compiled for MPI then executes an MPI_Barrier, else does nothing
Definition: speckley/src/SpeckleyDomain.h:115
virtual int getContinuousFunctionCode() const
returns a continuous FunctionSpace code
Definition: speckley/src/SpeckleyDomain.h:364
std::map< std::string, int > simap_t
Definition: speckley/src/SpeckleyDomain.h:44
int index_t
Definition: types.h:24
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:328
std::map< std::string, int > TagMap
Definition: Speckley.h:41
virtual int getReducedSolutionCode() const
returns a ReducedSolution FunctionSpace code
Definition: speckley/src/SpeckleyDomain.h:449
virtual int getReducedFunctionOnContactZeroCode() const
returns a FunctionOnContactZero code with reduced integration order
Definition: speckley/src/SpeckleyDomain.h:417
dim_t node
Definition: speckley/src/SpeckleyDomain.h:74
Give a short description of what AbstractTransportProblem does.
Definition: AbstractTransportProblem.h:45
virtual bool isValidTagName(const std::string &name) const
returns true if name is a defined tag name
Definition: speckley/src/SpeckleyDomain.h:208
index_t dim_t
Definition: types.h:27
virtual int getFunctionOnContactOneCode() const
returns a FunctionOnContactOne code
Definition: speckley/src/SpeckleyDomain.h:425
Base class for escript system matrices.
Definition: AbstractSystemMatrix.h:37
virtual int getApproximationOrder(int fsType) const
returns the approximation order used for a function space
Definition: speckley/src/SpeckleyDomain.h:352
virtual int getTag(const std::string &name) const
returns the tag key for tag name
Definition: speckley/src/SpeckleyDomain.h:195
int m_numDim
Definition: speckley/src/SpeckleyDomain.h:724
std::vector< int > multiplier
Definition: speckley/src/SpeckleyDomain.h:59
std::vector< dim_t > numValues
the number of values to read from file
Definition: speckley/src/SpeckleyDomain.h:56
boost::shared_ptr< AbstractSystemMatrix > ASM_ptr
Definition: AbstractSystemMatrix.h:170
Base class for all escript domains.
Definition: AbstractDomain.h:45
virtual int getReducedFunctionOnBoundaryCode() const
returns a function on boundary with reduced integration order FunctionSpace code
Definition: speckley/src/SpeckleyDomain.h:401
int tag
Definition: speckley/src/SpeckleyDomain.h:75
boost::shared_ptr< JMPI_ > JMPI
Definition: Esys_MPI.h:79
esysUtils::JMPI m_mpiInfo
Definition: speckley/src/SpeckleyDomain.h:726
virtual bool supportsContactElements() const
returns true if this domain supports contact elements, false otherwise
Definition: speckley/src/SpeckleyDomain.h:358
std::vector< DiracPoint > m_diracPoints
Definition: speckley/src/SpeckleyDomain.h:730
virtual int getDiracDeltaFunctionsCode() const
returns a DiracDeltaFunctions FunctionSpace code
Definition: speckley/src/SpeckleyDomain.h:457
virtual int getReducedFunctionOnContactOneCode() const
returns a FunctionOnContactOne code with reduced integration order
Definition: speckley/src/SpeckleyDomain.h:433