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