Nektar++
Geometry.h
Go to the documentation of this file.
1 ////////////////////////////////////////////////////////////////////////////////
2 //
3 // File: Geometry.h
4 //
5 // For more information, please see: http://www.nektar.info/
6 //
7 // The MIT License
8 //
9 // Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10 // Department of Aeronautics, Imperial College London (UK), and Scientific
11 // Computing and Imaging Institute, University of Utah (USA).
12 //
13 // Permission is hereby granted, free of charge, to any person obtaining a
14 // copy of this software and associated documentation files (the "Software"),
15 // to deal in the Software without restriction, including without limitation
16 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
17 // and/or sell copies of the Software, and to permit persons to whom the
18 // Software is furnished to do so, subject to the following conditions:
19 //
20 // The above copyright notice and this permission notice shall be included
21 // in all copies or substantial portions of the Software.
22 //
23 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29 // DEALINGS IN THE SOFTWARE.
30 //
31 // Description: This file contains the base class specification for the
32 // Geometry class.
33 //
34 //
35 ////////////////////////////////////////////////////////////////////////////////
36 
37 #ifndef NEKTAR_SPATIALDOMAINS_GEOMETRY_H
38 #define NEKTAR_SPATIALDOMAINS_GEOMETRY_H
39 
43 
44 #include <unordered_map>
45 #include <array>
46 
47 namespace Nektar
48 {
49 
50 namespace SpatialDomains
51 {
52 
53 class Geometry; // Forward declaration for typedef.
54 typedef std::shared_ptr<Geometry> GeometrySharedPtr;
55 typedef std::vector<GeometrySharedPtr> GeometryVector;
56 typedef std::unordered_set<GeometrySharedPtr> GeometrySet;
57 typedef std::shared_ptr<GeometryVector> GeometryVectorSharedPtr;
58 
59 class PointGeom;
60 typedef std::shared_ptr<PointGeom> PointGeomSharedPtr;
61 
62 class Geometry1D;
63 class Geometry2D;
64 typedef std::shared_ptr<Geometry1D> Geometry1DSharedPtr;
65 typedef std::shared_ptr<Geometry2D> Geometry2DSharedPtr;
66 
67 struct Curve;
68 typedef std::shared_ptr<Curve> CurveSharedPtr;
69 typedef std::unordered_map<int, CurveSharedPtr> CurveMap;
71 
72 /// \brief Less than operator to sort Geometry objects by global id when sorting
73 /// STL containers.
75  const std::shared_ptr<Geometry> &lhs,
76  const std::shared_ptr<Geometry> &rhs);
77 
79  const std::shared_ptr<Geometry> &lhs,
80  const std::shared_ptr<Geometry> &rhs);
81 
82 /// Base class for shape geometry information
83 class Geometry
84 {
85 public:
87  SPATIAL_DOMAINS_EXPORT Geometry(int coordim);
88 
90 
91  //---------------------------------------
92  // Helper functions
93  //---------------------------------------
94 
95  SPATIAL_DOMAINS_EXPORT inline int GetCoordim() const;
96  SPATIAL_DOMAINS_EXPORT inline void SetCoordim(int coordim);
97 
103 
104  //---------------------------------------
105  // Set and get ID
106  //---------------------------------------
107  SPATIAL_DOMAINS_EXPORT inline int GetGlobalID(void) const;
108  SPATIAL_DOMAINS_EXPORT inline void SetGlobalID(int globalid);
109 
110  //---------------------------------------
111  // Vertex, edge and face access
112  //---------------------------------------
113  SPATIAL_DOMAINS_EXPORT int GetVid(int i) const;
114  SPATIAL_DOMAINS_EXPORT int GetEid(int i) const;
115  SPATIAL_DOMAINS_EXPORT int GetFid(int i) const;
116  SPATIAL_DOMAINS_EXPORT inline int GetTid(int i) const;
121  const int i) const;
123  const int i) const;
124  SPATIAL_DOMAINS_EXPORT inline int GetNumVerts() const;
125  SPATIAL_DOMAINS_EXPORT inline int GetNumEdges() const;
126  SPATIAL_DOMAINS_EXPORT inline int GetNumFaces() const;
127  SPATIAL_DOMAINS_EXPORT inline int GetShapeDim() const;
128 
129  //---------------------------------------
130  // \chi mapping access
131  //---------------------------------------
134  SPATIAL_DOMAINS_EXPORT inline const
135  Array<OneD, const NekDouble> &GetCoeffs(const int i) const;
136  SPATIAL_DOMAINS_EXPORT inline void FillGeom();
137 
138  //---------------------------------------
139  // Point lookups
140  //---------------------------------------
141  SPATIAL_DOMAINS_EXPORT std::array<NekDouble, 6> GetBoundingBox();
142 
144  const Array<OneD, const NekDouble> &gloCoord,
145  NekDouble tol = 0.0);
147  const Array<OneD, const NekDouble> &gloCoord,
148  Array<OneD, NekDouble> &locCoord,
149  NekDouble tol);
151  const Array<OneD, const NekDouble> &gloCoord,
152  Array<OneD, NekDouble> &locCoord,
153  NekDouble tol,
154  NekDouble &dist);
156  const Array<OneD, const NekDouble> &coords,
157  Array<OneD, NekDouble> &Lcoords);
159  const int i, const Array<OneD, const NekDouble> &Lcoord);
161  const Array<OneD, const NekDouble> &gloCoord);
163  Array<OneD, NekDouble> &locCoord,
164  NekDouble tol);
165 
166  //---------------------------------------
167  // Misc. helper functions
168  //---------------------------------------
169  SPATIAL_DOMAINS_EXPORT inline int GetVertexEdgeMap(int i, int j) const;
170  SPATIAL_DOMAINS_EXPORT inline int GetVertexFaceMap(int i, int j) const;
171  SPATIAL_DOMAINS_EXPORT inline int GetEdgeFaceMap(int i, int j) const;
172  SPATIAL_DOMAINS_EXPORT inline int GetDir(const int i, const int j = 0) const;
173 
174 
175  SPATIAL_DOMAINS_EXPORT inline void Reset(CurveMap &curvedEdges,
176  CurveMap &curvedFaces);
177  SPATIAL_DOMAINS_EXPORT inline void Setup();
178 
179 protected:
181  GeomFactorsSharedPtr geomFactor);
183 
184  /// Coordinate dimension of this geometry object.
186  /// Geometric factors.
188  /// State of the geometric factors
190  /// \f$\chi\f$ mapping containing isoparametric transformation.
192  /// Enumeration to dictate whether coefficients are filled.
194  /// Wether or not the setup routines have been run
196  /// Type of geometry.
198  /// Type of shape.
200  /// Global ID
202  /// Array containing expansion coefficients of @p m_xmap
204  /// Array containing bounding box
206 
207  /// Handles generation of geometry factors.
208  void GenGeomFactors();
209 
210  //---------------------------------------
211  // Helper functions
212  //---------------------------------------
213  virtual PointGeomSharedPtr v_GetVertex(int i) const = 0;
214  virtual Geometry1DSharedPtr v_GetEdge(int i) const;
215  virtual Geometry2DSharedPtr v_GetFace(int i) const;
216  virtual StdRegions::Orientation v_GetEorient(const int i) const;
217  virtual StdRegions::Orientation v_GetForient(const int i) const;
218  virtual int v_GetNumVerts() const;
219  virtual int v_GetNumEdges() const;
220  virtual int v_GetNumFaces() const;
221  virtual int v_GetShapeDim() const;
222 
224  virtual void v_FillGeom();
225 
226  virtual bool v_ContainsPoint(const Array<OneD, const NekDouble> &gloCoord,
227  Array<OneD, NekDouble> &locCoord,
228  NekDouble tol,
229  NekDouble &dist);
230 
231  virtual NekDouble v_GetCoord(const int i,
232  const Array<OneD, const NekDouble> &Lcoord);
234  Array<OneD, NekDouble> &Lcoords);
235 
236  virtual int v_GetVertexEdgeMap(int i, int j) const;
237  virtual int v_GetVertexFaceMap(int i, int j) const;
238  virtual int v_GetEdgeFaceMap(int i, int j) const;
239  virtual int v_GetDir(const int faceidx, const int facedir) const;
240 
241  virtual void v_Reset(CurveMap &curvedEdges, CurveMap &curvedFaces);
242  virtual void v_Setup();
243  virtual void v_GenGeomFactors() = 0;
244 
245  inline void SetUpCoeffs(const int nCoeffs);
246 }; // class Geometry
247 
248 /**
249  * @brief Unary function that constructs a hash of a Geometry object, based on
250  * the vertex IDs.
251  */
252 struct GeometryHash : std::unary_function<GeometrySharedPtr, std::size_t>
253 {
254  std::size_t operator()(GeometrySharedPtr const& p) const
255  {
256  int i;
257  size_t seed = 0;
258  int nVert = p->GetNumVerts();
259  std::vector<unsigned int> ids(nVert);
260 
261  for (i = 0; i < nVert; ++i)
262  {
263  ids[i] = p->GetVid(i);
264  }
265  std::sort(ids.begin(), ids.end());
266  hash_range(seed, ids.begin(), ids.end());
267 
268  return seed;
269  }
270 };
271 
272 /**
273  * @brief Return the coordinate dimension of this object (i.e. the dimension of
274  * the space in which this object is embedded).
275  */
276 inline int Geometry::GetCoordim() const
277 {
278  return m_coordim;
279 }
280 
281 /**
282  * @brief Sets the coordinate dimension of this object (i.e. the dimension of
283  * the space in which this object is embedded).
284  */
285 inline void Geometry::SetCoordim(int dim)
286 {
287  m_coordim = dim;
288 }
289 
290 /**
291  * @brief Get the geometric factors for this object, generating them if
292  * required.
293  */
295 {
296  GenGeomFactors();
298 }
299 
300 /**
301  * @brief Get the geometric factors for this object.
302  */
304 {
305  return m_geomFactors;
306 }
307 
308 /**
309  * @brief Get the geometric shape type of this object.
310  */
312 {
313  return m_shapeType;
314 }
315 
316 /**
317  * @brief Get the ID of this object.
318  */
319 inline int Geometry::GetGlobalID(void) const
320 {
321  return m_globalID;
322 }
323 
324 /**
325  * @brief Set the ID of this object.
326  */
327 inline void Geometry::SetGlobalID(int globalid)
328 {
329  m_globalID = globalid;
330 }
331 
332 /**
333  * @brief Get the ID of trace @p i of this object.
334  *
335  * The trace element is the facet one dimension lower than the object; for
336  * example, a quadrilateral has four trace segments forming its boundary.
337  */
338 inline int Geometry::GetTid(int i) const
339 {
340  const int nDim = GetShapeDim();
341  return nDim == 1 ? GetVid(i) : nDim == 2 ? GetEid(i) :
342  nDim == 3 ? GetFid(i) : 0;
343 }
344 
345 /**
346  * @brief Returns vertex @p i of this object.
347  */
349 {
350  return v_GetVertex(i);
351 }
352 
353 /**
354  * @brief Returns edge @p i of this object.
355  */
357 {
358  return v_GetEdge(i);
359 }
360 
361 /**
362  * @brief Returns face @p i of this object.
363  */
365 {
366  return v_GetFace(i);
367 }
368 
369 /**
370  * @brief Returns the orientation of edge @p i with respect to the ordering of
371  * edges in the standard element.
372  */
374 {
375  return v_GetEorient(i);
376 }
377 
378 /**
379  * @brief Returns the orientation of face @p i with respect to the ordering of
380  * faces in the standard element.
381  */
383 {
384  return v_GetForient(i);
385 }
386 
387 /**
388  * @brief Get the number of vertices of this object.
389  */
390 inline int Geometry::GetNumVerts() const
391 {
392  return v_GetNumVerts();
393 }
394 
395 /**
396  * @brief Get the number of edges of this object.
397  */
398 inline int Geometry::GetNumEdges() const
399 {
400  return v_GetNumEdges();
401 }
402 
403 /**
404  * @brief Get the number of faces of this object.
405  */
406 inline int Geometry::GetNumFaces() const
407 {
408  return v_GetNumFaces();
409 }
410 
411 /**
412  * @brief Get the object's shape dimension.
413  *
414  * For example, a segment is one dimensional and quadrilateral is two
415  * dimensional.
416  */
417 inline int Geometry::GetShapeDim() const
418 {
419  return v_GetShapeDim();
420 }
421 
422 /**
423  * @brief Return the mapping object Geometry::m_xmap that represents the
424  * coordinate transformation from standard element to physical element.
425  */
427 {
428  return v_GetXmap();
429 }
430 
431 /**
432  * @brief Return the coefficients of the transformation Geometry::m_xmap in
433  * coordinate direction @p i.
434  */
436  const int i) const
437 {
438  return m_coeffs[i];
439 }
440 
441 /**
442  * @brief Populate the coordinate mapping Geometry::m_coeffs information from
443  * any children geometry elements.
444  *
445  * @see v_FillGeom()
446  */
447 inline void Geometry::FillGeom()
448 {
449  v_FillGeom();
450 }
451 
452 /**
453  * @brief Determine whether an element contains a particular Cartesian
454  * coordinate \f$(x,y,z)\f$.
455  *
456  * @see Geometry::ContainsPoint
457  */
459  const Array<OneD, const NekDouble> &gloCoord,
460  NekDouble tol)
461 {
462  Array<OneD,NekDouble> locCoord(GetCoordim(), 0.0);
463  NekDouble dist;
464  return v_ContainsPoint(gloCoord, locCoord, tol, dist);
465 }
466 
467 /**
468  * @brief Determine whether an element contains a particular Cartesian
469  * coordinate \f$(x,y,z)\f$.
470  *
471  * @see Geometry::ContainsPoint
472  */
474  const Array<OneD, const NekDouble> &gloCoord,
475  Array<OneD, NekDouble> &locCoord,
476  NekDouble tol)
477 {
478  NekDouble dist;
479  return v_ContainsPoint(gloCoord, locCoord, tol, dist);
480 }
481 
482 /**
483  * @brief Determine whether an element contains a particular Cartesian
484  * coordinate \f$\vec{x} = (x,y,z)\f$.
485  *
486  * For curvilinear and non-affine elements (i.e. where the Jacobian varies as a
487  * function of the standard element coordinates), this is a non-linear
488  * optimisation problem that requires the use of a Newton iteration. Note
489  * therefore that this can be an expensive operation.
490  *
491  * The parameter @p tol which is by default 0, can be used to expand the
492  * coordinate range of the standard element from \f$[-1,1]^d\f$ to
493  * \f$[-1-\epsilon,1+\epsilon\f$ to handle challenging edge cases. The function
494  * also returns the local coordinates corresponding to @p gloCoord that can be
495  * used to speed up subsequent searches.
496  *
497  * @param gloCoord The coordinate \f$ (x,y,z) \f$.
498  * @param locCoord On exit, this is the local coordinate \f$\vec{\xi}\f$ such
499  * that \f$\chi(\vec{\xi}) = \vec{x}\f$.
500  * @param tol The tolerance used to dictate the bounding box of the
501  * standard coordinates \f$\vec{\xi}\f$.
502  * @param dist On exit, returns the minimum distance between @p gloCoord
503  * and the quadrature points inside the element.
504  *
505  * @return `true` if the coordinate @p gloCoord is contained in the element;
506  * false otherwise.
507  *
508  * @see Geometry::GetLocCoords.
509  */
511  const Array<OneD, const NekDouble> &gloCoord,
512  Array<OneD, NekDouble> &locCoord,
513  NekDouble tol,
514  NekDouble &dist)
515 {
516  return v_ContainsPoint(gloCoord, locCoord, tol, dist);
517 }
518 
519 /**
520  * @brief Determine the local collapsed coordinates that correspond to a
521  * given Cartesian coordinate for this geometry object.
522  *
523  * For curvilinear and non-affine elements (i.e. where the Jacobian varies as a
524  * function of the standard element coordinates), this is a non-linear
525  * optimisation problem that requires the use of a Newton iteration. Note
526  * therefore that this can be an expensive operation.
527  *
528  * Note that, clearly, the provided Cartesian coordinate lie outside the
529  * element. The function therefore returns the minimum distance from some
530  * position in the element to . @p Lcoords will also be constrained to fit
531  * within the range \f$[-1,1]^d\f$ where \f$ d \f$ is the dimension of the
532  * element.
533  *
534  * @param coords Input Cartesian global coordinates
535  * @param Lcoords Corresponding local coordinates
536  *
537  * @return Distance between obtained coordinates and provided ones.
538  */
540  const Array<OneD, const NekDouble> &coords,
541  Array<OneD, NekDouble> &Lcoords)
542 {
543  return v_GetLocCoords(coords, Lcoords);
544 }
545 
546 /**
547  * @brief Given local collapsed coordinate @p Lcoord, return the value of
548  * physical coordinate in direction @p i.
549  */
550 inline NekDouble Geometry::GetCoord(const int i,
551  const Array<OneD, const NekDouble> &Lcoord)
552 {
553  return v_GetCoord(i, Lcoord);
554 }
555 
556 /**
557  * @brief Returns the standard element edge IDs that are connected to a given
558  * vertex.
559  *
560  * For example, on a prism, vertex 0 is connnected to edges 0, 3, and 4;
561  * `GetVertexEdgeMap(0,j)` would therefore return the values 0, 1 and 4
562  * respectively. We assume that @p j runs between 0 and 2 inclusive, which is
563  * true for every 3D element asides from the pyramid.
564  *
565  * This function is used in the construction of the low-energy preconditioner.
566  *
567  * @param i The vertex to query connectivity for.
568  * @param j The local edge index between 0 and 2 connected to this element.
569  *
570  * @todo Expand to work with pyramid elements.
571  * @see MultiRegions::PreconditionerLowEnergy
572  */
573 inline int Geometry::GetVertexEdgeMap(int i, int j) const
574 {
575  return v_GetVertexEdgeMap(i, j);
576 }
577 
578 /**
579  * @brief Returns the standard element face IDs that are connected to a given
580  * vertex.
581  *
582  * For example, on a hexahedron, vertex 0 is connnected to faces 0, 1, and 4;
583  * `GetVertexFaceMap(0,j)` would therefore return the values 0, 1 and 4
584  * respectively. We assume that @p j runs between 0 and 2 inclusive, which is
585  * true for every 3D element asides from the pyramid.
586  *
587  * This is used in the construction of the low-energy preconditioner.
588  *
589  * @param i The vertex to query connectivity for.
590  * @param j The local face index between 0 and 2 connected to this element.
591  *
592  * @todo Expand to work with pyramid elements.
593  * @see MultiRegions::PreconditionerLowEnergy
594  */
595 inline int Geometry::GetVertexFaceMap(int i, int j) const
596 {
597  return v_GetVertexFaceMap(i, j);
598 }
599 
600 /**
601  * @brief Returns the standard element edge IDs that are connected to a given
602  * face.
603  *
604  * For example, on a prism, edge 0 is connnected to faces 0 and 1;
605  * `GetEdgeFaceMap(0,j)` would therefore return the values 0 and 1
606  * respectively. We assume that @p j runs between 0 and 1 inclusive, since every
607  * face is connected to precisely two faces for all 3D elements.
608  *
609  * This function is used in the construction of the low-energy preconditioner.
610  *
611  * @param i The edge to query connectivity for.
612  * @param j The local face index between 0 and 1 connected to this element.
613  *
614  * @see MultiRegions::PreconditionerLowEnergy
615  */
616 inline int Geometry::GetEdgeFaceMap(int i, int j) const
617 {
618  return v_GetEdgeFaceMap(i, j);
619 }
620 
621 
622 /**
623  * @brief Returns the element coordinate direction corresponding to a given face
624  * coordinate direction
625  */
626 inline int Geometry::GetDir(const int faceidx, const int facedir) const
627 {
628  return v_GetDir(faceidx, facedir);
629 }
630 
631 /**
632  * @brief Reset this geometry object: unset the current state, zero
633  * Geometry::m_coeffs and remove allocated GeomFactors.
634  */
635 inline void Geometry::Reset(CurveMap &curvedEdges, CurveMap &curvedFaces)
636 {
637  v_Reset(curvedEdges, curvedFaces);
638 }
639 inline void Geometry::Setup()
640 {
641  v_Setup();
642 }
643 
644 /**
645  * @brief Generate the geometric factors (i.e. derivatives of \f$\chi\f$) and
646  * related metrics.
647  *
648  * @see SpatialDomains::GeomFactors
649  */
651 {
652  return v_GenGeomFactors();
653 }
654 
655 /**
656  * @brief Initialise the Geometry::m_coeffs array.
657  */
658 inline void Geometry::SetUpCoeffs(const int nCoeffs)
659 {
661 
662  for (int i = 0; i < m_coordim; ++i)
663  {
664  m_coeffs[i] = Array<OneD, NekDouble>(nCoeffs, 0.0);
665  }
666 }
667 
668 }
669 }
670 
671 #endif // NEKTAR_SPATIALDOMAINS_GEOMETRY_H
#define SPATIAL_DOMAINS_EXPORT
1D geometry information
Definition: Geometry1D.h:55
2D geometry information
Definition: Geometry2D.h:69
Base class for shape geometry information.
Definition: Geometry.h:84
virtual int v_GetNumEdges() const
Get the number of edges of this object.
Definition: Geometry.cpp:215
virtual StdRegions::Orientation v_GetForient(const int i) const
Returns the orientation of face i with respect to the ordering of faces in the standard element.
Definition: Geometry.cpp:204
NekDouble GetCoord(const int i, const Array< OneD, const NekDouble > &Lcoord)
Given local collapsed coordinate Lcoord, return the value of physical coordinate in direction i.
Definition: Geometry.h:550
LibUtilities::ShapeType GetShapeType(void)
Get the geometric shape type of this object.
Definition: Geometry.h:311
virtual NekDouble v_GetCoord(const int i, const Array< OneD, const NekDouble > &Lcoord)
Given local collapsed coordinate Lcoord, return the value of physical coordinate in direction i.
Definition: Geometry.cpp:322
Array< OneD, NekDouble > m_boundingBox
Array containing bounding box.
Definition: Geometry.h:205
GeomFactorsSharedPtr GetRefGeomFactors(const Array< OneD, const LibUtilities::BasisSharedPtr > &tbasis)
bool m_setupState
Wether or not the setup routines have been run.
Definition: Geometry.h:195
PointGeomSharedPtr GetVertex(int i) const
Returns vertex i of this object.
Definition: Geometry.h:348
GeomState m_state
Enumeration to dictate whether coefficients are filled.
Definition: Geometry.h:193
bool ClampLocCoords(Array< OneD, NekDouble > &locCoord, NekDouble tol)
Clamp local coords to be within standard regions [-1, 1]^dim.
Definition: Geometry.cpp:486
void SetUpCoeffs(const int nCoeffs)
Initialise the Geometry::m_coeffs array.
Definition: Geometry.h:658
NekDouble GetLocCoords(const Array< OneD, const NekDouble > &coords, Array< OneD, NekDouble > &Lcoords)
Determine the local collapsed coordinates that correspond to a given Cartesian coordinate for this ge...
Definition: Geometry.h:539
GeomType m_geomType
Type of geometry.
Definition: Geometry.h:197
Geometry1DSharedPtr GetEdge(int i) const
Returns edge i of this object.
Definition: Geometry.h:356
void SetGlobalID(int globalid)
Set the ID of this object.
Definition: Geometry.h:327
int GetNumFaces() const
Get the number of faces of this object.
Definition: Geometry.h:406
int GetShapeDim() const
Get the object's shape dimension.
Definition: Geometry.h:417
virtual void v_FillGeom()
Populate the coordinate mapping Geometry::m_coeffs information from any children geometry elements.
Definition: Geometry.cpp:346
Geometry()
Default constructor.
Definition: Geometry.cpp:54
virtual NekDouble v_GetLocCoords(const Array< OneD, const NekDouble > &coords, Array< OneD, NekDouble > &Lcoords)
Determine the local collapsed coordinates that correspond to a given Cartesian coordinate for this ge...
Definition: Geometry.cpp:334
virtual Geometry2DSharedPtr v_GetFace(int i) const
Returns face i of this object.
Definition: Geometry.cpp:172
int GetVid(int i) const
Get the ID of vertex i of this object.
Definition: Geometry.cpp:137
virtual int v_GetNumFaces() const
Get the number of faces of this object.
Definition: Geometry.cpp:223
virtual StdRegions::StdExpansionSharedPtr v_GetXmap() const
Return the mapping object Geometry::m_xmap that represents the coordinate transformation from standar...
Definition: Geometry.cpp:241
virtual void v_Reset(CurveMap &curvedEdges, CurveMap &curvedFaces)
Reset this geometry object: unset the current state, zero Geometry::m_coeffs and remove allocated Geo...
Definition: Geometry.cpp:355
const Array< OneD, const NekDouble > & GetCoeffs(const int i) const
Return the coefficients of the transformation Geometry::m_xmap in coordinate direction i.
Definition: Geometry.h:435
int GetGlobalID(void) const
Get the ID of this object.
Definition: Geometry.h:319
int GetCoordim() const
Return the coordinate dimension of this object (i.e. the dimension of the space in which this object ...
Definition: Geometry.h:276
int GetFid(int i) const
Get the ID of face i of this object.
Definition: Geometry.cpp:153
virtual int v_GetEdgeFaceMap(int i, int j) const
Returns the standard element edge IDs that are connected to a given face.
Definition: Geometry.cpp:300
LibUtilities::ShapeType m_shapeType
Type of shape.
Definition: Geometry.h:199
int GetDir(const int i, const int j=0) const
Returns the element coordinate direction corresponding to a given face coordinate direction.
Definition: Geometry.h:626
void FillGeom()
Populate the coordinate mapping Geometry::m_coeffs information from any children geometry elements.
Definition: Geometry.h:447
Array< OneD, Array< OneD, NekDouble > > m_coeffs
Array containing expansion coefficients of m_xmap.
Definition: Geometry.h:203
virtual ~Geometry()
Default destructor.
Definition: Geometry.cpp:72
int GetVertexEdgeMap(int i, int j) const
Returns the standard element edge IDs that are connected to a given vertex.
Definition: Geometry.h:573
GeomState m_geomFactorsState
State of the geometric factors.
Definition: Geometry.h:189
virtual PointGeomSharedPtr v_GetVertex(int i) const =0
GeomFactorsSharedPtr GetMetricInfo()
Get the geometric factors for this object.
Definition: Geometry.h:303
StdRegions::StdExpansionSharedPtr m_xmap
mapping containing isoparametric transformation.
Definition: Geometry.h:191
virtual Geometry1DSharedPtr v_GetEdge(int i) const
Returns edge i of this object.
Definition: Geometry.cpp:161
virtual StdRegions::Orientation v_GetEorient(const int i) const
Returns the orientation of edge i with respect to the ordering of edges in the standard element.
Definition: Geometry.cpp:193
StdRegions::StdExpansionSharedPtr GetXmap() const
Return the mapping object Geometry::m_xmap that represents the coordinate transformation from standar...
Definition: Geometry.h:426
void GenGeomFactors()
Handles generation of geometry factors.
Definition: Geometry.h:650
static GeomFactorsVector m_regGeomFactorsManager
Definition: Geometry.h:182
int GetNumEdges() const
Get the number of edges of this object.
Definition: Geometry.h:398
std::array< NekDouble, 6 > GetBoundingBox()
Generates the bounding box for the element.
Definition: Geometry.cpp:383
int GetNumVerts() const
Get the number of vertices of this object.
Definition: Geometry.h:390
Geometry2DSharedPtr GetFace(int i) const
Returns face i of this object.
Definition: Geometry.h:364
int GetVertexFaceMap(int i, int j) const
Returns the standard element face IDs that are connected to a given vertex.
Definition: Geometry.h:595
bool ContainsPoint(const Array< OneD, const NekDouble > &gloCoord, NekDouble tol=0.0)
Determine whether an element contains a particular Cartesian coordinate .
Definition: Geometry.h:458
int GetEdgeFaceMap(int i, int j) const
Returns the standard element edge IDs that are connected to a given face.
Definition: Geometry.h:616
bool MinMaxCheck(const Array< OneD, const NekDouble > &gloCoord)
Check if given global coord is within the BoundingBox of the element.
Definition: Geometry.cpp:462
virtual int v_GetVertexEdgeMap(int i, int j) const
Returns the standard element edge IDs that are connected to a given vertex.
Definition: Geometry.cpp:278
GeomFactorsSharedPtr GetGeomFactors()
Get the geometric factors for this object, generating them if required.
Definition: Geometry.h:294
static GeomFactorsSharedPtr ValidateRegGeomFactor(GeomFactorsSharedPtr geomFactor)
Check to see if a geometric factor has already been created that contains the same regular informatio...
Definition: Geometry.cpp:90
StdRegions::Orientation GetEorient(const int i) const
Returns the orientation of edge i with respect to the ordering of edges in the standard element.
Definition: Geometry.h:373
virtual int v_GetNumVerts() const
Get the number of vertices of this object.
Definition: Geometry.cpp:183
GeomFactorsSharedPtr m_geomFactors
Geometric factors.
Definition: Geometry.h:187
int m_coordim
Coordinate dimension of this geometry object.
Definition: Geometry.h:185
int GetTid(int i) const
Get the ID of trace i of this object.
Definition: Geometry.h:338
virtual int v_GetDir(const int faceidx, const int facedir) const
Returns the element coordinate direction corresponding to a given face coordinate direction.
Definition: Geometry.cpp:311
void SetCoordim(int coordim)
Sets the coordinate dimension of this object (i.e. the dimension of the space in which this object is...
Definition: Geometry.h:285
virtual int v_GetShapeDim() const
Get the object's shape dimension.
Definition: Geometry.cpp:231
void Reset(CurveMap &curvedEdges, CurveMap &curvedFaces)
Reset this geometry object: unset the current state, zero Geometry::m_coeffs and remove allocated Geo...
Definition: Geometry.h:635
virtual int v_GetVertexFaceMap(int i, int j) const
Returns the standard element face IDs that are connected to a given vertex.
Definition: Geometry.cpp:289
virtual bool v_ContainsPoint(const Array< OneD, const NekDouble > &gloCoord, Array< OneD, NekDouble > &locCoord, NekDouble tol, NekDouble &dist)
Definition: Geometry.cpp:251
int GetEid(int i) const
Get the ID of edge i of this object.
Definition: Geometry.cpp:145
StdRegions::Orientation GetForient(const int i) const
Returns the orientation of face i with respect to the ordering of faces in the standard element.
Definition: Geometry.h:382
bool GlobalIdEquality(const std::shared_ptr< Geometry > &lhs, const std::shared_ptr< Geometry > &rhs)
Definition: Geometry.cpp:128
std::shared_ptr< GeomFactors > GeomFactorsSharedPtr
Pointer to a GeomFactors object.
Definition: GeomFactors.h:62
std::shared_ptr< GeometryVector > GeometryVectorSharedPtr
Definition: Geometry.h:57
bool SortByGlobalId(const std::shared_ptr< Geometry > &lhs, const std::shared_ptr< Geometry > &rhs)
Less than operator to sort Geometry objects by global id when sorting STL containers.
Definition: Geometry.cpp:122
std::unordered_set< GeometrySharedPtr > GeometrySet
Definition: Geometry.h:56
std::shared_ptr< Curve > CurveSharedPtr
Definition: Curve.hpp:61
std::unordered_map< int, CurveSharedPtr > CurveMap
Definition: Curve.hpp:62
GeomType
Indicates the type of element geometry.
static CurveMap NullCurveMap
Definition: Geometry.h:70
std::shared_ptr< PointGeom > PointGeomSharedPtr
Definition: Geometry.h:59
std::shared_ptr< Geometry2D > Geometry2DSharedPtr
Definition: Geometry.h:65
GeomState
Indicates if the geometric information for an element has been populated.
std::shared_ptr< Geometry > GeometrySharedPtr
Definition: Geometry.h:53
std::shared_ptr< Geometry1D > Geometry1DSharedPtr
Definition: Geometry.h:63
std::vector< GeomFactorsSharedPtr > GeomFactorsVector
A vector of GeomFactor pointers.
Definition: GeomFactors.h:64
std::vector< GeometrySharedPtr > GeometryVector
Definition: Geometry.h:55
std::shared_ptr< StdExpansion > StdExpansionSharedPtr
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:1
std::size_t hash_range(Iter first, Iter last)
Definition: HashUtils.hpp:69
double NekDouble
Unary function that constructs a hash of a Geometry object, based on the vertex IDs.
Definition: Geometry.h:253
std::size_t operator()(GeometrySharedPtr const &p) const
Definition: Geometry.h:254