escript  Revision_
speckley/src/SpeckleyDomain.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 __Speckley_DOMAIN_H__
19 #define __Speckley_DOMAIN_H__
20 
21 #include <speckley/Speckley.h>
22 #include <speckley/SpeckleyException.h>
23 #include <speckley/AbstractAssembler.h>
24 #include <speckley/domainhelpers.h>
25 
26 #include <escript/AbstractContinuousDomain.h>
27 #include <escript/Data.h>
28 #include <escript/FunctionSpace.h>
29 #include <escript/SubWorld.h>
30 
31 #include <boost/python/tuple.hpp>
32 #include <boost/python/list.hpp>
33 
34 namespace speckley {
35 
38 };
39 
40 /* There is no particular significance to this type,
41 It is here as a typedef because a bug in clang++ prevents
42 that compiler from recognising it as a valid part of
43 a constant expression.
44 */
45 typedef std::map<std::string, int> simap_t;
46 
47 
53 {
55  std::vector<dim_t> first;
57  std::vector<dim_t> numValues;
60  std::vector<int> multiplier;
62  std::vector<int> reverse;
64  int byteOrder;
66  int dataType;
67 };
68 
73 struct DiracPoint
74 {
76  int tag;
77 };
78 
86 {
87 public:
93 
98  ~SpeckleyDomain();
99 
104  virtual escript::JMPI getMPI() const { return m_mpiInfo; }
105 
110  virtual int getMPISize() const { return m_mpiInfo->size; }
111 
116  virtual int getMPIRank() const { return m_mpiInfo->rank; }
117 
122  virtual void MPIBarrier() const {
123 #ifdef ESYS_MPI
124  MPI_Barrier(m_mpiInfo->comm);
125 #endif
126  }
127 
132  virtual bool onMasterProcessor() const { return getMPIRank()==0; }
133 
138  MPI_Comm getMPIComm() const { return m_mpiInfo->comm; }
139 
145  virtual bool isValidFunctionSpaceType(int fsType) const;
146 
151  virtual std::string functionSpaceTypeAsString(int fsType) const;
152 
157  virtual int getDim() const { return m_numDim; }
158 
162  virtual bool operator==(const escript::AbstractDomain& other) const;
163 
167  virtual bool operator!=(const escript::AbstractDomain& other) const {
168  return !(operator==(other));
169  }
170 
177  virtual std::pair<int,dim_t> getDataShape(int fsType) const;
178 
185  int getTagFromSampleNo(int fsType, dim_t sampleNo) const;
186 
193  virtual void setTagMap(const std::string& name, int tag) {
194  m_tagMap[name] = tag;
195  }
196 
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;
205  } else {
206  throw SpeckleyException("getTag: invalid tag name");
207  }
208  }
209 
215  virtual bool isValidTagName(const std::string& name) const {
216  return (m_tagMap.find(name)!=m_tagMap.end());
217  }
218 
223  virtual std::string showTagNames() const;
224 
230  virtual void setNewX(const escript::Data& arg);
231 
237  virtual void interpolateOnDomain(escript::Data& target,
238  const escript::Data& source) const;
239 
245  virtual bool probeInterpolationOnDomain(int fsType_source,
246  int fsType_target) const;
247 
255  virtual signed char preferredInterpolationOnDomain(int fsType_source,
256  int fsType_target) const;
257 
264  bool
265  commonFunctionSpace(const std::vector<int>& fs, int& resultcode) const;
266 
272  virtual void interpolateAcross(escript::Data& target,
273  const escript::Data& source) const = 0;
274 
280  int) const = 0;
281 
286  virtual escript::Data getX() const;
287 
288 #ifdef ESYS_HAVE_BOOST_NUMPY
293  virtual boost::python::numpy::ndarray getNumpyX() const;
294 
298  virtual boost::python::numpy::ndarray getConnectivityInfo() const;
299 #endif
300 
305  virtual escript::Data getNormal() const;
306 
310  virtual escript::Data getSize() const;
311 
317  virtual void setToX(escript::Data& arg) const;
318 
325  virtual void setToGradient(escript::Data& out,
326  const escript::Data& in) const;
327 
333  virtual void setTags(int fsType, int newTag, const escript::Data& mask) const;
334 
340  virtual bool isCellOriented(int fsType) const;
341 
348  virtual StatusType getStatus() const { return m_status; }
349 
354  virtual int getNumberOfTagsInUse(int fsType) const;
355 
360  virtual const int* borrowListOfTagsInUse(int fsType) const;
361 
366  virtual bool canTag(int fsType) const;
367 
372  virtual int getApproximationOrder(int fsType) const { return 1; }
373 
378  virtual bool supportsContactElements() const { return false; }
379 
384  virtual int getContinuousFunctionCode() const { return Nodes; }
385 
390  virtual int getReducedContinuousFunctionCode() const {
391  throw SpeckleyException("Speckley does not support reduced functionspaces");
392  }
393 
398  virtual int getFunctionCode() const { return Elements; }
399 
404  virtual int getReducedFunctionCode() const {
405  return ReducedElements;
406  }
407 
412  virtual int getFunctionOnBoundaryCode() const {
413  throw SpeckleyException("Speckley does not support face elements");
414  }
415 
421  virtual int getReducedFunctionOnBoundaryCode() const {
422  throw SpeckleyException("Speckley does not support face elements");
423  }
424 
429  virtual int getFunctionOnContactZeroCode() const {
430  throw SpeckleyException("Speckley does not support contact elements");
431  }
432 
438  throw SpeckleyException("Speckley does not support contact elements");
439  }
440 
445  virtual int getFunctionOnContactOneCode() const {
446  throw SpeckleyException("Speckley does not support contact elements");
447  }
448 
454  throw SpeckleyException("Speckley does not support contact elements");
455  }
456 
461  virtual int getSolutionCode() const {
462  return DegreesOfFreedom;
463  }
464 
469  virtual int getReducedSolutionCode() const {
470  throw SpeckleyException("Speckley does not support reduced function spaces");
471  }
472 
477  virtual int getDiracDeltaFunctionsCode() const { return Points; }
478 
487  virtual int getSystemMatrixTypeId(const boost::python::object& options) const;
488 
498  virtual int getTransportTypeId(int solver, int preconditioner, int package,
499  bool symmetry) const;
500 
506  virtual void setToIntegrals(std::vector<real_t>& integrals,
507  const escript::Data& arg) const;
508  virtual void setToIntegrals(std::vector<cplx_t>& integrals,
509  const escript::Data& arg) const;
510 
511 
517  virtual void addToSystem(escript::AbstractSystemMatrix& mat,
518  escript::Data& rhs, const DataMap& data,
519  Assembler_ptr assembler) const;
520 
525  virtual void addToSystemFromPython(escript::AbstractSystemMatrix& mat,
526  escript::Data& rhs, const boost::python::list& data,
527  Assembler_ptr assembler) const;
528 
534  virtual void addToRHS(escript::Data& rhs, const DataMap& data,
535  Assembler_ptr assembler) const;
536 
541  virtual void addToRHSFromPython(escript::Data& rhs,
542  const boost::python::list& data,
543  Assembler_ptr assembler) const;
544 
550  virtual void addPDEToTransportProblem(escript::AbstractTransportProblem& tp,
551  escript::Data& source, const DataMap& data,
552  Assembler_ptr assembler) const;
557  void addPDEToTransportProblemFromPython(
559  escript::Data& source, const boost::python::list& data,
560  Assembler_ptr assembler) const;
561 
566  virtual escript::ASM_ptr newSystemMatrix(int row_blocksize,
567  const escript::FunctionSpace& row_functionspace,
568  int column_blocksize,
569  const escript::FunctionSpace& column_functionspace, int type) const;
570 
575  virtual escript::ATP_ptr newTransportProblem(int blocksize,
576  const escript::FunctionSpace& functionspace,
577  int type) const;
578 
584  virtual void Print_Mesh_Info(bool full=false) const;
585 
586 
587  /************************************************************************/
588 
594  virtual void write(const std::string& filename) const = 0;
595 
600  virtual std::string getDescription() const = 0;
601 
607  void dump(const std::string& filename) const = 0;
608 
614  const index_t* borrowSampleReferenceIDs(int fsType) const = 0;
615 
622  virtual void setToNormal(escript::Data& out) const = 0;
623 
629  virtual void setToSize(escript::Data& out) const = 0;
630 
635  virtual void readNcGrid(escript::Data& out, std::string filename,
636  std::string varname, const ReaderParameters& params) const = 0;
637 
642  virtual void readBinaryGrid(escript::Data& out, std::string filename,
643  const ReaderParameters& params) const = 0;
644 
650  std::string filename, const ReaderParameters& params) const = 0;
651 
656  virtual void writeBinaryGrid(const escript::Data& in, std::string filename,
657  int byteOrder, int dataType) const = 0;
658 
663  virtual bool ownSample(int fsType, index_t id) const = 0;
664 
669  virtual dim_t getNumDataPointsGlobal() const = 0;
670 
675  virtual const dim_t* getNumNodesPerDim() const = 0;
676 
681  virtual const dim_t* getNumElementsPerDim() const = 0;
682 
688  virtual const dim_t* getNumFacesPerBoundary() const = 0;
689 
694  virtual IndexVector getNodeDistribution() const = 0;
695 
700  virtual const int* getNumSubdivisionsPerDim() const = 0;
701 
706  virtual double getLocalCoordinate(dim_t index, int dim) const = 0;
707 
712  virtual boost::python::tuple getGridParameters() const = 0;
713 
719  virtual bool supportsFilter(const boost::python::tuple& t) const;
720 
724  virtual Assembler_ptr createAssembler(const std::string type,
725  const DataMap& options) const {
726  throw SpeckleyException("Domain does not support custom assemblers");
727  }
728 
729  Assembler_ptr createAssemblerFromPython(const std::string type,
730  const boost::python::list& options) const;
731 
736  virtual const double *getLength() const = 0;
737 
742  inline int getOrder() const { return m_order;}
743 
744 protected:
745  int m_numDim;
749  mutable std::vector<int> m_nodeTags, m_nodeTagsInUse;
750  mutable std::vector<int> m_elementTags, m_elementTagsInUse;
751  std::vector<DiracPoint> m_diracPoints;
752  IndexVector m_diracPointNodeIDs; //for borrowSampleID
755  int m_order;
756 
758  template<typename Scalar>
759  void copyData(escript::Data& out, const escript::Data& in) const;
760 
761  // this is const because setTags is const
762  void updateTagsInUse(int fsType) const;
763 
764  void addToSystemMatrix(escript::AbstractSystemMatrix* mat,
765  const IndexVector& nodes, dim_t numEq,
766  const DoubleVector& array) const;
767 
768  void addPoints(const std::vector<double>& coords,
769  const std::vector<int>& tags);
770 
772  template<typename Scalar>
773  void multiplyData(escript::Data& out, const escript::Data& in) const;
774 
775  /***********************************************************************/
776 
778  virtual dim_t getNumNodes() const = 0;
779 
781  virtual dim_t getNumElements() const = 0;
782 
784  virtual dim_t getNumDOF() const = 0;
785 
787  virtual void assembleCoordinates(escript::Data& arg) const = 0;
788 
790  virtual void assembleGradient(escript::Data& out,
791  const escript::Data& in) const = 0;
792 
794  virtual void assembleIntegrate(std::vector<real_t>& integrals,
795  const escript::Data& arg) const = 0;
796  virtual void assembleIntegrate(std::vector<cplx_t>& integrals,
797  const escript::Data& arg) const = 0;
798 
801  const escript::Data& in,
802  bool reduced) const = 0;
803 
806  const escript::Data& in) const = 0;
807 
808  virtual dim_t getDofOfNode(dim_t node) const = 0;
809 
811  virtual void reduceElements(escript::Data& out, const escript::Data& in) const = 0;
812 
813 #ifdef ESYS_MPI
815  virtual void balanceNeighbours(escript::Data& data, bool average) const = 0;
816 #endif
817 
818 private:
820  void assemblePDE(escript::AbstractSystemMatrix* mat, escript::Data& rhs,
821  const DataMap& coefs, Assembler_ptr assembler) const;
822 
825  void assemblePDEBoundary(escript::AbstractSystemMatrix* mat,
826  escript::Data& rhs, const DataMap& coefs,
827  Assembler_ptr assembler) const;
828 
829  void assemblePDEDiracWrap(escript::AbstractSystemMatrix* mat,
830  escript::Data& rhs, const DataMap& coefs,
831  Assembler_ptr assembler) const;
832 
833  void assemblePDEDirac(escript::AbstractSystemMatrix* mat,
834  escript::Data& rhs, const DataMap& coefs,
835  Assembler_ptr assembler) const;
836  void assembleComplexPDEDirac(escript::AbstractSystemMatrix* mat,
837  escript::Data& rhs, const DataMap& coefs,
838  Assembler_ptr assembler) const;
839 
840  template<typename Scalar>
841  void setToIntegralsWorker(std::vector<Scalar>& integrals,
842  const escript::Data& arg) const;
843 
845  virtual dim_t findNode(const double *coords) const = 0;
846 };
847 
848 } // end of namespace speckley
849 
850 #endif // __Speckley_DOMAIN_H__
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 &params) 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 &params) 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 &params) 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