18#ifndef __RIPLEY_RECTANGLE_H__
19#define __RIPLEY_RECTANGLE_H__
22#define _USE_MATH_DEFINES
26#include <ripley/RipleyDomain.h>
48 Rectangle(dim_t n0, dim_t n1,
double x0,
double y0,
double x1,
double y1,
50 const std::vector<double>& points = std::vector<double>(),
51 const std::vector<int>& tags = std::vector<int>(),
66 virtual std::string getDescription()
const;
78 virtual void write(
const std::string& filename)
const;
85 void dump(
const std::string& filename)
const;
105 std::string filename,
106 int byteOrder,
int dataType)
const;
113 const dim_t* borrowSampleReferenceIDs(
int fsType)
const;
119 virtual bool ownSample(
int fsType, index_t
id)
const;
140 virtual dim_t getNumDataPointsGlobal()
const;
147 virtual void Print_Mesh_Info(
const bool full=
false)
const;
184 virtual double getLocalCoordinate(index_t index,
int dim)
const;
190 virtual boost::python::tuple getGridParameters()
const;
198 const boost::python::tuple& filter)
const;
204 virtual Assembler_ptr createAssembler(std::string type,
224 virtual RankVector getOwnerVector(
int fsType)
const;
227 virtual dim_t getNumNodes()
const;
228 virtual dim_t getNumElements()
const;
229 virtual dim_t getNumFaceElements()
const;
230 virtual dim_t getNumDOF()
const;
231 virtual dim_t getNumDOFInAxis(
unsigned axis)
const;
232 virtual index_t getFirstInDim(
unsigned axis)
const;
234 virtual IndexVector getDiagonalIndices(
bool upperOnly)
const;
238 virtual void assembleIntegrate(std::vector<real_t>& integrals,
240 virtual void assembleIntegrate(std::vector<cplx_t>& integrals,
242 virtual std::vector<IndexVector> getConnections(
bool includeShared=
false)
const;
244#ifdef ESYS_HAVE_TRILINOS
245 virtual esys_trilinos::const_TrilinosGraph_ptr getTrilinosGraph()
const;
250 bool reducedRowOrder,
bool reducedColOrder)
const;
259 virtual dim_t getDofOfNode(dim_t node)
const;
261 virtual void populateSampleIds();
262 virtual void populateDofMap();
264 template<
typename Scalar>
268 template<
typename Scalar>
270 const std::vector<Scalar>& EM_S,
const std::vector<Scalar>& EM_F,
271 bool addS,
bool addF, index_t firstNode,
int nEq=1,
int nComp=1)
const;
273 template<
typename ValueType>
274 void readBinaryGridImpl(
escript::Data& out,
const std::string& filename,
277#ifdef ESYS_HAVE_BOOST_IO
278 template<
typename ValueType>
283 template<
typename ValueType>
285 const std::string& filename,
int byteOrder)
const;
287 virtual dim_t findNode(
const double *coords)
const;
291 long seed,
const boost::python::tuple& filter)
const;
321 dim_t m_faceCount[4];
344#ifdef ESYS_HAVE_TRILINOS
346 mutable esys_trilinos::const_TrilinosGraph_ptr m_graph;
349 template<
typename Scalar>
350 void assembleIntegrateImpl(std::vector<Scalar>& integrals,
const escript::Data& arg)
const;
352 template <
typename S>
355 template <
typename S>
358 bool reduced,
S sentinel)
const;
374 ESYS_ASSERT(dim>=0 && dim<2,
"'dim' out of bounds");
381 return boost::python::make_tuple(
383 boost::python::make_tuple(
m_dx[0],
m_dx[1]),
384 boost::python::make_tuple(
m_gNE[0],
m_gNE[1]));
#define ESYS_ASSERT(a, b)
EsysAssert is a MACRO that will throw an exception if the boolean condition specified is false.
Definition Assert.h:79
#define S(_J_, _I_)
Definition ShapeFunctions.cpp:122
Base class for all escript domains.
Definition AbstractDomain.h:51
Base class for escript system matrices.
Definition AbstractSystemMatrix.h:44
Data represents a collection of datapoints.
Definition Data.h:64
Definition FunctionSpace.h:36
Definition ripley/src/DefaultAssembler2D.h:26
Definition LameAssembler2D.h:26
Rectangle is the 2-dimensional implementation of a RipleyDomain.
Definition ripley/src/Rectangle.h:35
virtual const dim_t * getNumElementsPerDim() const
returns the number of elements per MPI rank in each dimension
Definition ripley/src/Rectangle.h:159
virtual dim_t getNumNodes() const
returns the number of nodes per MPI rank
Definition ripley/src/Rectangle.h:401
virtual dim_t getNumFaceElements() const
returns the number of face elements on current MPI rank
Definition ripley/src/Rectangle.h:413
dim_t m_faceCount[4]
number of face elements per edge (left, right, bottom, top)
Definition ripley/src/Rectangle.h:321
virtual const int * getNumSubdivisionsPerDim() const
returns the number of spatial subdivisions in each dimension
Definition ripley/src/Rectangle.h:178
virtual dim_t getNumDOF() const
returns the number of degrees of freedom per MPI rank
Definition ripley/src/Rectangle.h:388
double m_origin[2]
origin of domain
Definition ripley/src/Rectangle.h:297
dim_t m_NN[2]
number of nodes for this rank in each dimension
Definition ripley/src/Rectangle.h:315
dim_t m_offset[2]
first node on this rank is at (offset0,offset1) in global mesh
Definition ripley/src/Rectangle.h:318
const double * getLength() const
returns the lengths of the domain
Definition ripley/src/Rectangle.h:211
IndexVector m_faceOffset
Definition ripley/src/Rectangle.h:325
IndexVector m_dofId
vector of sample reference identifiers
Definition ripley/src/Rectangle.h:328
virtual escript::Data randomFill(const escript::DataTypes::ShapeType &shape, const escript::FunctionSpace &what, long seed, const boost::python::tuple &filter) const
returns a Data object filled with random data passed through filter.
virtual dim_t getNumDOFInAxis(unsigned axis) const
Definition ripley/src/Rectangle.h:394
int m_NX[2]
number of spatial subdivisions
Definition ripley/src/Rectangle.h:306
virtual dim_t getNumElements() const
returns the number of elements per MPI rank
Definition ripley/src/Rectangle.h:407
const double * getElementLength() const
returns the lengths of an element
Definition ripley/src/Rectangle.h:217
virtual const dim_t * getNumNodesPerDim() const
returns the number of nodes per MPI rank in each dimension
Definition ripley/src/Rectangle.h:153
virtual dim_t getNumDataPointsGlobal() const
returns the number of data points summed across all MPI processes
Definition ripley/src/Rectangle.h:367
IndexVector m_nodeId
Definition ripley/src/Rectangle.h:329
virtual IndexVector getNodeDistribution() const
returns the node distribution vector
Definition ripley/src/Rectangle.h:172
dim_t m_NE[2]
number of elements for this rank in each dimension including shared
Definition ripley/src/Rectangle.h:309
double m_dx[2]
grid spacings / cell sizes of domain
Definition ripley/src/Rectangle.h:303
virtual const dim_t * getNumFacesPerBoundary() const
returns the number of face elements in the order (left,right,bottom,top) on current MPI rank
Definition ripley/src/Rectangle.h:166
IndexVector m_dofMap
Definition ripley/src/Rectangle.h:337
IndexVector m_elementId
Definition ripley/src/Rectangle.h:330
virtual dim_t getDofOfNode(dim_t node) const
Definition ripley/src/Rectangle.h:362
virtual index_t getFirstInDim(unsigned axis) const
Definition ripley/src/Rectangle.h:419
virtual double getLocalCoordinate(index_t index, int dim) const
returns the index'th coordinate value in given dimension for this rank
Definition ripley/src/Rectangle.h:372
dim_t m_gNE[2]
total number of elements in each dimension
Definition ripley/src/Rectangle.h:294
virtual boost::python::tuple getGridParameters() const
returns the tuple (origin, spacing, number_of_elements)
Definition ripley/src/Rectangle.h:379
IndexVector m_faceId
Definition ripley/src/Rectangle.h:331
IndexVector m_nodeDistribution
Definition ripley/src/Rectangle.h:334
RipleyDomain extends the AbstractContinuousDomain interface for the Ripley library and is the base cl...
Definition ripley/src/RipleyDomain.h:103
int m_numDim
Definition ripley/src/RipleyDomain.h:771
Definition ripley/src/WaveAssembler2D.h:26
std::vector< int > ShapeType
The shape of a single datapoint.
Definition DataTypes.h:44
boost::shared_ptr< SubWorld > SubWorld_ptr
Definition SubWorld.h:147
boost::shared_ptr< SystemMatrixPattern > SystemMatrixPattern_ptr
Definition SystemMatrixPattern.h:41
Definition ripley/src/AbstractAssembler.h:26
std::vector< index_t > IndexVector
Definition Ripley.h:44
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 ripleycpp.cpp:63
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 ripleycpp.cpp:117
std::map< std::string, int > TagMap
Definition Ripley.h:47
std::vector< int > RankVector
Definition Ripley.h:46
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 ripleycpp.cpp:88
std::map< std::string, escript::Data > DataMap
Definition ripley/src/domainhelpers.h:25
#define RIPLEY_DLL_API
Definition ripley/src/system_dep.h:21
Structure that wraps parameters for the grid reading routines.
Definition ripley/src/RipleyDomain.h:70