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