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