escript  Revision_
speckley/src/Brick.h
Go to the documentation of this file.
1 
2 /*****************************************************************************
3 *
4 * Copyright (c) 2003-2016 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_BRICK_H__
18 #define __Speckley_BRICK_H__
19 
20 #include <speckley/SpeckleyDomain.h>
21 
22 namespace speckley {
23 
24 #ifdef USE_RIPLEY
25 class RipleyCoupler; //forward declaration of coupler to avoid circles
26 #endif
27 
33 {
34  friend class DefaultAssembler3D;
35  friend class WaveAssembler3D;
36 public:
37 
45  Brick(int order, dim_t n0, dim_t n1, dim_t n2, double x0, double y0, double z0, double x1,
46  double y1, double z1, int d0=-1, int d1=-1, int d2=-1,
47  const std::vector<double>& points = std::vector<double>(),
48  const std::vector<int>& tags = std::vector<int>(),
49  const TagMap& tagnamestonums = TagMap(),
51 
56  ~Brick();
57 
62  virtual std::string getDescription() const;
63 
67  virtual bool operator==(const escript::AbstractDomain& other) const;
68 
74  virtual void write(const std::string& filename) const;
75 
81  void dump(const std::string& filename) const;
82 
85  virtual void readNcGrid(escript::Data& out, std::string filename,
86  std::string varname, const ReaderParameters& params) const;
87 
90  virtual void readBinaryGrid(escript::Data& out, std::string filename,
91  const ReaderParameters& params) const;
92 
93  virtual void readBinaryGridFromZipped(escript::Data& out,
94  std::string filename, const ReaderParameters& params) const;
95 
98  virtual void writeBinaryGrid(const escript::Data& in,
99  std::string filename,
100  int byteOrder, int dataType) const;
101 
107  const dim_t* borrowSampleReferenceIDs(int fsType) const;
108 
113  virtual bool ownSample(int fsType, index_t id) const;
114 
121  virtual void setToNormal(escript::Data& out) const;
122 
128  virtual void setToSize(escript::Data& out) const;
129 
134  virtual dim_t getNumDataPointsGlobal() const;
135 
141  virtual void Print_Mesh_Info(const bool full=false) const;
142 
147  virtual const dim_t* getNumNodesPerDim() const { return m_NN; }
148 
153  virtual const dim_t* getNumElementsPerDim() const { return m_NE; }
154 
160  virtual const dim_t* getNumFacesPerBoundary() const { return m_faceCount; }
161 
166  virtual IndexVector getNodeDistribution() const { return m_nodeDistribution; }
167 
172  virtual const int* getNumSubdivisionsPerDim() const { return m_NX; }
173 
178  virtual double getLocalCoordinate(index_t index, int dim) const;
179 
184  virtual boost::python::tuple getGridParameters() const;
185 
190  virtual escript::Data randomFill(const escript::DataTypes::ShapeType& shape,
191  const escript::FunctionSpace& what,
192  long seed,
193  const boost::python::tuple& filter) const;
194 
200  virtual void interpolateAcross(escript::Data& target,
201  const escript::Data& source) const;
202 
207  virtual bool probeInterpolationAcross(int, const escript::AbstractDomain&,
208  int) const;
209 
214  const double *getLength() const { return m_length; }
215 
216 protected:
217  virtual dim_t getNumNodes() const;
218  virtual dim_t getNumElements() const;
219  virtual dim_t getNumDOF() const;
220  virtual void assembleCoordinates(escript::Data& arg) const;
221  virtual void assembleGradient(escript::Data& out,
222  const escript::Data& in) const;
223  virtual void assembleIntegrate(DoubleVector& integrals,
224  const escript::Data& arg) const;
225  virtual void interpolateNodesOnElements(escript::Data& out,
226  const escript::Data& in,
227  bool reduced) const;
228  virtual void interpolateElementsOnNodes(escript::Data& out,
229  const escript::Data& in) const;
230  virtual dim_t getDofOfNode(dim_t node) const;
231  Assembler_ptr createAssembler(std::string type, const DataMap& constants) const;
232  virtual void reduceElements(escript::Data& out, const escript::Data& in) const;
233 #ifdef ESYS_MPI
234  virtual void balanceNeighbours(escript::Data& data, bool average) const;
235 #endif
236 
237 private:
238  void gradient_order2(escript::Data&, const escript::Data&) const;
239  void gradient_order3(escript::Data&, const escript::Data&) const;
240  void gradient_order4(escript::Data&, const escript::Data&) const;
241  void gradient_order5(escript::Data&, const escript::Data&) const;
242  void gradient_order6(escript::Data&, const escript::Data&) const;
243  void gradient_order7(escript::Data&, const escript::Data&) const;
244  void gradient_order8(escript::Data&, const escript::Data&) const;
245  void gradient_order9(escript::Data&, const escript::Data&) const;
246  void gradient_order10(escript::Data&, const escript::Data&) const;
247 
248  void reduction_order2(const escript::Data&, escript::Data&) const;
249  void reduction_order3(const escript::Data&, escript::Data&) const;
250  void reduction_order4(const escript::Data&, escript::Data&) const;
251  void reduction_order5(const escript::Data&, escript::Data&) const;
252  void reduction_order6(const escript::Data&, escript::Data&) const;
253  void reduction_order7(const escript::Data&, escript::Data&) const;
254  void reduction_order8(const escript::Data&, escript::Data&) const;
255  void reduction_order9(const escript::Data&, escript::Data&) const;
256  void reduction_order10(const escript::Data&, escript::Data&) const;
257 
258  void integral_order2(std::vector<double>&, const escript::Data&) const;
259  void integral_order3(std::vector<double>&, const escript::Data&) const;
260  void integral_order4(std::vector<double>&, const escript::Data&) const;
261  void integral_order5(std::vector<double>&, const escript::Data&) const;
262  void integral_order6(std::vector<double>&, const escript::Data&) const;
263  void integral_order7(std::vector<double>&, const escript::Data&) const;
264  void integral_order8(std::vector<double>&, const escript::Data&) const;
265  void integral_order9(std::vector<double>&, const escript::Data&) const;
266  void integral_order10(std::vector<double>&, const escript::Data&) const;
267 #ifdef ESYS_MPI
268  void setCornerNeighbours();
269  void shareEdges(escript::Data& out, int rx, int ry, int rz) const;
270  void shareFaces(escript::Data& out, int rx, int ry, int rz) const;
271  void shareCorners(escript::Data& out) const;
272 #endif
273  /* \brief
274  interpolates the non-corner point values of an element
275  from the corner values
276  */
277  void interpolateFromCorners(escript::Data& out) const;
278 
279  void populateSampleIds();
280  void addToMatrixAndRHS(escript::AbstractSystemMatrix* S, escript::Data& F,
281  const DoubleVector& EM_S, const DoubleVector& EM_F,
282  bool addS, bool addF, index_t firstNode, int nEq=1, int nComp=1) const;
283 
284  template<typename ValueType>
285  void readBinaryGridImpl(escript::Data& out, const std::string& filename,
286  const ReaderParameters& params) const;
287 
288 #ifdef ESYS_HAVE_BOOST_IO
289  template<typename ValueType>
290  void readBinaryGridZippedImpl(escript::Data& out,
291  const std::string& filename, const ReaderParameters& params) const;
292 #endif
293 
294  template<typename ValueType>
295  void writeBinaryGridImpl(const escript::Data& in,
296  const std::string& filename, int byteOrder) const;
297 
298  dim_t findNode(const double *coords) const;
299 
300 
301  escript::Data randomFillWorker(const escript::DataTypes::ShapeType& shape,
302  long seed, const boost::python::tuple& filter) const;
303 
304 #ifdef ESYS_MPI
305  int neighbour_ranks[8];
307 
309  bool neighbour_exists[8];
310 #endif
311 
313  dim_t m_gNE[3];
314 
316  double m_origin[3];
317 
319  double m_length[3];
320 
322  double m_dx[3];
323 
325  int m_NX[3];
326 
328  dim_t m_NE[3];
329 
331  dim_t m_NN[3];
332 
334  dim_t m_offset[3];
335 
337  dim_t m_faceCount[6];
338 
339 
344 
345  // vector with first node id on each rank
347 
348 #ifdef USE_RIPLEY
349  mutable RipleyCoupler *coupler;
350 #endif
351 };
352 
354 inline dim_t Brick::getDofOfNode(dim_t node) const
355 {
356  return m_nodeId[node];
357 }
358 
360 {
361  return (m_gNE[0]*m_order+1)*(m_gNE[1]*m_order+1)*(m_gNE[2]*m_order+1);
362 }
363 
364 inline double Brick::getLocalCoordinate(index_t index, int dim) const
365 {
366  ESYS_ASSERT(dim>=0 && dim<m_numDim, "'dim' out of bounds");
367  ESYS_ASSERT(index>=0 && index<m_NN[dim], "'index' out of bounds");
368  return m_origin[dim] //origin
369  + m_dx[dim]*(m_offset[dim] + index/m_order //elements
370  + point_locations[m_order-2][index%m_order]); //quads
371 }
372 
373 inline boost::python::tuple Brick::getGridParameters() const
374 {
375  return boost::python::make_tuple(
376  boost::python::make_tuple(m_origin[0], m_origin[1], m_origin[2]),
377  boost::python::make_tuple(m_dx[0], m_dx[1], m_dx[2]),
378  boost::python::make_tuple(m_gNE[0], m_gNE[1], m_gNE[2]));
379 }
380 
381 //protected
382 inline dim_t Brick::getNumDOF() const //global points
383 {
384  return getNumNodes();
385 }
386 
387 //protected
388 inline dim_t Brick::getNumNodes() const //points per rank
389 {
390  return m_NN[0] * m_NN[1] * m_NN[2];
391 }
392 
393 //protected
395 {
396  return m_NE[0]*m_NE[1]*m_NE[2];
397 }
398 
399 } // end of namespace speckley
400 
401 #endif // __Speckley_BRICK_H__
402 
Definition: FunctionSpace.h:34
Definition: AbstractAssembler.cpp:18
SpeckleyDomain extends the AbstractContinuousDomain interface for the Speckley library and is the bas...
Definition: speckley/src/SpeckleyDomain.h:84
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:85
bool probeInterpolationAcross(int fsType_source, const escript::AbstractDomain &domain, int fsType_target, int dim)
Definition: CrossDomainCoupler.cpp:31
IndexVector m_elementId
Definition: speckley/src/Brick.h:343
virtual dim_t getDofOfNode(dim_t node) const
Definition: speckley/src/Brick.h:354
std::vector< real_t > DoubleVector
Definition: Speckley.h:43
IndexVector m_dofId
vector of sample reference identifiers
Definition: speckley/src/Brick.h:341
virtual dim_t getNumDOF() const
returns the number of degrees of freedom per MPI rank
Definition: speckley/src/Brick.h:382
Definition: speckley/src/DefaultAssembler3D.h:31
virtual IndexVector getNodeDistribution() const
returns the node distribution vector
Definition: speckley/src/Brick.h:166
virtual const int * getNumSubdivisionsPerDim() const
returns the number of spatial subdivisions in each dimension
Definition: speckley/src/Brick.h:172
Definition: CrossDomainCoupler.h:28
std::map< std::string, escript::Data > DataMap
Definition: speckley/src/domainhelpers.h:24
virtual dim_t getNumDataPointsGlobal() const
returns the number of data points summed across all MPI processes
Definition: speckley/src/Brick.h:359
std::vector< int > ShapeType
The shape of a single datapoint.
Definition: DataTypes.h:42
#define Speckley_DLL_API
Definition: speckley/src/system_dep.h:22
const double * getLength() const
returns the lengths of the domain
Definition: speckley/src/Brick.h:214
Structure that wraps parameters for the grid reading routines.
Definition: speckley/src/SpeckleyDomain.h:51
int index_t
type for array/matrix indices used both globally and on each rank
Definition: DataTypes.h:59
boost::shared_ptr< SubWorld > SubWorld_ptr
Definition: SubWorld.h:146
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:114
Data represents a collection of datapoints.
Definition: Data.h:63
IndexVector m_nodeDistribution
Definition: speckley/src/Brick.h:346
std::vector< index_t > IndexVector
Definition: Speckley.h:42
virtual dim_t getNumElements() const
returns the number of elements per MPI rank
Definition: speckley/src/Brick.h:394
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:60
std::map< std::string, int > TagMap
Definition: Speckley.h:44
virtual boost::python::tuple getGridParameters() const
returns the tuple (origin, spacing, number_of_elements)
Definition: speckley/src/Brick.h:373
Brick is the 3-dimensional implementation of a SpeckleyDomain.
Definition: speckley/src/Brick.h:32
Base class for escript system matrices.
Definition: AbstractSystemMatrix.h:42
const double point_locations[][11]
Definition: Speckley.h:59
Definition: speckley/src/WaveAssembler3D.h:31
virtual dim_t getNumNodes() const
returns the number of nodes per MPI rank
Definition: speckley/src/Brick.h:388
virtual double getLocalCoordinate(index_t index, int dim) const
returns the index&#39;th coordinate value in given dimension for this rank
Definition: speckley/src/Brick.h:364
Base class for all escript domains.
Definition: AbstractDomain.h:45
#define ESYS_ASSERT(a, b)
EsysAssert is a MACRO that will throw an exception if the boolean condition specified is false...
Definition: Assert.h:78
virtual const dim_t * getNumNodesPerDim() const
returns the number of nodes per MPI rank in each dimension
Definition: speckley/src/Brick.h:147
IndexVector m_nodeId
Definition: speckley/src/Brick.h:342
#define S(_J_, _I_)
Definition: ShapeFunctions.cpp:121
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:160
index_t dim_t
Definition: DataTypes.h:64
virtual const dim_t * getNumElementsPerDim() const
returns the number of elements per MPI rank in each dimension
Definition: speckley/src/Brick.h:153