escript  Revision_
ripley/src/Rectangle.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_RECTANGLE_H__
19 #define __RIPLEY_RECTANGLE_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 
28 namespace ripley {
29 
35 {
36  template<class Scalar> friend class DefaultAssembler2D;
37  friend class WaveAssembler2D;
38  friend class LameAssembler2D;
39 public:
40 
48  Rectangle(dim_t n0, dim_t n1, double x0, double y0, double x1, double y1,
49  int d0=-1, int d1=-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  );
55 
60  ~Rectangle();
61 
66  virtual std::string getDescription() const;
67 
71  virtual bool operator==(const escript::AbstractDomain& other) const;
72 
78  virtual void write(const std::string& filename) const;
79 
85  void dump(const std::string& filename) const;
86 
89  virtual void readNcGrid(escript::Data& out, std::string filename,
90  std::string varname, const ReaderParameters& params) const;
91 
94  virtual void readBinaryGrid(escript::Data& out, std::string filename,
95  const ReaderParameters& params) const;
96 
99  virtual void readBinaryGridFromZipped(escript::Data& out, std::string filename,
100  const ReaderParameters& params) const;
101 
104  virtual void writeBinaryGrid(const escript::Data& in,
105  std::string filename,
106  int byteOrder, int dataType) const;
107 
113  const dim_t* borrowSampleReferenceIDs(int fsType) const;
114 
119  virtual bool ownSample(int fsType, index_t id) const;
120 
127  virtual void setToNormal(escript::Data& out) const;
128 
134  virtual void setToSize(escript::Data& out) const;
135 
140  virtual dim_t getNumDataPointsGlobal() const;
141 
147  virtual void Print_Mesh_Info(const bool full=false) const;
148 
153  virtual const dim_t* getNumNodesPerDim() const { return m_NN; }
154 
159  virtual const dim_t* getNumElementsPerDim() const { return m_NE; }
160 
166  virtual const dim_t* getNumFacesPerBoundary() const { return m_faceCount; }
167 
172  virtual IndexVector getNodeDistribution() const { return m_nodeDistribution; }
173 
178  virtual const int* getNumSubdivisionsPerDim() const { return m_NX; }
179 
184  virtual double getLocalCoordinate(index_t index, int dim) const;
185 
190  virtual boost::python::tuple getGridParameters() const;
191 
197  const escript::FunctionSpace& what, long seed,
198  const boost::python::tuple& filter) const;
199 
204  virtual Assembler_ptr createAssembler(std::string type,
205  const DataMap& options) const;
206 
211  const double *getLength() const { return m_length; }
212 
217  const double *getElementLength() const { return m_dx; }
218 
224  virtual RankVector getOwnerVector(int fsType) const;
225 
226 protected:
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;
233 
234  virtual IndexVector getDiagonalIndices(bool upperOnly) const;
235  virtual void assembleCoordinates(escript::Data& arg) const;
236  virtual void assembleGradient(escript::Data& out,
237  const escript::Data& in) const;
238  virtual void assembleIntegrate(std::vector<real_t>& integrals,
239  const escript::Data& arg) const;
240  virtual void assembleIntegrate(std::vector<cplx_t>& integrals,
241  const escript::Data& arg) const;
242  virtual std::vector<IndexVector> getConnections(bool includeShared=false) const;
243 
244 #ifdef ESYS_HAVE_TRILINOS
245  virtual esys_trilinos::const_TrilinosGraph_ptr getTrilinosGraph() const;
246 #endif
247 
248 #ifdef ESYS_HAVE_PASO
249  virtual paso::SystemMatrixPattern_ptr getPasoMatrixPattern(
250  bool reducedRowOrder, bool reducedColOrder) const;
251 #endif
252  virtual void interpolateNodesOnElements(escript::Data& out,
253  const escript::Data& in, bool reduced) const;
254  virtual void interpolateNodesOnFaces(escript::Data& out,
255  const escript::Data& in,
256  bool reduced) const;
257 
258  virtual void nodesToDOF(escript::Data& out, const escript::Data& in) const;
259  virtual dim_t getDofOfNode(dim_t node) const;
260 
261  virtual void populateSampleIds();
262  virtual void populateDofMap();
263 
264  template<typename Scalar>
265  void assembleGradientImpl(escript::Data& out,
266  const escript::Data& in) const;
267 
268  template<typename Scalar>
269  void addToMatrixAndRHS(escript::AbstractSystemMatrix* S, escript::Data& F,
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;
272 
273  template<typename ValueType>
274  void readBinaryGridImpl(escript::Data& out, const std::string& filename,
275  const ReaderParameters& params) const;
276 
277 #ifdef ESYS_HAVE_BOOST_IO
278  template<typename ValueType>
279  void readBinaryGridZippedImpl(escript::Data& out,
280  const std::string& filename, const ReaderParameters& params) const;
281 #endif
282 
283  template<typename ValueType>
284  void writeBinaryGridImpl(const escript::Data& in,
285  const std::string& filename, int byteOrder) const;
286 
287  virtual dim_t findNode(const double *coords) const;
288 
289 
290  escript::Data randomFillWorker(const escript::DataTypes::ShapeType& shape,
291  long seed, const boost::python::tuple& filter) const;
292 
294  dim_t m_gNE[2];
295 
297  double m_origin[2];
298 
300  double m_length[2];
301 
303  double m_dx[2];
304 
306  int m_NX[2];
307 
309  dim_t m_NE[2];
310 
312  dim_t m_ownNE[2];
313 
315  dim_t m_NN[2];
316 
318  dim_t m_offset[2];
319 
321  dim_t m_faceCount[4];
322 
326 
332 
333  // vector with first node id on each rank
335 
336  // vector that maps each node to a DOF index (used for the coupler)
338 
339 #ifdef ESYS_HAVE_PASO
340  // the Paso System Matrix pattern
341  mutable paso::SystemMatrixPattern_ptr m_pattern;
342 #endif
343 
344 #ifdef ESYS_HAVE_TRILINOS
346  mutable esys_trilinos::const_TrilinosGraph_ptr m_graph;
347 #endif
348 private:
349  template<typename Scalar>
350  void assembleIntegrateImpl(std::vector<Scalar>& integrals, const escript::Data& arg) const;
351 
352  template <typename S>
353  void interpolateNodesOnElementsWorker(escript::Data& out,
354  const escript::Data& in, bool reduced, S sentinel) const;
355  template <typename S>
356  void interpolateNodesOnFacesWorker(escript::Data& out,
357  const escript::Data& in,
358  bool reduced, S sentinel) const;
359 };
360 
363 {
364  return m_dofMap[node];
365 }
366 
368 {
369  return (m_gNE[0]+1)*(m_gNE[1]+1);
370 }
371 
372 inline double Rectangle::getLocalCoordinate(index_t index, int dim) const
373 {
374  ESYS_ASSERT(dim>=0 && dim<2, "'dim' out of bounds");
375  ESYS_ASSERT(index>=0 && index<m_NN[dim], "'index' out of bounds");
376  return m_origin[dim]+m_dx[dim]*(m_offset[dim]+index);
377 }
378 
379 inline boost::python::tuple Rectangle::getGridParameters() const
380 {
381  return boost::python::make_tuple(
382  boost::python::make_tuple(m_origin[0], m_origin[1]),
383  boost::python::make_tuple(m_dx[0], m_dx[1]),
384  boost::python::make_tuple(m_gNE[0], m_gNE[1]));
385 }
386 
387 //protected
389 {
390  return (m_gNE[0]+1)/m_NX[0]*(m_gNE[1]+1)/m_NX[1];
391 }
392 
393 //protected
394 inline dim_t Rectangle::getNumDOFInAxis(unsigned axis) const
395 {
396  ESYS_ASSERT(axis < m_numDim, "Invalid axis");
397  return (m_gNE[axis]+1)/m_NX[axis];
398 }
399 
400 //protected
402 {
403  return m_NN[0]*m_NN[1];
404 }
405 
406 //protected
408 {
409  return m_NE[0]*m_NE[1];
410 }
411 
412 //protected
414 {
415  return m_faceCount[0] + m_faceCount[1] + m_faceCount[2] + m_faceCount[3];
416 }
417 
418 //protected
419 inline index_t Rectangle::getFirstInDim(unsigned axis) const
420 {
421  return m_offset[axis] == 0 ? 0 : 1;
422 }
423 
424 } // end of namespace ripley
425 
426 #endif // __RIPLEY_RECTANGLE_H__
427 
#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
const double * getLength() const
returns the lengths of the domain
Definition: ripley/src/Rectangle.h:211
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 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
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
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
const double * getElementLength() const
returns the lengths of an element
Definition: ripley/src/Rectangle.h:217
virtual IndexVector getNodeDistribution() const
returns the node distribution vector
Definition: ripley/src/Rectangle.h:172
virtual const dim_t * getNumNodesPerDim() const
returns the number of nodes per MPI rank in each dimension
Definition: ripley/src/Rectangle.h:153
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
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 const int * getNumSubdivisionsPerDim() const
returns the number of spatial subdivisions in each dimension
Definition: ripley/src/Rectangle.h:178
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
virtual const dim_t * getNumElementsPerDim() const
returns the number of elements per MPI rank in each dimension
Definition: ripley/src/Rectangle.h:159
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
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
index_t dim_t
Definition: DataTypes.h:66
int index_t
type for array/matrix indices used both globally and on each rank
Definition: DataTypes.h:61
boost::shared_ptr< SubWorld > SubWorld_ptr
Definition: SubWorld.h:147
boost::shared_ptr< SystemMatrixPattern > SystemMatrixPattern_ptr
Definition: SystemMatrixPattern.h:40
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