escript  Revision_
speckley/src/SpeckleyDomain.h
Go to the documentation of this file.
1 
2 /*****************************************************************************
3 *
4 * Copyright (c) 2003-2018 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 
35 enum assembler_t {
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 
51 struct ReaderParameters
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 
423  virtual int getReducedFunctionOnContactZeroCode() const {
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 
439  virtual int getReducedFunctionOnContactOneCode() const {
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(std::vector<real_t>& integrals,
493  const escript::Data& arg) const;
494  virtual void setToIntegrals(std::vector<cplx_t>& integrals,
495  const escript::Data& arg) const;
496 
497 
503  virtual void addToSystem(escript::AbstractSystemMatrix& mat,
504  escript::Data& rhs, const DataMap& data,
505  Assembler_ptr assembler) const;
506 
511  virtual void addToSystemFromPython(escript::AbstractSystemMatrix& mat,
512  escript::Data& rhs, const boost::python::list& data,
513  Assembler_ptr assembler) const;
514 
520  virtual void addToRHS(escript::Data& rhs, const DataMap& data,
521  Assembler_ptr assembler) const;
522 
527  virtual void addToRHSFromPython(escript::Data& rhs,
528  const boost::python::list& data,
529  Assembler_ptr assembler) const;
530 
536  virtual void addPDEToTransportProblem(escript::AbstractTransportProblem& tp,
537  escript::Data& source, const DataMap& data,
538  Assembler_ptr assembler) const;
543  void addPDEToTransportProblemFromPython(
545  escript::Data& source, const boost::python::list& data,
546  Assembler_ptr assembler) const;
547 
552  virtual escript::ASM_ptr newSystemMatrix(int row_blocksize,
553  const escript::FunctionSpace& row_functionspace,
554  int column_blocksize,
555  const escript::FunctionSpace& column_functionspace, int type) const;
556 
561  virtual escript::ATP_ptr newTransportProblem(int blocksize,
562  const escript::FunctionSpace& functionspace,
563  int type) const;
564 
570  virtual void Print_Mesh_Info(bool full=false) const;
571 
572 
573  /************************************************************************/
574 
580  virtual void write(const std::string& filename) const = 0;
581 
586  virtual std::string getDescription() const = 0;
587 
593  void dump(const std::string& filename) const = 0;
594 
600  const index_t* borrowSampleReferenceIDs(int fsType) const = 0;
601 
608  virtual void setToNormal(escript::Data& out) const = 0;
609 
615  virtual void setToSize(escript::Data& out) const = 0;
616 
621  virtual void readNcGrid(escript::Data& out, std::string filename,
622  std::string varname, const ReaderParameters& params) const = 0;
623 
628  virtual void readBinaryGrid(escript::Data& out, std::string filename,
629  const ReaderParameters& params) const = 0;
630 
635  virtual void readBinaryGridFromZipped(escript::Data& out,
636  std::string filename, const ReaderParameters& params) const = 0;
637 
642  virtual void writeBinaryGrid(const escript::Data& in, std::string filename,
643  int byteOrder, int dataType) const = 0;
644 
649  virtual bool ownSample(int fsType, index_t id) const = 0;
650 
655  virtual dim_t getNumDataPointsGlobal() const = 0;
656 
661  virtual const dim_t* getNumNodesPerDim() const = 0;
662 
667  virtual const dim_t* getNumElementsPerDim() const = 0;
668 
674  virtual const dim_t* getNumFacesPerBoundary() const = 0;
675 
680  virtual IndexVector getNodeDistribution() const = 0;
681 
686  virtual const int* getNumSubdivisionsPerDim() const = 0;
687 
692  virtual double getLocalCoordinate(dim_t index, int dim) const = 0;
693 
698  virtual boost::python::tuple getGridParameters() const = 0;
699 
705  virtual bool supportsFilter(const boost::python::tuple& t) const;
706 
710  virtual Assembler_ptr createAssembler(const std::string type,
711  const DataMap& options) const {
712  throw SpeckleyException("Domain does not support custom assemblers");
713  }
714 
715  Assembler_ptr createAssemblerFromPython(const std::string type,
716  const boost::python::list& options) const;
717 
722  virtual const double *getLength() const = 0;
723 
728  inline int getOrder() const { return m_order;}
729 
730 protected:
731  int m_numDim;
732  StatusType m_status;
733  escript::JMPI m_mpiInfo;
734  TagMap m_tagMap;
735  mutable std::vector<int> m_nodeTags, m_nodeTagsInUse;
736  mutable std::vector<int> m_elementTags, m_elementTagsInUse;
737  std::vector<DiracPoint> m_diracPoints;
738  IndexVector m_diracPointNodeIDs; //for borrowSampleID
739  assembler_t assembler_type;
741  int m_order;
742 
744  template<typename Scalar>
745  void copyData(escript::Data& out, const escript::Data& in) const;
746 
747  // this is const because setTags is const
748  void updateTagsInUse(int fsType) const;
749 
750  void addToSystemMatrix(escript::AbstractSystemMatrix* mat,
751  const IndexVector& nodes, dim_t numEq,
752  const DoubleVector& array) const;
753 
754  void addPoints(const std::vector<double>& coords,
755  const std::vector<int>& tags);
756 
758  template<typename Scalar>
759  void multiplyData(escript::Data& out, const escript::Data& in) const;
760 
761  /***********************************************************************/
762 
764  virtual dim_t getNumNodes() const = 0;
765 
767  virtual dim_t getNumElements() const = 0;
768 
770  virtual dim_t getNumDOF() const = 0;
771 
773  virtual void assembleCoordinates(escript::Data& arg) const = 0;
774 
776  virtual void assembleGradient(escript::Data& out,
777  const escript::Data& in) const = 0;
778 
780  virtual void assembleIntegrate(std::vector<real_t>& integrals,
781  const escript::Data& arg) const = 0;
782  virtual void assembleIntegrate(std::vector<cplx_t>& integrals,
783  const escript::Data& arg) const = 0;
784 
786  virtual void interpolateNodesOnElements(escript::Data& out,
787  const escript::Data& in,
788  bool reduced) const = 0;
789 
791  virtual void interpolateElementsOnNodes(escript::Data& out,
792  const escript::Data& in) const = 0;
793 
794  virtual dim_t getDofOfNode(dim_t node) const = 0;
795 
797  virtual void reduceElements(escript::Data& out, const escript::Data& in) const = 0;
798 
799 #ifdef ESYS_MPI
800  virtual void balanceNeighbours(escript::Data& data, bool average) const = 0;
802 #endif
803 
804 private:
806  void assemblePDE(escript::AbstractSystemMatrix* mat, escript::Data& rhs,
807  const DataMap& coefs, Assembler_ptr assembler) const;
808 
811  void assemblePDEBoundary(escript::AbstractSystemMatrix* mat,
812  escript::Data& rhs, const DataMap& coefs,
813  Assembler_ptr assembler) const;
814 
815  void assemblePDEDirac(escript::AbstractSystemMatrix* mat,
816  escript::Data& rhs, const DataMap& coefs,
817  Assembler_ptr assembler) const;
818 
819  template<typename Scalar>
820  void setToIntegralsWorker(std::vector<Scalar>& integrals,
821  const escript::Data& arg) const;
822 
824  virtual dim_t findNode(const double *coords) const = 0;
825 };
826 
827 } // end of namespace speckley
828 
829 #endif // __Speckley_DOMAIN_H__
830 
speckley::DATATYPE_INT32
Definition: speckley/src/system_dep.h:59
speckley::simap_t
std::map< std::string, int > simap_t
Definition: speckley/src/SpeckleyDomain.h:55
speckley::TagMap
std::map< std::string, int > TagMap
Definition: Speckley.h:57
weipa::SpeckleyDomain::initFromFile
virtual bool initFromFile(const std::string &filename)
Reads the domain from a dump file.
Definition: weipa/src/SpeckleyDomain.cpp:67
escript::Data::getNumDataPointsPerSample
int getNumDataPointsPerSample() const
Return the number of data points per sample.
Definition: Data.h:530
ripley::isNotEmpty
bool isNotEmpty(const std::string target, const DataMap &mapping)
Definition: ripley/src/domainhelpers.h:50
escript::Data::isComplex
bool isComplex() const
True if components of this data are stored as complex.
Definition: Data.cpp:1163
speckley::readNcGrid
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:113
escript::FunctionSpace::getSize
escript::Data getSize() const
Returns the sample size (e.g. the diameter of elements, radius of particles).
Definition: FunctionSpace.cpp:245
speckley::SpeckleyDomain::showTagNames
virtual std::string showTagNames() const
returns all tag names in a single string separated by commas
Definition: speckley/src/SpeckleyDomain.cpp:146
dudley::Points
Definition: Dudley.h:68
ripley::DataMap
std::map< std::string, escript::Data > DataMap
Definition: ripley/src/domainhelpers.h:35
speckley::DoubleVector
std::vector< real_t > DoubleVector
Definition: Speckley.h:56
REGISTER_ESCRIPT_EXCEPTION_TRANSLATORS
#define REGISTER_ESCRIPT_EXCEPTION_TRANSLATORS
Definition: ExceptionTranslators.h:24
getTag
int getTag(unsigned char sourcex, unsigned char sourcey, unsigned char sourcez, unsigned char targetx, unsigned char targety, unsigned char targetz)
Definition: blocktools.cpp:412
speckley::assembler_t
assembler_t
Definition: speckley/src/SpeckleyDomain.h:46
speckley::SpeckleyDomain::getDescription
virtual std::string getDescription() const =0
returns a description for this domain
escript::Data::requireWrite
void requireWrite()
Ensures data is ready for write access. This means that the data will be resolved if lazy and will be...
Definition: Data.cpp:1242
weipa::SpeckleyDomain
Represents a full Speckley domain including nodes and elements.
Definition: weipa/src/SpeckleyDomain.h:46
weipa::SpeckleyDomain::getCenteringForFunctionSpace
virtual Centering getCenteringForFunctionSpace(int fsCode) const
Returns whether data on given function space is node or cell centered.
Definition: weipa/src/SpeckleyDomain.cpp:73
speckley::SpeckleyDomain::getNormal
virtual escript::Data getNormal() const
returns boundary normals at the quadrature point on the face elements
Definition: speckley/src/SpeckleyDomain.cpp:436
dudley::Nodes
Definition: Dudley.h:63
speckley::DATATYPE_FLOAT64
Definition: speckley/src/system_dep.h:61
INDEX2
#define INDEX2(_X1_, _X2_, _N1_)
Definition: index.h:21
speckley::SpeckleyDomain::addToRHSFromPython
virtual void addToRHSFromPython(escript::Data &rhs, const boost::python::list &data, Assembler_ptr assembler) const
a wrapper for addToRHS that allows calling from Python
Definition: speckley/src/SpeckleyDomain.cpp:746
escript::Data::getSampleDataRO
const DataTypes::real_t * getSampleDataRO(DataTypes::RealVectorType::size_type sampleNo, DataTypes::real_t dummy=0) const
Return the sample data for the given sample no. Please do not use this unless you NEED to access samp...
Definition: Data.h:1975
speckley::_rectangle
escript::Domain_ptr _rectangle(int order, double _n0, double _n1, const object &l0, const object &l1, int d0, int d1, const object &objpoints, const object &objtags, escript::SubWorld_ptr world)
Definition: speckleycpp.cpp:233
weipa::VisItControl::initialized
bool initialized
Definition: VisItControl.cpp:60
escript::continuousFunction
FunctionSpace continuousFunction(const AbstractDomain &domain)
Create function space objects.
Definition: FunctionSpaceFactory.cpp:41
weipa::SpeckleyDomain::siloPath
std::string siloPath
Definition: weipa/src/SpeckleyDomain.h:88
escript::Data::getDataPointShape
const DataTypes::ShapeType & getDataPointShape() const
Return a reference to the data point shape.
Definition: Data.h:691
speckley::SpeckleyDomain::onMasterProcessor
virtual bool onMasterProcessor() const
returns true if on MPI processor 0, else false
Definition: speckley/src/SpeckleyDomain.h:142
speckley::SpeckleyDomain::createAssemblerFromPython
Assembler_ptr createAssemblerFromPython(const std::string type, const boost::python::list &options) const
Definition: speckley/src/SpeckleyDomain.cpp:738
speckley::SpeckleyDomain::dump
void dump(const std::string &filename) const =0
dumps the mesh to a file with the given name
speckley::SpeckleyDomain::getOrder
int getOrder() const
returns the order of the domain
Definition: speckley/src/SpeckleyDomain.h:739
speckley::ReaderParameters::byteOrder
int byteOrder
byte order in the file (used by binary reader only)
Definition: speckley/src/SpeckleyDomain.h:74
speckley::readBinaryGridFromZipped
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:84
speckley::SpeckleyDomain::SpeckleyDomain
SpeckleyDomain(dim_t dim, int order, escript::SubWorld_ptr p=escript::SubWorld_ptr())
Constructor with number of dimensions. Allocates MPI info structure.
Definition: speckley/src/SpeckleyDomain.cpp:59
weipa::SpeckleyDomain::getMeshByName
virtual NodeData_ptr getMeshByName(const std::string &name) const
Returns the node mesh with given name.
Definition: weipa/src/SpeckleyDomain.cpp:198
weipa::SpeckleyNodes_ptr
boost::shared_ptr< SpeckleyNodes > SpeckleyNodes_ptr
Definition: SpeckleyNodes.h:29
escript::Domain_ptr
boost::shared_ptr< AbstractDomain > Domain_ptr
Definition: AbstractDomain.h:47
escript::AbstractTransportProblem
Give a short description of what AbstractTransportProblem does.
Definition: AbstractTransportProblem.h:54
finley::ReducedNodes
Definition: Finley.h:77
dudley::DegreesOfFreedom
Definition: Dudley.h:62
speckley::SpeckleyDomain::Print_Mesh_Info
virtual void Print_Mesh_Info(bool full=false) const
writes information about the mesh to standard output
Definition: speckley/src/SpeckleyDomain.cpp:687
escript::Data::getDataPointSize
int getDataPointSize() const
Return the size of the data point. It is the product of the data point shape dimensions.
Definition: Data.cpp:1363
weipa::SpeckleyDomain::writeToSilo
virtual bool writeToSilo(DBfile *dbfile, const std::string &pathInSilo, const StringVec &labels, const StringVec &units, bool writeMeshData)
Writes the domain to a Silo file.
Definition: weipa/src/SpeckleyDomain.cpp:220
speckley::ReaderParameters::first
std::vector< dim_t > first
the (global) offset into the data object to start writing into
Definition: speckley/src/SpeckleyDomain.h:65
speckley::SpeckleyDomain::newSystemMatrix
virtual escript::ASM_ptr newSystemMatrix(int row_blocksize, const escript::FunctionSpace &row_functionspace, int column_blocksize, const escript::FunctionSpace &column_functionspace, int type) const
creates a stiffness matrix and initializes it with zeros
Definition: speckley/src/SpeckleyDomain.cpp:714
speckley::SpeckleyDomain::getGridParameters
virtual boost::python::tuple getGridParameters() const =0
returns the tuple (origin, spacing, number_of_elements)
ripley::unpackData
const escript::Data unpackData(const std::string target, const DataMap &mapping)
Definition: ripley/src/domainhelpers.h:39
weipa::SpeckleyDomain::getMeshNames
virtual StringVec getMeshNames() const
Returns the names of all meshes within this domain.
Definition: weipa/src/SpeckleyDomain.cpp:124
escript::FunctionSpace
Definition: FunctionSpace.h:45
escript::makeInfo
JMPI makeInfo(MPI_Comm comm, bool owncom)
Definition: EsysMPI.cpp:39
weipa::SpeckleyDomain::faces
SpeckleyElements_ptr faces
Definition: weipa/src/SpeckleyDomain.h:87
escript::Data::getFunctionSpace
const FunctionSpace & getFunctionSpace() const
Return the function space.
Definition: Data.h:461
dudley::Elements
Definition: Dudley.h:64
weipa::DataVar_ptr
boost::shared_ptr< DataVar > DataVar_ptr
Definition: weipa.h:63
weipa
Definition: DataVar.cpp:49
speckley::probeInterpolationAcross
bool probeInterpolationAcross(int fsType_source, const escript::AbstractDomain &domain, int fsType_target, int dim)
Definition: CrossDomainCoupler.cpp:30
escript::FunctionSpace::getDomain
const_Domain_ptr getDomain() const
Returns the function space domain.
Definition: FunctionSpace.cpp:101
MPI_INT
#define MPI_INT
Definition: EsysMPI.h:44
speckley::BYTEORDER_NATIVE
Definition: speckley/src/system_dep.h:46
speckley
Definition: AbstractAssembler.cpp:17
weipa::ZONE_CENTERED
Definition: DomainChunk.h:31
escript::Data::getNumSamples
int getNumSamples() const
Return the number of samples.
Definition: Data.h:519
paso::util::copy
void copy(dim_t N, double *out, const double *in)
out = in
Definition: PasoUtil.h:110
speckley::SpeckleyDomain::setTagMap
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:203
speckley::readBinaryGrid
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:59
weipa::SpeckleyDomain::initFromEscript
virtual bool initFromEscript(const escript::AbstractDomain *domain)
Initialises the domain using an escript domain instance.
Definition: weipa/src/SpeckleyDomain.cpp:48
MPI_COMM_WORLD
#define MPI_COMM_WORLD
Definition: EsysMPI.h:46
escript::Data
Data represents a collection of datapoints.
Definition: Data.h:62
paso::util::l2
double l2(dim_t n, const double *x, escript::JMPI mpiinfo)
returns the global L2 norm of x
Definition: PasoUtil.cpp:524
weipa::SpeckleyDomain::getDataVarByName
virtual DataVar_ptr getDataVarByName(const std::string &name) const
Creates and returns a variable with domain data.
Definition: weipa/src/SpeckleyDomain.cpp:155
escript::FunctionSpace::getX
escript::Data getX() const
Returns the spatial locations of the data points.
Definition: FunctionSpace.cpp:227
escript::DataTypes::dim_t
index_t dim_t
Definition: DataTypes.h:87
speckley::tupleListToMap
void tupleListToMap(DataMap &mapping, const bp::list &list)
Definition: speckley/src/SpeckleyDomain.cpp:44
speckley::ReaderParameters::reverse
std::vector< int > reverse
if non-zero, values are written from last index to first index
Definition: speckley/src/SpeckleyDomain.h:72
escript::JMPI
boost::shared_ptr< JMPI_ > JMPI
Definition: EsysMPI.h:70
speckley::SpeckleyDomain::readBinaryGridFromZipped
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
finley::ReducedDegreesOfFreedom
Definition: Finley.h:75
escript::ATP_ptr
boost::shared_ptr< AbstractTransportProblem > ATP_ptr
Definition: AbstractTransportProblem.h:171
escript::AbstractSystemMatrix
Base class for escript system matrices.
Definition: AbstractSystemMatrix.h:53
weipa::SpeckleyElements
Stores and manipulates one type of speckley mesh elements (cells, faces).
Definition: SpeckleyElements.h:43
speckley::DATATYPE_FLOAT32
Definition: speckley/src/system_dep.h:60
speckley::SpeckleyDomain::getX
virtual escript::Data getX() const
returns locations in the SEM nodes
Definition: speckley/src/SpeckleyDomain.cpp:431
escript::Data::getSampleDataRW
DataTypes::real_t * getSampleDataRW(DataTypes::RealVectorType::size_type sampleNo, DataTypes::real_t dummy=0)
Return the sample data for the given sample no. Please do not use this unless you NEED to access samp...
Definition: Data.h:1940
speckley::SpeckleyDomain::getMPISize
virtual int getMPISize() const
returns the number of processors used for this domain
Definition: speckley/src/SpeckleyDomain.h:120
weipa::ElementData_ptr
boost::shared_ptr< ElementData > ElementData_ptr
Definition: weipa.h:65
speckley::Brick
Brick is the 3-dimensional implementation of a SpeckleyDomain.
Definition: speckley/src/Brick.h:43
weipa::SpeckleyDomain::reorderGhostZones
virtual void reorderGhostZones(int ownIndex)
Reorders elements so that 'ghost' elements (i.e. those that do not belong to ownIndex) appear last.
Definition: weipa/src/SpeckleyDomain.cpp:210
speckley::SpeckleyDomain::writeBinaryGrid
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
weipa::NODE_CENTERED
Definition: DomainChunk.h:30
speckley::SpeckleyDomain
SpeckleyDomain extends the AbstractContinuousDomain interface for the Speckley library and is the bas...
Definition: speckley/src/SpeckleyDomain.h:95
speckley::ReducedElements
Definition: Speckley.h:75
speckley::SpeckleyDomain::MPIBarrier
virtual void MPIBarrier() const
if compiled for MPI then executes an MPI_Barrier, else does nothing
Definition: speckley/src/SpeckleyDomain.h:132
escript::Data::isEmpty
bool isEmpty() const
Definition: Data.cpp:1135
ripley::DoubleVector
std::vector< real_t > DoubleVector
Definition: Ripley.h:56
weipa::SpeckleyDomain::initialized
bool initialized
Definition: weipa/src/SpeckleyDomain.h:84
speckley::FaceElements
Definition: Speckley.h:76
speckley::SpeckleyDomain::isValidTagName
virtual bool isValidTagName(const std::string &name) const
returns true if name is a defined tag name
Definition: speckley/src/SpeckleyDomain.h:225
weipa::SpeckleyElements_ptr
boost::shared_ptr< SpeckleyElements > SpeckleyElements_ptr
Definition: SpeckleyElements.h:31
speckley::Rectangle
Rectangle is the 2-dimensional implementation of a SpeckleyDomain.
Definition: speckley/src/Rectangle.h:43
speckley::SpeckleyDomain::getTag
virtual int getTag(const std::string &name) const
returns the tag key for tag name
Definition: speckley/src/SpeckleyDomain.h:212
speckley::DegreesOfFreedom
Definition: Speckley.h:70
dudley::ReducedElements
Definition: Dudley.h:65
speckley::Elements
Definition: Speckley.h:74
weipa::SpeckleyDomain::getElementsForFunctionSpace
virtual ElementData_ptr getElementsForFunctionSpace(int fsCode) const
Returns the element data for given function space code.
Definition: weipa/src/SpeckleyDomain.cpp:95
weipa::SpeckleyDomain::getVarNames
virtual StringVec getVarNames() const
Returns the names of all 'special' domain variables.
Definition: weipa/src/SpeckleyDomain.cpp:140
speckley::ReducedNodes
Definition: Speckley.h:73
escript::SubWorld_ptr
boost::shared_ptr< SubWorld > SubWorld_ptr
Definition: SubWorld.h:157
speckley::SpeckleyDomain::getDataShape
virtual std::pair< int, dim_t > getDataShape(int fsType) const
returns the number of data points per sample, and the number of samples as a pair.
Definition: speckley/src/SpeckleyDomain.cpp:120
MPI_MIN
#define MPI_MIN
Definition: EsysMPI.h:51
weipa::SpeckleyDomain::cells
SpeckleyElements_ptr cells
Definition: weipa/src/SpeckleyDomain.h:86
speckley::BYTEORDER_LITTLE_ENDIAN
Definition: speckley/src/system_dep.h:54
weipa::SpeckleyDomain::getElementsByName
virtual ElementData_ptr getElementsByName(const std::string &name) const
Returns element data with given name.
Definition: weipa/src/SpeckleyDomain.cpp:187
weipa::SpeckleyDomain::nodes
SpeckleyNodes_ptr nodes
Definition: weipa/src/SpeckleyDomain.h:85
speckley::BYTEORDER_BIG_ENDIAN
Definition: speckley/src/system_dep.h:55
speckley::ReaderParameters::dataType
int dataType
data type in the file (used by binary reader only)
Definition: speckley/src/SpeckleyDomain.h:76
BOOST_PYTHON_MODULE
BOOST_PYTHON_MODULE(speckleycpp)
Definition: speckleycpp.cpp:317
dudley::FaceElements
Definition: Dudley.h:66
escript::RuntimeErrorTranslator
void RuntimeErrorTranslator(const EsysException &e)
Function which translates an EsysException into a python RuntimeError.
Definition: ExceptionTranslators.cpp:47
weipa::SpeckleyDomain::removeGhostZones
virtual void removeGhostZones(int ownIndex)
Removes 'ghost' elements and nodes.
Definition: weipa/src/SpeckleyDomain.cpp:215
escript::DataTypes::index_t
int index_t
type for array/matrix indices used both globally and on each rank
Definition: DataTypes.h:82
escript::Data::isExpanded
bool isExpanded() const
Return true if this Data is expanded.
Definition: Data.cpp:1114
speckley::ReaderParameters::multiplier
std::vector< int > multiplier
Definition: speckley/src/SpeckleyDomain.h:70
dudley::ReducedFaceElements
Definition: Dudley.h:67
speckley::SpeckleyDomain::readBinaryGrid
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
weipa::NodeData_ptr
boost::shared_ptr< NodeData > NodeData_ptr
Definition: weipa.h:67
weipa::StringVec
std::vector< std::string > StringVec
Definition: weipa.h:59
escript::FunctionSpace::getTypeCode
int getTypeCode() const
Returns the function space type code.
Definition: FunctionSpace.cpp:93
speckley::DEFAULT_ASSEMBLER
Definition: speckley/src/SpeckleyDomain.h:59
escript::FunctionSpace::getDim
int getDim() const
Return the number of spatial dimensions of the underlying domain.
Definition: FunctionSpace.h:186
speckley::SpeckleyDomain::getMPIRank
virtual int getMPIRank() const
returns the MPI rank of this processor
Definition: speckley/src/SpeckleyDomain.h:126
escript::function
FunctionSpace function(const AbstractDomain &domain)
Return a function FunctionSpace.
Definition: FunctionSpaceFactory.cpp:53
weipa::SpeckleyDomain::getNodes
virtual NodeData_ptr getNodes() const
Returns a pointer to the full nodes.
Definition: weipa/src/SpeckleyDomain.h:79
speckley::DiracPoint::node
dim_t node
Definition: speckley/src/SpeckleyDomain.h:85
escript::Vector
Data Vector(double value, const FunctionSpace &what, bool expanded)
Return a Data object containing vector data-points. ie: rank 1 data-points.
Definition: DataFactory.cpp:96
speckley::ReaderParameters
Structure that wraps parameters for the grid reading routines.
Definition: speckley/src/SpeckleyDomain.h:62
speckley::SpeckleyDomain::getDim
virtual int getDim() const
returns the number of spatial dimensions of the domain
Definition: speckley/src/SpeckleyDomain.h:167
weipa::SpeckleyDomain::SpeckleyDomain
SpeckleyDomain()
Definition: weipa/src/SpeckleyDomain.cpp:34
speckley::Points
Definition: Speckley.h:78
MPI_Comm
int MPI_Comm
Definition: EsysMPI.h:40
escript::Data::actsExpanded
bool actsExpanded() const
Return true if this Data is expanded or resolves to expanded. That is, if it has a separate value for...
Definition: Data.cpp:1121
escript::AbstractDomain
Base class for all escript domains.
Definition: AbstractDomain.h:56
speckley::IndexVector
std::vector< index_t > IndexVector
Definition: Speckley.h:55
speckley::ReaderParameters::numValues
std::vector< dim_t > numValues
the number of values to read from file
Definition: speckley/src/SpeckleyDomain.h:67
Speckley_DLL_API
#define Speckley_DLL_API
Definition: speckley/src/system_dep.h:21
speckley::SpeckleyDomain::getNumDataPointsGlobal
virtual dim_t getNumDataPointsGlobal() const =0
returns the number of data points summed across all MPI processes
escript::ASM_ptr
boost::shared_ptr< AbstractSystemMatrix > ASM_ptr
Definition: AbstractSystemMatrix.h:43
speckley::DiracPoint
A struct to contain a dirac point's information.
Definition: speckley/src/SpeckleyDomain.h:83
weipa::Centering
Centering
Definition: DomainChunk.h:29
escript::AbstractSystemMatrix::getRowBlockSize
int getRowBlockSize() const
returns the row block size
Definition: AbstractSystemMatrix.h:124
escript::Scalar
Data Scalar(double value, const FunctionSpace &what, bool expanded)
A collection of factory functions for creating Data objects which contain data points of various shap...
Definition: DataFactory.cpp:60
speckley::SpeckleyDomain::getSize
virtual escript::Data getSize() const
returns the element size
Definition: speckley/src/SpeckleyDomain.cpp:441
speckley::SpeckleyDomain::readNcGrid
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
speckley::DataMap
std::map< std::string, escript::Data > DataMap
Definition: speckley/src/domainhelpers.h:35
escript::Data::expand
void expand()
Whatever the current Data type make this into a DataExpanded.
Definition: Data.cpp:1183
escript::AbstractSystemMatrix::getColumnBlockSize
int getColumnBlockSize() const
returns the column block size
Definition: AbstractSystemMatrix.h:135
weipa::IntVec
std::vector< int > IntVec
Definition: weipa.h:58
weipa::SpeckleyDomain::getMeshForFunctionSpace
virtual NodeData_ptr getMeshForFunctionSpace(int fsCode) const
Returns the node mesh for given function space code.
Definition: weipa/src/SpeckleyDomain.cpp:79
speckley::SpeckleyDomain::getSystemMatrixTypeId
virtual int getSystemMatrixTypeId(const boost::python::object &options) const
returns the identifier of the matrix type to be used for the global stiffness matrix when a particula...
Definition: speckley/src/SpeckleyDomain.cpp:703
escript::ValueError
An exception class that signals an invalid argument value.
Definition: EsysException.h:99
speckley::DiracPoint::tag
int tag
Definition: speckley/src/SpeckleyDomain.h:86
weipa::DataVar
A class that provides functionality to read an escript data object from a dump file or an escript::Da...
Definition: DataVar.h:34
escript::AbstractContinuousDomain
AbstractContinuousDomain, base class for continuous domains.
Definition: AbstractContinuousDomain.h:56
speckley::SpeckleyException
SpeckleyException exception class.
Definition: SpeckleyException.h:40
speckley::_brick
escript::Domain_ptr _brick(int order, double _n0, double _n1, double _n2, const object &l0, const object &l1, const object &l2, int d0, int d1, int d2, const object &objpoints, const object &objtags, escript::SubWorld_ptr world)
Definition: speckleycpp.cpp:141
escript::DataTypes::IndexVector
std::vector< index_t > IndexVector
Definition: DataTypes.h:85
weipa::SpeckleyNodes
Stores and manipulates speckley mesh nodes.
Definition: SpeckleyNodes.h:37
escript::AbstractContinuousDomain::addPDEToTransportProblem
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:165
speckley::SpeckleyDomain::addToSystemFromPython
virtual void addToSystemFromPython(escript::AbstractSystemMatrix &mat, escript::Data &rhs, const boost::python::list &data, Assembler_ptr assembler) const
a wrapper for addToSystem that allows calling from Python
Definition: speckley/src/SpeckleyDomain.cpp:728
speckley::Nodes
Definition: Speckley.h:72
speckley::extractPyArray
std::vector< T > extractPyArray(const object &obj, const std::string &name, int expectedLength=0)
Definition: speckleycpp.cpp:34