escript Revision_
speckley/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 __Speckley_BRICK_H__
19#define __Speckley_BRICK_H__
20
21#include <speckley/SpeckleyDomain.h>
22
23namespace speckley {
24
25#ifdef USE_RIPLEY
26class RipleyCoupler; //forward declaration of coupler to avoid circles
27#endif
28
34{
35 friend class DefaultAssembler3D;
36 friend class WaveAssembler3D;
37public:
38
46 Brick(int order, dim_t n0, dim_t n1, dim_t n2, double x0, double y0, double z0, double x1,
47 double y1, double z1, int d0=-1, int d1=-1, int d2=-1,
48 const std::vector<double>& points = std::vector<double>(),
49 const std::vector<int>& tags = std::vector<int>(),
50 const TagMap& tagnamestonums = TagMap(),
52
57 ~Brick();
58
63 virtual std::string getDescription() const;
64
68 virtual bool operator==(const escript::AbstractDomain& other) const;
69
75 virtual void write(const std::string& filename) const;
76
82 void dump(const std::string& filename) const;
83
86 virtual void readNcGrid(escript::Data& out, std::string filename,
87 std::string varname, const ReaderParameters& params) const;
88
91 virtual void readBinaryGrid(escript::Data& out, std::string filename,
92 const ReaderParameters& params) const;
93
94 virtual void readBinaryGridFromZipped(escript::Data& out,
95 std::string filename, const ReaderParameters& params) const;
96
99 virtual void writeBinaryGrid(const escript::Data& in,
100 std::string filename,
101 int byteOrder, int dataType) const;
102
108 const dim_t* borrowSampleReferenceIDs(int fsType) const;
109
114 virtual bool ownSample(int fsType, index_t id) const;
115
122 virtual void setToNormal(escript::Data& out) const;
123
129 virtual void setToSize(escript::Data& out) const;
130
135 virtual dim_t getNumDataPointsGlobal() const;
136
142 virtual void Print_Mesh_Info(const bool full=false) const;
143
148 virtual const dim_t* getNumNodesPerDim() const { return m_NN; }
149
154 virtual const dim_t* getNumElementsPerDim() const { return m_NE; }
155
161 virtual const dim_t* getNumFacesPerBoundary() const { return m_faceCount; }
162
167 virtual IndexVector getNodeDistribution() const { return m_nodeDistribution; }
168
173 virtual const int* getNumSubdivisionsPerDim() const { return m_NX; }
174
179 virtual double getLocalCoordinate(index_t index, int dim) const;
180
185 virtual boost::python::tuple getGridParameters() const;
186
191 virtual escript::Data randomFill(const escript::DataTypes::ShapeType& shape,
192 const escript::FunctionSpace& what,
193 long seed,
194 const boost::python::tuple& filter) const;
195
201 virtual void interpolateAcross(escript::Data& target,
202 const escript::Data& source) const;
203
208 virtual bool probeInterpolationAcross(int, const escript::AbstractDomain&,
209 int) const;
210
215 const double *getLength() const { return m_length; }
216
217protected:
218 virtual dim_t getNumNodes() const;
219 virtual dim_t getNumElements() const;
220 virtual dim_t getNumDOF() const;
221 virtual void assembleCoordinates(escript::Data& arg) const;
222 virtual void assembleGradient(escript::Data& out,
223 const escript::Data& in) const;
224 virtual void assembleIntegrate(std::vector<real_t>& integrals,
225 const escript::Data& arg) const;
226 virtual void assembleIntegrate(std::vector<cplx_t>& integrals,
227 const escript::Data& arg) const;
228 virtual void interpolateNodesOnElements(escript::Data& out,
229 const escript::Data& in,
230 bool reduced) const;
231 virtual void interpolateElementsOnNodes(escript::Data& out,
232 const escript::Data& in) const;
233 virtual dim_t getDofOfNode(dim_t node) const;
234 Assembler_ptr createAssembler(std::string type, const DataMap& constants) const;
235 virtual void reduceElements(escript::Data& out, const escript::Data& in) const;
236#ifdef ESYS_MPI
237 virtual void balanceNeighbours(escript::Data& data, bool average) const;
238#endif
239
240private:
241 template<typename Scalar>
242 void gradient_order2(escript::Data&, const escript::Data&) const;
243 template<typename Scalar>
244 void gradient_order3(escript::Data&, const escript::Data&) const;
245 template<typename Scalar>
246 void gradient_order4(escript::Data&, const escript::Data&) const;
247 template<typename Scalar>
248 void gradient_order5(escript::Data&, const escript::Data&) const;
249 template<typename Scalar>
250 void gradient_order6(escript::Data&, const escript::Data&) const;
251 template<typename Scalar>
252 void gradient_order7(escript::Data&, const escript::Data&) const;
253 template<typename Scalar>
254 void gradient_order8(escript::Data&, const escript::Data&) const;
255 template<typename Scalar>
256 void gradient_order9(escript::Data&, const escript::Data&) const;
257 template<typename Scalar>
258 void gradient_order10(escript::Data&, const escript::Data&) const;
259
260 template<typename Scalar>
261 void reduction_order2(const escript::Data&, escript::Data&) const;
262 template<typename Scalar>
263 void reduction_order3(const escript::Data&, escript::Data&) const;
264 template<typename Scalar>
265 void reduction_order4(const escript::Data&, escript::Data&) const;
266 template<typename Scalar>
267 void reduction_order5(const escript::Data&, escript::Data&) const;
268 template<typename Scalar>
269 void reduction_order6(const escript::Data&, escript::Data&) const;
270 template<typename Scalar>
271 void reduction_order7(const escript::Data&, escript::Data&) const;
272 template<typename Scalar>
273 void reduction_order8(const escript::Data&, escript::Data&) const;
274 template<typename Scalar>
275 void reduction_order9(const escript::Data&, escript::Data&) const;
276 template<typename Scalar>
277 void reduction_order10(const escript::Data&, escript::Data&) const;
278
279 template<typename Scalar>
280 void integral_order2(std::vector<Scalar>&, const escript::Data&) const;
281 template<typename Scalar>
282 void integral_order3(std::vector<Scalar>&, const escript::Data&) const;
283 template<typename Scalar>
284 void integral_order4(std::vector<Scalar>&, const escript::Data&) const;
285 template<typename Scalar>
286 void integral_order5(std::vector<Scalar>&, const escript::Data&) const;
287 template<typename Scalar>
288 void integral_order6(std::vector<Scalar>&, const escript::Data&) const;
289 template<typename Scalar>
290 void integral_order7(std::vector<Scalar>&, const escript::Data&) const;
291 template<typename Scalar>
292 void integral_order8(std::vector<Scalar>&, const escript::Data&) const;
293 template<typename Scalar>
294 void integral_order9(std::vector<Scalar>&, const escript::Data&) const;
295 template<typename Scalar>
296 void integral_order10(std::vector<Scalar>&, const escript::Data&) const;
297
298 template<typename Scalar>
299 void assembleIntegrateWorker(std::vector<Scalar>& integrals,
300 const escript::Data& arg) const;
301
302 template<typename Scalar>
303 void interpolateNodesOnElementsWorker(escript::Data& out,
304 const escript::Data& in,
305 bool reduced) const;
306#ifdef ESYS_MPI
307 void setCornerNeighbours();
308 void shareEdges(escript::Data& out, int rx, int ry, int rz) const;
309 void shareFaces(escript::Data& out, int rx, int ry, int rz) const;
310 void shareCorners(escript::Data& out) const;
311#endif
312 /* \brief
313 interpolates the non-corner point values of an element
314 from the corner values
315 */
316 void interpolateFromCorners(escript::Data& out) const;
317
318 void populateSampleIds();
319
320 template<typename ValueType>
321 void readBinaryGridImpl(escript::Data& out, const std::string& filename,
322 const ReaderParameters& params) const;
323
324#ifdef ESYS_HAVE_BOOST_IO
325 template<typename ValueType>
326 void readBinaryGridZippedImpl(escript::Data& out,
327 const std::string& filename, const ReaderParameters& params) const;
328#endif
329
330 template<typename ValueType>
331 void writeBinaryGridImpl(const escript::Data& in,
332 const std::string& filename, int byteOrder) const;
333
334 dim_t findNode(const double *coords) const;
335
336
337 escript::Data randomFillWorker(const escript::DataTypes::ShapeType& shape,
338 long seed, const boost::python::tuple& filter) const;
339
340#ifdef ESYS_MPI
342 int neighbour_ranks[8];
343
345 bool neighbour_exists[8];
346#endif
347
349 dim_t m_gNE[3];
350
352 double m_origin[3];
353
355 double m_length[3];
356
358 double m_dx[3];
359
361 int m_NX[3];
362
364 dim_t m_NE[3];
365
367 dim_t m_NN[3];
368
370 dim_t m_offset[3];
371
373 dim_t m_faceCount[6];
374
375
380
381 // vector with first node id on each rank
383
384#ifdef USE_RIPLEY
385 mutable RipleyCoupler *coupler;
386#endif
387};
388
390inline dim_t Brick::getDofOfNode(dim_t node) const
391{
392 return m_nodeId[node];
393}
394
396{
397 return (m_gNE[0]*m_order+1)*(m_gNE[1]*m_order+1)*(m_gNE[2]*m_order+1);
398}
399
400inline double Brick::getLocalCoordinate(index_t index, int dim) const
401{
402 ESYS_ASSERT(dim>=0 && dim<m_numDim, "'dim' out of bounds");
403 ESYS_ASSERT(index>=0 && index<m_NN[dim], "'index' out of bounds");
404 return m_origin[dim] //origin
405 + m_dx[dim]*(m_offset[dim] + index/m_order //elements
406 + point_locations[m_order-2][index%m_order]); //quads
407}
408
409inline boost::python::tuple Brick::getGridParameters() const
410{
411 return boost::python::make_tuple(
412 boost::python::make_tuple(m_origin[0], m_origin[1], m_origin[2]),
413 boost::python::make_tuple(m_dx[0], m_dx[1], m_dx[2]),
414 boost::python::make_tuple(m_gNE[0], m_gNE[1], m_gNE[2]));
415}
416
417//protected
418inline dim_t Brick::getNumDOF() const //global points
419{
420 return getNumNodes();
421}
422
423//protected
424inline dim_t Brick::getNumNodes() const //points per rank
425{
426 return m_NN[0] * m_NN[1] * m_NN[2];
427}
428
429//protected
430inline dim_t Brick::getNumElements() const
431{
432 return m_NE[0]*m_NE[1]*m_NE[2];
433}
434
435} // end of namespace speckley
436
437#endif // __Speckley_BRICK_H__
438
#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
Base class for all escript domains.
Definition AbstractDomain.h:51
Data represents a collection of datapoints.
Definition Data.h:64
Definition FunctionSpace.h:36
Brick is the 3-dimensional implementation of a SpeckleyDomain.
Definition speckley/src/Brick.h:34
virtual dim_t getNumDataPointsGlobal() const
returns the number of data points summed across all MPI processes
Definition speckley/src/Brick.h:395
virtual const dim_t * getNumElementsPerDim() const
returns the number of elements per MPI rank in each dimension
Definition speckley/src/Brick.h:154
double m_dx[3]
grid spacings / cell sizes of domain
Definition speckley/src/Brick.h:358
const double * getLength() const
returns the lengths of the domain
Definition speckley/src/Brick.h:215
virtual double getLocalCoordinate(index_t index, int dim) const
returns the index'th coordinate value in given dimension for this rank
Definition speckley/src/Brick.h:400
virtual boost::python::tuple getGridParameters() const
returns the tuple (origin, spacing, number_of_elements)
Definition speckley/src/Brick.h:409
virtual dim_t getNumNodes() const
returns the number of nodes per MPI rank
Definition speckley/src/Brick.h:424
dim_t m_NE[3]
number of elements for this rank in each dimension including shared
Definition speckley/src/Brick.h:364
IndexVector m_nodeDistribution
Definition speckley/src/Brick.h:382
virtual dim_t getDofOfNode(dim_t node) const
Definition speckley/src/Brick.h:390
virtual const int * getNumSubdivisionsPerDim() const
returns the number of spatial subdivisions in each dimension
Definition speckley/src/Brick.h:173
IndexVector m_nodeId
Definition speckley/src/Brick.h:378
virtual dim_t getNumElements() const
returns the number of elements per MPI rank
Definition speckley/src/Brick.h:430
double m_origin[3]
origin of domain
Definition speckley/src/Brick.h:352
dim_t m_gNE[3]
total number of elements in each dimension
Definition speckley/src/Brick.h:349
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 speckley/src/Brick.h:161
IndexVector m_elementId
Definition speckley/src/Brick.h:379
virtual IndexVector getNodeDistribution() const
returns the node distribution vector
Definition speckley/src/Brick.h:167
dim_t m_offset[3]
first node on this rank is at (offset0,offset1,offset2) in global mesh
Definition speckley/src/Brick.h:370
IndexVector m_dofId
vector of sample reference identifiers
Definition speckley/src/Brick.h:377
virtual dim_t getNumDOF() const
returns the number of degrees of freedom per MPI rank
Definition speckley/src/Brick.h:418
dim_t m_NN[3]
number of nodes for this rank in each dimension
Definition speckley/src/Brick.h:367
virtual const dim_t * getNumNodesPerDim() const
returns the number of nodes per MPI rank in each dimension
Definition speckley/src/Brick.h:148
Definition speckley/src/DefaultAssembler3D.h:33
Definition CrossDomainCoupler.h:29
SpeckleyDomain extends the AbstractContinuousDomain interface for the Speckley library and is the bas...
Definition speckley/src/SpeckleyDomain.h:86
int m_numDim
Definition speckley/src/SpeckleyDomain.h:745
int m_order
element order (will be m_order + 1 quad points in each axis)
Definition speckley/src/SpeckleyDomain.h:755
Definition speckley/src/WaveAssembler3D.h:33
std::vector< int > ShapeType
The shape of a single datapoint.
Definition DataTypes.h:44
boost::shared_ptr< SubWorld > SubWorld_ptr
Definition SubWorld.h:147
Definition AbstractAssembler.cpp:19
bool probeInterpolationAcross(int fsType_source, const escript::AbstractDomain &domain, int fsType_target, int dim)
Definition CrossDomainCoupler.cpp:32
std::map< std::string, int > TagMap
Definition Speckley.h:46
std::vector< index_t > IndexVector
Definition Speckley.h:44
const double point_locations[][11]
Definition Speckley.h:71
std::map< std::string, escript::Data > DataMap
Definition speckley/src/domainhelpers.h:25
#define Speckley_DLL_API
Definition speckley/src/system_dep.h:23
Structure that wraps parameters for the grid reading routines.
Definition speckley/src/SpeckleyDomain.h:53