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 <array>
45 #include <unordered_map>
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, const std::shared_ptr<Geometry> &rhs);
76 
78  const std::shared_ptr<Geometry> &lhs, const std::shared_ptr<Geometry> &rhs);
79 
80 /// Base class for shape geometry information
81 class Geometry
82 {
83 public:
85  SPATIAL_DOMAINS_EXPORT Geometry(int coordim);
86 
88 
89  //---------------------------------------
90  // Helper functions
91  //---------------------------------------
92 
93  SPATIAL_DOMAINS_EXPORT inline int GetCoordim() const;
94  SPATIAL_DOMAINS_EXPORT inline void SetCoordim(int coordim);
95 
101 
102  //---------------------------------------
103  // Set and get ID
104  //---------------------------------------
105  SPATIAL_DOMAINS_EXPORT inline int GetGlobalID(void) const;
106  SPATIAL_DOMAINS_EXPORT inline void SetGlobalID(int globalid);
107 
108  //---------------------------------------
109  // Vertex, edge and face access
110  //---------------------------------------
111  SPATIAL_DOMAINS_EXPORT int GetVid(int i) const;
112  SPATIAL_DOMAINS_EXPORT int GetEid(int i) const;
113  SPATIAL_DOMAINS_EXPORT int GetFid(int i) const;
114  SPATIAL_DOMAINS_EXPORT inline int GetTid(int i) const;
119  const int i) const;
121  const int i) const;
122  SPATIAL_DOMAINS_EXPORT inline int GetNumVerts() const;
123  SPATIAL_DOMAINS_EXPORT inline int GetNumEdges() const;
124  SPATIAL_DOMAINS_EXPORT inline int GetNumFaces() const;
125  SPATIAL_DOMAINS_EXPORT inline int GetShapeDim() const;
126 
127  //---------------------------------------
128  // \chi mapping access
129  //---------------------------------------
131  const;
133  const int i) const;
134  SPATIAL_DOMAINS_EXPORT inline void FillGeom();
135 
136  //---------------------------------------
137  // Point lookups
138  //---------------------------------------
139  SPATIAL_DOMAINS_EXPORT std::array<NekDouble, 6> GetBoundingBox();
141 
143  const Array<OneD, const NekDouble> &gloCoord, NekDouble tol = 0.0);
145  const Array<OneD, const NekDouble> &gloCoord,
146  Array<OneD, NekDouble> &locCoord, NekDouble tol);
148  const Array<OneD, const NekDouble> &gloCoord,
149  Array<OneD, NekDouble> &locCoord, NekDouble tol, NekDouble &dist);
151  const Array<OneD, const NekDouble> &coords,
152  Array<OneD, NekDouble> &Lcoords);
154  const int i, const Array<OneD, const NekDouble> &Lcoord);
156  const Array<OneD, const NekDouble> &gloCoord);
158  Array<OneD, NekDouble> &locCoord,
159  NekDouble tol = std::numeric_limits<NekDouble>::epsilon());
162 
163  //---------------------------------------
164  // Misc. helper functions
165  //---------------------------------------
166  SPATIAL_DOMAINS_EXPORT inline int GetVertexEdgeMap(int i, int j) const;
167  SPATIAL_DOMAINS_EXPORT inline int GetVertexFaceMap(int i, int j) const;
168  SPATIAL_DOMAINS_EXPORT inline int GetEdgeFaceMap(int i, int j) const;
170  int j) const;
171  SPATIAL_DOMAINS_EXPORT inline int GetDir(const int i,
172  const int j = 0) const;
173 
174  SPATIAL_DOMAINS_EXPORT inline void Reset(CurveMap &curvedEdges,
175  CurveMap &curvedFaces);
176  SPATIAL_DOMAINS_EXPORT inline void Setup();
177 
178 protected:
180  GeomFactorsSharedPtr geomFactor);
182 
183  /// Coordinate dimension of this geometry object.
185  /// Geometric factors.
187  /// State of the geometric factors
189  /// \f$\chi\f$ mapping containing isoparametric transformation.
191  /// Enumeration to dictate whether coefficients are filled.
193  /// Wether or not the setup routines have been run
195  /// Type of geometry.
197  /// Type of shape.
199  /// Global ID
201  /// Array containing expansion coefficients of @p m_xmap
203  /// Array containing bounding box
205 
206  /// Handles generation of geometry factors.
207  void GenGeomFactors();
208 
209  //---------------------------------------
210  // Helper functions
211  //---------------------------------------
212  virtual PointGeomSharedPtr v_GetVertex(int i) const = 0;
213  virtual Geometry1DSharedPtr v_GetEdge(int i) const;
214  virtual Geometry2DSharedPtr v_GetFace(int i) const;
215  virtual StdRegions::Orientation v_GetEorient(const int i) const;
216  virtual StdRegions::Orientation v_GetForient(const int i) const;
217  virtual int v_GetNumVerts() const;
218  virtual int v_GetNumEdges() const;
219  virtual int v_GetNumFaces() const;
220  virtual int v_GetShapeDim() const;
221 
223  virtual void v_FillGeom();
224 
225  virtual bool v_ContainsPoint(const Array<OneD, const NekDouble> &gloCoord,
226  Array<OneD, NekDouble> &locCoord,
227  NekDouble tol, NekDouble &dist);
228 
229  virtual NekDouble v_GetCoord(const int i,
230  const Array<OneD, const NekDouble> &Lcoord);
232  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_GetEdgeNormalToFaceVert(const int i, const int j) const;
240  virtual int v_GetDir(const int faceidx, const int facedir) const;
241 
242  virtual void v_Reset(CurveMap &curvedEdges, CurveMap &curvedFaces);
243  virtual void v_Setup();
244  virtual void v_GenGeomFactors() = 0;
245 
246  inline void SetUpCoeffs(const int nCoeffs);
247 }; // class Geometry
248 
249 /**
250  * @brief Unary function that constructs a hash of a Geometry object, based on
251  * the vertex IDs.
252  */
254 {
255  std::size_t operator()(GeometrySharedPtr const &p) const
256  {
257  int i;
258  size_t seed = 0;
259  int nVert = p->GetNumVerts();
260  std::vector<unsigned int> ids(nVert);
261 
262  for (i = 0; i < nVert; ++i)
263  {
264  ids[i] = p->GetVid(i);
265  }
266  std::sort(ids.begin(), ids.end());
267  hash_range(seed, ids.begin(), ids.end());
268 
269  return seed;
270  }
271 };
272 
273 /**
274  * @brief Return the coordinate dimension of this object (i.e. the dimension of
275  * the space in which this object is embedded).
276  */
277 inline int Geometry::GetCoordim() const
278 {
279  return m_coordim;
280 }
281 
282 /**
283  * @brief Sets the coordinate dimension of this object (i.e. the dimension of
284  * the space in which this object is embedded).
285  */
286 inline void Geometry::SetCoordim(int dim)
287 {
288  m_coordim = dim;
289 }
290 
291 /**
292  * @brief Get the geometric factors for this object, generating them if
293  * required.
294  */
296 {
297  GenGeomFactors();
299 }
300 
301 /**
302  * @brief Get the geometric factors for this object.
303  */
305 {
306  return m_geomFactors;
307 }
308 
309 /**
310  * @brief Get the geometric shape type of this object.
311  */
313 {
314  return m_shapeType;
315 }
316 
317 /**
318  * @brief Get the ID of this object.
319  */
320 inline int Geometry::GetGlobalID(void) const
321 {
322  return m_globalID;
323 }
324 
325 /**
326  * @brief Set the ID of this object.
327  */
328 inline void Geometry::SetGlobalID(int globalid)
329 {
330  m_globalID = globalid;
331 }
332 
333 /**
334  * @brief Get the ID of trace @p i of this object.
335  *
336  * The trace element is the facet one dimension lower than the object; for
337  * example, a quadrilateral has four trace segments forming its boundary.
338  */
339 inline int Geometry::GetTid(int i) const
340 {
341  const int nDim = GetShapeDim();
342  return nDim == 1 ? GetVid(i)
343  : nDim == 2 ? GetEid(i)
344  : nDim == 3 ? GetFid(i)
345  : 0;
346 }
347 
348 /**
349  * @brief Returns vertex @p i of this object.
350  */
352 {
353  return v_GetVertex(i);
354 }
355 
356 /**
357  * @brief Returns edge @p i of this object.
358  */
360 {
361  return v_GetEdge(i);
362 }
363 
364 /**
365  * @brief Returns face @p i of this object.
366  */
368 {
369  return v_GetFace(i);
370 }
371 
372 /**
373  * @brief Returns the orientation of edge @p i with respect to the ordering of
374  * edges in the standard element.
375  */
377 {
378  return v_GetEorient(i);
379 }
380 
381 /**
382  * @brief Returns the orientation of face @p i with respect to the ordering of
383  * faces in the standard element.
384  */
386 {
387  return v_GetForient(i);
388 }
389 
390 /**
391  * @brief Get the number of vertices of this object.
392  */
393 inline int Geometry::GetNumVerts() const
394 {
395  return v_GetNumVerts();
396 }
397 
398 /**
399  * @brief Get the number of edges of this object.
400  */
401 inline int Geometry::GetNumEdges() const
402 {
403  return v_GetNumEdges();
404 }
405 
406 /**
407  * @brief Get the number of faces of this object.
408  */
409 inline int Geometry::GetNumFaces() const
410 {
411  return v_GetNumFaces();
412 }
413 
414 /**
415  * @brief Get the object's shape dimension.
416  *
417  * For example, a segment is one dimensional and quadrilateral is two
418  * dimensional.
419  */
420 inline int Geometry::GetShapeDim() const
421 {
422  return v_GetShapeDim();
423 }
424 
425 /**
426  * @brief Return the mapping object Geometry::m_xmap that represents the
427  * coordinate transformation from standard element to physical element.
428  */
430 {
431  return v_GetXmap();
432 }
433 
434 /**
435  * @brief Return the coefficients of the transformation Geometry::m_xmap in
436  * coordinate direction @p i.
437  */
439  const int i) const
440 {
441  return m_coeffs[i];
442 }
443 
444 /**
445  * @brief Populate the coordinate mapping Geometry::m_coeffs information from
446  * any children geometry elements.
447  *
448  * @see v_FillGeom()
449  */
450 inline void Geometry::FillGeom()
451 {
452  v_FillGeom();
453 }
454 
455 /**
456  * @brief Determine whether an element contains a particular Cartesian
457  * coordinate \f$(x,y,z)\f$.
458  *
459  * @see Geometry::ContainsPoint
460  */
462  const Array<OneD, const NekDouble> &gloCoord, NekDouble tol)
463 {
464  Array<OneD, NekDouble> locCoord(GetCoordim(), 0.0);
465  NekDouble dist;
466  return v_ContainsPoint(gloCoord, locCoord, tol, dist);
467 }
468 
469 /**
470  * @brief Determine whether an element contains a particular Cartesian
471  * coordinate \f$(x,y,z)\f$.
472  *
473  * @see Geometry::ContainsPoint
474  */
476  const Array<OneD, const NekDouble> &gloCoord,
477  Array<OneD, NekDouble> &locCoord, NekDouble tol)
478 {
479  NekDouble dist;
480  return v_ContainsPoint(gloCoord, locCoord, tol, dist);
481 }
482 
483 /**
484  * @brief Determine whether an element contains a particular Cartesian
485  * coordinate \f$\vec{x} = (x,y,z)\f$.
486  *
487  * For curvilinear and non-affine elements (i.e. where the Jacobian varies as a
488  * function of the standard element coordinates), this is a non-linear
489  * optimisation problem that requires the use of a Newton iteration. Note
490  * therefore that this can be an expensive operation.
491  *
492  * The parameter @p tol which is by default 0, can be used to expand the
493  * coordinate range of the standard element from \f$[-1,1]^d\f$ to
494  * \f$[-1-\epsilon,1+\epsilon\f$ to handle challenging edge cases. The function
495  * also returns the local coordinates corresponding to @p gloCoord that can be
496  * used to speed up subsequent searches.
497  *
498  * @param gloCoord The coordinate \f$ (x,y,z) \f$.
499  * @param locCoord On exit, this is the local coordinate \f$\vec{\xi}\f$ such
500  * that \f$\chi(\vec{\xi}) = \vec{x}\f$.
501  * @param tol The tolerance used to dictate the bounding box of the
502  * standard coordinates \f$\vec{\xi}\f$.
503  * @param dist On exit, returns the minimum distance between @p gloCoord
504  * and the quadrature points inside the element.
505  *
506  * @return `true` if the coordinate @p gloCoord is contained in the element;
507  * false otherwise.
508  *
509  * @see Geometry::GetLocCoords.
510  */
512  const Array<OneD, const NekDouble> &gloCoord,
513  Array<OneD, NekDouble> &locCoord, NekDouble tol, NekDouble &dist)
514 {
515  return v_ContainsPoint(gloCoord, locCoord, tol, dist);
516 }
517 
518 /**
519  * @brief Determine the local collapsed coordinates that correspond to a
520  * given Cartesian coordinate for this geometry object.
521  *
522  * For curvilinear and non-affine elements (i.e. where the Jacobian varies as a
523  * function of the standard element coordinates), this is a non-linear
524  * optimisation problem that requires the use of a Newton iteration. Note
525  * therefore that this can be an expensive operation.
526  *
527  * Note that, clearly, the provided Cartesian coordinate lie outside the
528  * element. The function therefore returns the minimum distance from some
529  * position in the element to . @p Lcoords will also be constrained to fit
530  * within the range \f$[-1,1]^d\f$ where \f$ d \f$ is the dimension of the
531  * element.
532  *
533  * @param coords Input Cartesian global coordinates
534  * @param Lcoords Corresponding local coordinates
535  *
536  * @return Distance between obtained coordinates and provided ones.
537  */
540 {
541  return v_GetLocCoords(coords, Lcoords);
542 }
543 
544 /**
545  * @brief Given local collapsed coordinate @p Lcoord, return the value of
546  * physical coordinate in direction @p i.
547  */
548 inline NekDouble Geometry::GetCoord(const int i,
549  const Array<OneD, const NekDouble> &Lcoord)
550 {
551  return v_GetCoord(i, Lcoord);
552 }
553 
556 {
557  return v_FindDistance(xs, xi);
558 }
559 
560 /**
561  * @brief Returns the standard element edge IDs that are connected to a given
562  * vertex.
563  *
564  * For example, on a prism, vertex 0 is connnected to edges 0, 3, and 4;
565  * `GetVertexEdgeMap(0,j)` would therefore return the values 0, 1 and 4
566  * respectively. We assume that @p j runs between 0 and 2 inclusive, which is
567  * true for every 3D element asides from the pyramid.
568  *
569  * This function is used in the construction of the low-energy preconditioner.
570  *
571  * @param i The vertex to query connectivity for.
572  * @param j The local edge index between 0 and 2 connected to this element.
573  *
574  * @todo Expand to work with pyramid elements.
575  * @see MultiRegions::PreconditionerLowEnergy
576  */
577 inline int Geometry::GetVertexEdgeMap(int i, int j) const
578 {
579  return v_GetVertexEdgeMap(i, j);
580 }
581 
582 /**
583  * @brief Returns the standard element face IDs that are connected to a given
584  * vertex.
585  *
586  * For example, on a hexahedron, vertex 0 is connnected to faces 0, 1, and 4;
587  * `GetVertexFaceMap(0,j)` would therefore return the values 0, 1 and 4
588  * respectively. We assume that @p j runs between 0 and 2 inclusive, which is
589  * true for every 3D element asides from the pyramid.
590  *
591  * This is used in the construction of the low-energy preconditioner.
592  *
593  * @param i The vertex to query connectivity for.
594  * @param j The local face index between 0 and 2 connected to this element.
595  *
596  * @todo Expand to work with pyramid elements.
597  * @see MultiRegions::PreconditionerLowEnergy
598  */
599 inline int Geometry::GetVertexFaceMap(int i, int j) const
600 {
601  return v_GetVertexFaceMap(i, j);
602 }
603 
604 /**
605  * @brief Returns the standard element edge IDs that are connected to a given
606  * face.
607  *
608  * For example, on a prism, edge 0 is connnected to faces 0 and 1;
609  * `GetEdgeFaceMap(0,j)` would therefore return the values 0 and 1
610  * respectively. We assume that @p j runs between 0 and 1 inclusive, since every
611  * face is connected to precisely two faces for all 3D elements.
612  *
613  * This function is used in the construction of the low-energy preconditioner.
614  *
615  * @param i The edge to query connectivity for.
616  * @param j The local face index between 0 and 1 connected to this element.
617  *
618  * @see MultiRegions::PreconditionerLowEnergy
619  */
620 inline int Geometry::GetEdgeFaceMap(int i, int j) const
621 {
622  return v_GetEdgeFaceMap(i, j);
623 }
624 
625 /**
626  * @brief Returns the standard lement edge IDs that are normal to a given face
627  * vertex.
628  *
629  * For example, on a hexahedron, on face 0 at vertices 0,1,2,3 the
630  * edges normal to that face are 4,5,6,7, ; so
631  * `GetEdgeNormalToFaceVert(0,j)` would therefore return the values 4,
632  * 5, 6 and 7 respectively. We assume that @p j runs between 0 and 3
633  * inclusive on a quadrilateral face and between 0 and 2 inclusive on
634  * a triangular face.
635  *
636  * This is used to help set up a length scale normal to an face
637  *
638  * @param i The face to query for the normal edge
639  * @param j The local vertex index between 0 and nverts on this face
640  *
641  */
642 inline int Geometry::GetEdgeNormalToFaceVert(int i, int j) const
643 {
644  return v_GetEdgeNormalToFaceVert(i, j);
645 }
646 
647 /**
648  * @brief Returns the element coordinate direction corresponding to a given face
649  * coordinate direction
650  */
651 inline int Geometry::GetDir(const int faceidx, const int facedir) const
652 {
653  return v_GetDir(faceidx, facedir);
654 }
655 
656 /**
657  * @brief Reset this geometry object: unset the current state, zero
658  * Geometry::m_coeffs and remove allocated GeomFactors.
659  */
660 inline void Geometry::Reset(CurveMap &curvedEdges, CurveMap &curvedFaces)
661 {
662  v_Reset(curvedEdges, curvedFaces);
663 }
664 inline void Geometry::Setup()
665 {
666  v_Setup();
667 }
668 
669 /**
670  * @brief Generate the geometric factors (i.e. derivatives of \f$\chi\f$) and
671  * related metrics.
672  *
673  * @see SpatialDomains::GeomFactors
674  */
676 {
677  return v_GenGeomFactors();
678 }
679 
680 /**
681  * @brief Initialise the Geometry::m_coeffs array.
682  */
683 inline void Geometry::SetUpCoeffs(const int nCoeffs)
684 {
686 
687  for (int i = 0; i < m_coordim; ++i)
688  {
689  m_coeffs[i] = Array<OneD, NekDouble>(nCoeffs, 0.0);
690  }
691 }
692 
693 } // namespace SpatialDomains
694 } // namespace Nektar
695 
696 #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:82
virtual int v_GetNumEdges() const
Get the number of edges of this object.
Definition: Geometry.cpp:217
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:206
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:548
LibUtilities::ShapeType GetShapeType(void)
Get the geometric shape type of this object.
Definition: Geometry.h:312
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:343
Array< OneD, NekDouble > m_boundingBox
Array containing bounding box.
Definition: Geometry.h:204
GeomFactorsSharedPtr GetRefGeomFactors(const Array< OneD, const LibUtilities::BasisSharedPtr > &tbasis)
bool m_setupState
Wether or not the setup routines have been run.
Definition: Geometry.h:194
PointGeomSharedPtr GetVertex(int i) const
Returns vertex i of this object.
Definition: Geometry.h:351
GeomState m_state
Enumeration to dictate whether coefficients are filled.
Definition: Geometry.h:192
void SetUpCoeffs(const int nCoeffs)
Initialise the Geometry::m_coeffs array.
Definition: Geometry.h:683
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:538
GeomType m_geomType
Type of geometry.
Definition: Geometry.h:196
virtual int v_GetEdgeNormalToFaceVert(const int i, const int j) const
Returns the standard lement edge IDs that are normal to a given face vertex.
Definition: Geometry.cpp:321
Geometry1DSharedPtr GetEdge(int i) const
Returns edge i of this object.
Definition: Geometry.h:359
void SetGlobalID(int globalid)
Set the ID of this object.
Definition: Geometry.h:328
int GetNumFaces() const
Get the number of faces of this object.
Definition: Geometry.h:409
int GetShapeDim() const
Get the object's shape dimension.
Definition: Geometry.h:420
virtual void v_FillGeom()
Populate the coordinate mapping Geometry::m_coeffs information from any children geometry elements.
Definition: Geometry.cpp:367
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:355
virtual Geometry2DSharedPtr v_GetFace(int i) const
Returns face i of this object.
Definition: Geometry.cpp:174
int GetVid(int i) const
Get the ID of vertex i of this object.
Definition: Geometry.cpp:139
virtual int v_GetNumFaces() const
Get the number of faces of this object.
Definition: Geometry.cpp:225
virtual StdRegions::StdExpansionSharedPtr v_GetXmap() const
Return the mapping object Geometry::m_xmap that represents the coordinate transformation from standar...
Definition: Geometry.cpp:243
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:376
bool ClampLocCoords(Array< OneD, NekDouble > &locCoord, NekDouble tol=std::numeric_limits< NekDouble >::epsilon())
Clamp local coords to be within standard regions [-1, 1]^dim.
Definition: Geometry.cpp:510
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:438
int GetGlobalID(void) const
Get the ID of this object.
Definition: Geometry.h:320
int GetCoordim() const
Return the coordinate dimension of this object (i.e. the dimension of the space in which this object ...
Definition: Geometry.h:277
int GetFid(int i) const
Get the ID of face i of this object.
Definition: Geometry.cpp:155
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:310
LibUtilities::ShapeType m_shapeType
Type of shape.
Definition: Geometry.h:198
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:651
void FillGeom()
Populate the coordinate mapping Geometry::m_coeffs information from any children geometry elements.
Definition: Geometry.h:450
Array< OneD, Array< OneD, NekDouble > > m_coeffs
Array containing expansion coefficients of m_xmap.
Definition: Geometry.h:202
virtual ~Geometry()
Default destructor.
Definition: Geometry.cpp:74
int GetVertexEdgeMap(int i, int j) const
Returns the standard element edge IDs that are connected to a given vertex.
Definition: Geometry.h:577
GeomState m_geomFactorsState
State of the geometric factors.
Definition: Geometry.h:188
virtual PointGeomSharedPtr v_GetVertex(int i) const =0
GeomFactorsSharedPtr GetMetricInfo()
Get the geometric factors for this object.
Definition: Geometry.h:304
StdRegions::StdExpansionSharedPtr m_xmap
mapping containing isoparametric transformation.
Definition: Geometry.h:190
virtual Geometry1DSharedPtr v_GetEdge(int i) const
Returns edge i of this object.
Definition: Geometry.cpp:163
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:195
StdRegions::StdExpansionSharedPtr GetXmap() const
Return the mapping object Geometry::m_xmap that represents the coordinate transformation from standar...
Definition: Geometry.h:429
NekDouble FindDistance(const Array< OneD, const NekDouble > &xs, Array< OneD, NekDouble > &xi)
Definition: Geometry.h:554
void GenGeomFactors()
Handles generation of geometry factors.
Definition: Geometry.h:675
static GeomFactorsVector m_regGeomFactorsManager
Definition: Geometry.h:181
int GetNumEdges() const
Get the number of edges of this object.
Definition: Geometry.h:401
std::array< NekDouble, 6 > GetBoundingBox()
Generates the bounding box for the element.
Definition: Geometry.cpp:403
int GetNumVerts() const
Get the number of vertices of this object.
Definition: Geometry.h:393
Geometry2DSharedPtr GetFace(int i) const
Returns face i of this object.
Definition: Geometry.h:367
int GetVertexFaceMap(int i, int j) const
Returns the standard element face IDs that are connected to a given vertex.
Definition: Geometry.h:599
bool ContainsPoint(const Array< OneD, const NekDouble > &gloCoord, NekDouble tol=0.0)
Determine whether an element contains a particular Cartesian coordinate .
Definition: Geometry.h:461
int GetEdgeFaceMap(int i, int j) const
Returns the standard element edge IDs that are connected to a given face.
Definition: Geometry.h:620
bool MinMaxCheck(const Array< OneD, const NekDouble > &gloCoord)
Check if given global coord is within the BoundingBox of the element.
Definition: Geometry.cpp:487
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:288
int GetEdgeNormalToFaceVert(int i, int j) const
Returns the standard lement edge IDs that are normal to a given face vertex.
Definition: Geometry.h:642
virtual NekDouble v_FindDistance(const Array< OneD, const NekDouble > &xs, Array< OneD, NekDouble > &xi)
Definition: Geometry.cpp:276
GeomFactorsSharedPtr GetGeomFactors()
Get the geometric factors for this object, generating them if required.
Definition: Geometry.h:295
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:92
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:376
virtual int v_GetNumVerts() const
Get the number of vertices of this object.
Definition: Geometry.cpp:185
GeomFactorsSharedPtr m_geomFactors
Geometric factors.
Definition: Geometry.h:186
int m_coordim
Coordinate dimension of this geometry object.
Definition: Geometry.h:184
int GetTid(int i) const
Get the ID of trace i of this object.
Definition: Geometry.h:339
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:332
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:286
virtual int v_GetShapeDim() const
Get the object's shape dimension.
Definition: Geometry.cpp:233
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:660
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:299
virtual bool v_ContainsPoint(const Array< OneD, const NekDouble > &gloCoord, Array< OneD, NekDouble > &locCoord, NekDouble tol, NekDouble &dist)
Definition: Geometry.cpp:253
int GetEid(int i) const
Get the ID of edge i of this object.
Definition: Geometry.cpp:147
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:385
std::vector< GeomFactorsSharedPtr > GeomFactorsVector
A vector of GeomFactor pointers.
Definition: GeomFactors.h:64
bool GlobalIdEquality(const std::shared_ptr< Geometry > &lhs, const std::shared_ptr< Geometry > &rhs)
Definition: Geometry.cpp:130
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:124
std::unordered_set< GeometrySharedPtr > GeometrySet
Definition: Geometry.h:56
std::shared_ptr< Curve > CurveSharedPtr
Definition: Curve.hpp:60
std::unordered_map< int, CurveSharedPtr > CurveMap
Definition: Curve.hpp:61
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< 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:2
std::size_t hash_range(Iter first, Iter last)
Definition: HashUtils.hpp:68
double NekDouble
Unary function that constructs a hash of a Geometry object, based on the vertex IDs.
Definition: Geometry.h:254
std::size_t operator()(GeometrySharedPtr const &p) const
Definition: Geometry.h:255