escript Revision_
ripley/src/Brick.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 __RIPLEY_BRICK_H__
19#define __RIPLEY_BRICK_H__
20
21#ifdef _WIN32 // for M_1_PI
22#define _USE_MATH_DEFINES
23#include <math.h>
24#endif
25
26#include <ripley/RipleyDomain.h>
27
28namespace ripley {
29
35{
36 template<class Scalar> friend class DefaultAssembler3D;
37 friend class WaveAssembler3D;
38 friend class LameAssembler3D;
39public:
40
48 Brick(dim_t n0, dim_t n1, dim_t n2, double x0, double y0, double z0,
49 double x1, double y1, double z1, int d0=-1, int d1=-1, int d2=-1,
50 const std::vector<double>& points = std::vector<double>(),
51 const std::vector<int>& tags = std::vector<int>(),
52 const TagMap& tagnamestonums = TagMap(),
54
59 ~Brick();
60
65 virtual std::string getDescription() const;
66
70 virtual bool operator==(const escript::AbstractDomain& other) const;
71
77 virtual void write(const std::string& filename) const;
78
84 void dump(const std::string& filename) const;
85
88 virtual void readNcGrid(escript::Data& out, std::string filename,
89 std::string varname, const ReaderParameters& params) const;
90
93 virtual void readBinaryGrid(escript::Data& out, std::string filename,
94 const ReaderParameters& params) const;
95
96 virtual void readBinaryGridFromZipped(escript::Data& out, std::string filename,
97 const ReaderParameters& params) const;
98
101 virtual void writeBinaryGrid(const escript::Data& in,
102 std::string filename,
103 int byteOrder, int dataType) const;
104
110 const dim_t* borrowSampleReferenceIDs(int fsType) const;
111
116 virtual bool ownSample(int fsType, index_t id) const;
117
124 virtual void setToNormal(escript::Data& out) const;
125
131 virtual void setToSize(escript::Data& out) const;
132
137 virtual dim_t getNumDataPointsGlobal() const;
138
144 virtual void Print_Mesh_Info(const bool full=false) const;
145
150 virtual const dim_t* getNumNodesPerDim() const { return m_NN; }
151
156 virtual const dim_t* getNumElementsPerDim() const { return m_NE; }
157
163 virtual const dim_t* getNumFacesPerBoundary() const { return m_faceCount; }
164
169 virtual IndexVector getNodeDistribution() const { return m_nodeDistribution; }
170
175 virtual const int* getNumSubdivisionsPerDim() const { return m_NX; }
176
181 virtual double getLocalCoordinate(index_t index, int dim) const;
182
187 virtual boost::python::tuple getGridParameters() const;
188
194 const escript::FunctionSpace& what,
195 long seed,
196 const boost::python::tuple& filter) const;
197
198 virtual Assembler_ptr createAssembler(std::string type,
199 const DataMap& options) const;
204 const double *getLength() const { return m_length; }
205
210 const double *getElementLength() const { return m_dx; }
211
217 virtual RankVector getOwnerVector(int fsType) const;
218
219protected:
220 virtual dim_t getNumNodes() const;
221 virtual dim_t getNumElements() const;
222 virtual dim_t getNumFaceElements() const;
223 virtual dim_t getNumDOF() const;
224 virtual dim_t getNumDOFInAxis(unsigned axis) const;
225 virtual dim_t getFirstInDim(unsigned axis) const;
226
227 virtual IndexVector getDiagonalIndices(bool upperOnly) const;
228 virtual void assembleCoordinates(escript::Data& arg) const;
229 virtual void assembleGradient(escript::Data& out,
230 const escript::Data& in) const;
231 virtual void assembleIntegrate(std::vector<real_t>& integrals,
232 const escript::Data& arg) const;
233 virtual void assembleIntegrate(std::vector<cplx_t>& integrals,
234 const escript::Data& arg) const;
235 virtual std::vector<IndexVector> getConnections(bool includeShared=false) const;
236#ifdef ESYS_HAVE_TRILINOS
237 virtual esys_trilinos::const_TrilinosGraph_ptr getTrilinosGraph() const;
238#endif
239
240#ifdef ESYS_HAVE_PASO
241 virtual paso::SystemMatrixPattern_ptr getPasoMatrixPattern(
242 bool reducedRowOrder, bool reducedColOrder) const;
243#endif
244 virtual void interpolateNodesOnElements(escript::Data& out,
245 const escript::Data& in, bool reduced) const;
246 virtual void interpolateNodesOnFaces(escript::Data& out,
247 const escript::Data& in,
248 bool reduced) const;
249 virtual void nodesToDOF(escript::Data& out, const escript::Data& in) const;
250 virtual dim_t getDofOfNode(dim_t node) const;
251
252 void populateSampleIds();
253 void populateDofMap();
254
255 template<typename Scalar>
256 void assembleGradientImpl(escript::Data& out,
257 const escript::Data& in) const;
258
259 template<typename Scalar>
260 void addToMatrixAndRHS(escript::AbstractSystemMatrix* S, escript::Data& F,
261 const std::vector<Scalar>& EM_S, const std::vector<Scalar>& EM_F,
262 bool addS, bool addF, index_t firstNode, int nEq=1, int nComp=1) const;
263
264 template<typename ValueType>
265 void readBinaryGridImpl(escript::Data& out, const std::string& filename,
266 const ReaderParameters& params) const;
267#ifdef ESYS_HAVE_BOOST_IO
268 template<typename ValueType>
269 void readBinaryGridZippedImpl(escript::Data& out,
270 const std::string& filename,
271 const ReaderParameters& params) const;
272#endif
273 template<typename ValueType>
274 void writeBinaryGridImpl(const escript::Data& in,
275 const std::string& filename, int byteOrder) const;
276
277 dim_t findNode(const double *coords) const;
278
280 const escript::DataTypes::ShapeType& shape, long seed,
281 const boost::python::tuple& filter) const;
282
284 dim_t m_gNE[3];
285
287 double m_origin[3];
288
290 double m_length[3];
291
293 double m_dx[3];
294
296 int m_NX[3];
297
299 dim_t m_NE[3];
300
302 dim_t m_ownNE[3];
303
305 dim_t m_NN[3];
306
308 dim_t m_offset[3];
309
311 dim_t m_faceCount[6];
312
317
323
324 // vector with first node id on each rank
326
327 // vector that maps each node to a DOF index
329
330#ifdef ESYS_HAVE_PASO
331 // the Paso System Matrix pattern
332 mutable paso::SystemMatrixPattern_ptr m_pattern;
333#endif
334
335#ifdef ESYS_HAVE_TRILINOS
337 mutable esys_trilinos::const_TrilinosGraph_ptr m_graph;
338#endif
339private:
340 template<typename Scalar>
341 void assembleIntegrateImpl(std::vector<Scalar>& integrals, const escript::Data& arg) const;
342
343 template <typename S>
344 void interpolateNodesOnElementsWorker(escript::Data& out,
345 const escript::Data& in, bool reduced, S sentinel) const;
346 template <typename S>
347 void interpolateNodesOnFacesWorker(escript::Data& out,
348 const escript::Data& in,
349 bool reduced, S sentinel) const;
350};
351
353inline dim_t Brick::getDofOfNode(dim_t node) const
354{
355 return m_dofMap[node];
356}
357
359{
360 return (m_gNE[0]+1)*(m_gNE[1]+1)*(m_gNE[2]+1);
361}
362
363inline double Brick::getLocalCoordinate(index_t index, int dim) const
364{
365 ESYS_ASSERT(dim>=0 && dim<3, "'dim' out of bounds");
366 ESYS_ASSERT(index>=0 && index<m_NN[dim], "'index' out of bounds");
367 return m_origin[dim]+m_dx[dim]*(m_offset[dim]+index);
368}
369
370inline boost::python::tuple Brick::getGridParameters() const
371{
372 return boost::python::make_tuple(
373 boost::python::make_tuple(m_origin[0], m_origin[1], m_origin[2]),
374 boost::python::make_tuple(m_dx[0], m_dx[1], m_dx[2]),
375 boost::python::make_tuple(m_gNE[0], m_gNE[1], m_gNE[2]));
376}
377
378//protected
379inline dim_t Brick::getNumDOF() const
380{
381 return (m_gNE[0]+1)/m_NX[0]*(m_gNE[1]+1)/m_NX[1]*(m_gNE[2]+1)/m_NX[2];
382}
383
384//protected
385inline dim_t Brick::getNumDOFInAxis(unsigned axis) const
386{
387 ESYS_ASSERT(axis < m_numDim, "Invalid axis");
388 return (m_gNE[axis]+1)/m_NX[axis];
389}
390
391//protected
392inline dim_t Brick::getNumNodes() const
393{
394 return m_NN[0]*m_NN[1]*m_NN[2];
395}
396
397//protected
398inline dim_t Brick::getNumElements() const
399{
400 return m_NE[0]*m_NE[1]*m_NE[2];
401}
402
403//protected
404inline dim_t Brick::getNumFaceElements() const
405{
406 return m_faceCount[0] + m_faceCount[1] + m_faceCount[2]
407 + m_faceCount[3] + m_faceCount[4] + m_faceCount[5];
408}
409
410//protected
411inline index_t Brick::getFirstInDim(unsigned axis) const
412{
413 return m_offset[axis] == 0 ? 0 : 1;
414}
415
416} // end of namespace ripley
417
418#endif // __RIPLEY_BRICK_H__
419
#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
Brick is the 3-dimensional implementation of a RipleyDomain.
Definition ripley/src/Brick.h:35
virtual dim_t getNumElements() const
returns the number of elements per MPI rank
Definition ripley/src/Brick.h:398
int m_NX[3]
number of spatial subdivisions
Definition ripley/src/Brick.h:296
IndexVector m_nodeId
Definition ripley/src/Brick.h:320
virtual const dim_t * getNumNodesPerDim() const
returns the number of nodes per MPI rank in each dimension
Definition ripley/src/Brick.h:150
IndexVector m_faceId
Definition ripley/src/Brick.h:322
IndexVector m_dofId
vector of sample reference identifiers
Definition ripley/src/Brick.h:319
dim_t m_offset[3]
first node on this rank is at (offset0,offset1,offset2) in global mesh
Definition ripley/src/Brick.h:308
double m_origin[3]
origin of domain
Definition ripley/src/Brick.h:287
virtual const dim_t * getNumElementsPerDim() const
returns the number of elements per MPI rank in each dimension
Definition ripley/src/Brick.h:156
virtual IndexVector getNodeDistribution() const
returns the node distribution vector
Definition ripley/src/Brick.h:169
virtual escript::Data randomFillWorker(const escript::DataTypes::ShapeType &shape, long seed, const boost::python::tuple &filter) const
double m_dx[3]
grid spacings / cell sizes of domain
Definition ripley/src/Brick.h:293
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 getNumDataPointsGlobal() const
returns the number of data points summed across all MPI processes
Definition ripley/src/Brick.h:358
IndexVector m_faceOffset
Definition ripley/src/Brick.h:316
virtual boost::python::tuple getGridParameters() const
returns the tuple (origin, spacing, number_of_elements)
Definition ripley/src/Brick.h:370
const double * getElementLength() const
returns the lengths of an element
Definition ripley/src/Brick.h:210
IndexVector m_dofMap
Definition ripley/src/Brick.h:328
dim_t m_NN[3]
number of nodes for this rank in each dimension
Definition ripley/src/Brick.h:305
virtual dim_t getNumNodes() const
returns the number of nodes per MPI rank
Definition ripley/src/Brick.h:392
virtual dim_t getFirstInDim(unsigned axis) const
Definition ripley/src/Brick.h:411
dim_t m_gNE[3]
total number of elements in each dimension
Definition ripley/src/Brick.h:284
virtual const dim_t * getNumFacesPerBoundary() const
returns the number of face elements in the order (left,right,bottom,top,front,back) on current MPI ra...
Definition ripley/src/Brick.h:163
virtual const int * getNumSubdivisionsPerDim() const
returns the number of spatial subdivisions in each dimension
Definition ripley/src/Brick.h:175
dim_t m_NE[3]
number of elements for this rank in each dimension including shared
Definition ripley/src/Brick.h:299
dim_t m_faceCount[6]
number of face elements per edge (left, right, bottom, top, front, back)
Definition ripley/src/Brick.h:311
virtual dim_t getNumDOF() const
returns the number of degrees of freedom per MPI rank
Definition ripley/src/Brick.h:379
IndexVector m_nodeDistribution
Definition ripley/src/Brick.h:325
IndexVector m_elementId
Definition ripley/src/Brick.h:321
const double * getLength() const
returns the lengths of the domain
Definition ripley/src/Brick.h:204
virtual dim_t getNumFaceElements() const
returns the number of face elements on current MPI rank
Definition ripley/src/Brick.h:404
virtual double getLocalCoordinate(index_t index, int dim) const
returns the index'th coordinate value in given dimension for this rank
Definition ripley/src/Brick.h:363
virtual dim_t getDofOfNode(dim_t node) const
Definition ripley/src/Brick.h:353
virtual dim_t getNumDOFInAxis(unsigned axis) const
Definition ripley/src/Brick.h:385
Definition ripley/src/DefaultAssembler3D.h:26
Definition LameAssembler3D.h:26
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/WaveAssembler3D.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