Nektar++
Public Member Functions | Protected Member Functions | Static Protected Member Functions | Protected Attributes | Static Protected Attributes | List of all members
Nektar::SpatialDomains::Geometry Class Referenceabstract

Base class for shape geometry information. More...

#include <Geometry.h>

Inheritance diagram for Nektar::SpatialDomains::Geometry:
[legend]

Public Member Functions

 Geometry ()
 Default constructor. More...
 
 Geometry (int coordim)
 Constructor when supplied a coordinate dimension. More...
 
virtual ~Geometry ()=default
 
int GetCoordim () const
 Return the coordinate dimension of this object (i.e. the dimension of the space in which this object is embedded). More...
 
void SetCoordim (int coordim)
 Sets the coordinate dimension of this object (i.e. the dimension of the space in which this object is embedded). More...
 
GeomFactorsSharedPtr GetGeomFactors ()
 Get the geometric factors for this object, generating them if required. More...
 
GeomFactorsSharedPtr GetRefGeomFactors (const Array< OneD, const LibUtilities::BasisSharedPtr > &tbasis)
 
GeomFactorsSharedPtr GetMetricInfo ()
 Get the geometric factors for this object. More...
 
LibUtilities::ShapeType GetShapeType (void)
 Get the geometric shape type of this object. More...
 
int GetGlobalID (void) const
 Get the ID of this object. More...
 
void SetGlobalID (int globalid)
 Set the ID of this object. More...
 
int GetVid (int i) const
 Get the ID of vertex i of this object. More...
 
int GetEid (int i) const
 Get the ID of edge i of this object. More...
 
int GetFid (int i) const
 Get the ID of face i of this object. More...
 
int GetTid (int i) const
 Get the ID of trace i of this object. More...
 
PointGeomSharedPtr GetVertex (int i) const
 Returns vertex i of this object. More...
 
Geometry1DSharedPtr GetEdge (int i) const
 Returns edge i of this object. More...
 
Geometry2DSharedPtr GetFace (int i) const
 Returns face i of this object. More...
 
StdRegions::Orientation GetEorient (const int i) const
 Returns the orientation of edge i with respect to the ordering of edges in the standard element. More...
 
StdRegions::Orientation GetForient (const int i) const
 Returns the orientation of face i with respect to the ordering of faces in the standard element. More...
 
int GetNumVerts () const
 Get the number of vertices of this object. More...
 
int GetNumEdges () const
 Get the number of edges of this object. More...
 
int GetNumFaces () const
 Get the number of faces of this object. More...
 
int GetShapeDim () const
 Get the object's shape dimension. More...
 
StdRegions::StdExpansionSharedPtr GetXmap () const
 Return the mapping object Geometry::m_xmap that represents the coordinate transformation from standard element to physical element. More...
 
const Array< OneD, const NekDouble > & GetCoeffs (const int i) const
 Return the coefficients of the transformation Geometry::m_xmap in coordinate direction i. More...
 
void FillGeom ()
 Populate the coordinate mapping Geometry::m_coeffs information from any children geometry elements. More...
 
std::array< NekDouble, 6 > GetBoundingBox ()
 Generates the bounding box for the element. More...
 
void ClearBoundingBox ()
 
bool ContainsPoint (const Array< OneD, const NekDouble > &gloCoord, NekDouble tol=0.0)
 Determine whether an element contains a particular Cartesian coordinate \((x,y,z)\). More...
 
bool ContainsPoint (const Array< OneD, const NekDouble > &gloCoord, Array< OneD, NekDouble > &locCoord, NekDouble tol)
 Determine whether an element contains a particular Cartesian coordinate \((x,y,z)\). More...
 
bool ContainsPoint (const Array< OneD, const NekDouble > &gloCoord, Array< OneD, NekDouble > &locCoord, NekDouble tol, NekDouble &dist)
 Determine whether an element contains a particular Cartesian coordinate \(\vec{x} = (x,y,z)\). More...
 
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 geometry object. More...
 
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. More...
 
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 elements, this check is unnecessary. More...
 
bool MinMaxCheck (const Array< OneD, const NekDouble > &gloCoord)
 Check if given global coord is within the BoundingBox of the element. More...
 
bool ClampLocCoords (Array< OneD, NekDouble > &locCoord, NekDouble tol=std::numeric_limits< NekDouble >::epsilon())
 Clamp local coords to be within standard regions [-1, 1]^dim. More...
 
NekDouble FindDistance (const Array< OneD, const NekDouble > &xs, Array< OneD, NekDouble > &xi)
 
int GetVertexEdgeMap (int i, int j) const
 Returns the standard element edge IDs that are connected to a given vertex. More...
 
int GetVertexFaceMap (int i, int j) const
 Returns the standard element face IDs that are connected to a given vertex. More...
 
int GetEdgeFaceMap (int i, int j) const
 Returns the standard element edge IDs that are connected to a given face. More...
 
int GetEdgeNormalToFaceVert (int i, int j) const
 Returns the standard lement edge IDs that are normal to a given face vertex. More...
 
int GetDir (const int i, const int j=0) const
 Returns the element coordinate direction corresponding to a given face coordinate direction. More...
 
void Reset (CurveMap &curvedEdges, CurveMap &curvedFaces)
 Reset this geometry object: unset the current state, zero Geometry::m_coeffs and remove allocated GeomFactors. More...
 
void ResetNonRecursive (CurveMap &curvedEdges, CurveMap &curvedFaces)
 Reset this geometry object non-recursively: unset the current state, zero Geometry::m_coeffs and remove allocated GeomFactors. More...
 
void Setup ()
 
void GenGeomFactors ()
 Handles generation of geometry factors. More...
 

Protected Member Functions

virtual PointGeomSharedPtr v_GetVertex (int i) const =0
 
virtual Geometry1DSharedPtr v_GetEdge (int i) const
 Returns edge i of this object. More...
 
virtual Geometry2DSharedPtr v_GetFace (int i) const
 Returns face i of this object. More...
 
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. More...
 
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. More...
 
virtual int v_GetNumVerts () const
 Get the number of vertices of this object. More...
 
virtual int v_GetNumEdges () const
 Get the number of edges of this object. More...
 
virtual int v_GetNumFaces () const
 Get the number of faces of this object. More...
 
virtual int v_GetShapeDim () const
 Get the object's shape dimension. More...
 
virtual StdRegions::StdExpansionSharedPtr v_GetXmap () const
 Return the mapping object Geometry::m_xmap that represents the coordinate transformation from standard element to physical element. More...
 
virtual void v_FillGeom ()
 Populate the coordinate mapping Geometry::m_coeffs information from any children geometry elements. More...
 
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 \(\vec{x} = (x,y,z)\). More...
 
virtual int v_AllLeftCheck (const Array< OneD, const NekDouble > &gloCoord)
 
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. More...
 
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 geometry object. More...
 
virtual NekDouble v_FindDistance (const Array< OneD, const NekDouble > &xs, Array< OneD, NekDouble > &xi)
 
virtual int v_GetVertexEdgeMap (int i, int j) const
 Returns the standard element edge IDs that are connected to a given vertex. More...
 
virtual int v_GetVertexFaceMap (int i, int j) const
 Returns the standard element face IDs that are connected to a given vertex. More...
 
virtual int v_GetEdgeFaceMap (int i, int j) const
 Returns the standard element edge IDs that are connected to a given face. More...
 
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. More...
 
virtual int v_GetDir (const int faceidx, const int facedir) const
 Returns the element coordinate direction corresponding to a given face coordinate direction. More...
 
virtual void v_Reset (CurveMap &curvedEdges, CurveMap &curvedFaces)
 Reset this geometry object: unset the current state, zero Geometry::m_coeffs and remove allocated GeomFactors. More...
 
virtual void v_Setup ()
 
virtual void v_GenGeomFactors ()=0
 
void SetUpCoeffs (const int nCoeffs)
 Initialise the Geometry::m_coeffs array. More...
 
virtual void v_CalculateInverseIsoParam ()
 

Static Protected Member Functions

static GeomFactorsSharedPtr ValidateRegGeomFactor (GeomFactorsSharedPtr geomFactor)
 Check to see if a geometric factor has already been created that contains the same regular information. More...
 

Protected Attributes

int m_coordim
 Coordinate dimension of this geometry object. More...
 
GeomFactorsSharedPtr m_geomFactors
 Geometric factors. More...
 
GeomState m_geomFactorsState
 State of the geometric factors. More...
 
StdRegions::StdExpansionSharedPtr m_xmap
 \(\chi\) mapping containing isoparametric transformation. More...
 
GeomState m_state
 Enumeration to dictate whether coefficients are filled. More...
 
bool m_setupState
 Wether or not the setup routines have been run. More...
 
GeomType m_geomType
 Type of geometry. More...
 
LibUtilities::ShapeType m_shapeType
 Type of shape. More...
 
int m_globalID
 Global ID. More...
 
Array< OneD, Array< OneD, NekDouble > > m_coeffs
 Array containing expansion coefficients of m_xmap. More...
 
Array< OneD, NekDoublem_boundingBox
 Array containing bounding box. More...
 
Array< OneD, Array< OneD, NekDouble > > m_isoParameter
 
Array< OneD, Array< OneD, NekDouble > > m_invIsoParam
 
int m_straightEdge
 

Static Protected Attributes

static GeomFactorsVector m_regGeomFactorsManager
 

Detailed Description

Base class for shape geometry information.

Definition at line 78 of file Geometry.h.

Constructor & Destructor Documentation

◆ Geometry() [1/2]

Nektar::SpatialDomains::Geometry::Geometry ( )

Default constructor.

Definition at line 50 of file Geometry.cpp.

54{
55}
bool m_setupState
Wether or not the setup routines have been run.
Definition: Geometry.h:198
GeomState m_state
Enumeration to dictate whether coefficients are filled.
Definition: Geometry.h:196
LibUtilities::ShapeType m_shapeType
Type of shape.
Definition: Geometry.h:202
GeomState m_geomFactorsState
State of the geometric factors.
Definition: Geometry.h:192
int m_coordim
Coordinate dimension of this geometry object.
Definition: Geometry.h:188
@ eNotFilled
Geometric information has not been generated.

◆ Geometry() [2/2]

Nektar::SpatialDomains::Geometry::Geometry ( int  coordim)

Constructor when supplied a coordinate dimension.

Definition at line 60 of file Geometry.cpp.

◆ ~Geometry()

virtual Nektar::SpatialDomains::Geometry::~Geometry ( )
virtualdefault

Member Function Documentation

◆ ClampLocCoords()

bool Nektar::SpatialDomains::Geometry::ClampLocCoords ( Array< OneD, NekDouble > &  locCoord,
NekDouble  tol = std::numeric_limits<NekDouble>::epsilon() 
)

Clamp local coords to be within standard regions [-1, 1]^dim.

Parameters
LcoordsCorresponding local coordinates

Definition at line 530 of file Geometry.cpp.

531{
532 // Validation checks
533 ASSERTL1(locCoord.size() >= GetShapeDim(),
534 "Expects local coordinates to be same or "
535 "larger than shape dimension.");
536
537 // If out of range clamp locCoord to be within [-1,1]^dim
538 // since any larger value will be very oscillatory if
539 // called by 'returnNearestElmt' option in
540 // ExpList::GetExpIndex
541 bool clamp = false;
542 for (int i = 0; i < GetShapeDim(); ++i)
543 {
544 if (!std::isfinite(locCoord[i]))
545 {
546 locCoord[i] = 0.;
547 clamp = true;
548 }
549 else if (locCoord[i] < -(1. + tol))
550 {
551 locCoord[i] = -(1. + tol);
552 clamp = true;
553 }
554 else if (locCoord[i] > (1. + tol))
555 {
556 locCoord[i] = 1. + tol;
557 clamp = true;
558 }
559 }
560 return clamp;
561}
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:242
int GetShapeDim() const
Get the object's shape dimension.
Definition: Geometry.h:426

References ASSERTL1, and GetShapeDim().

Referenced by Nektar::SpatialDomains::Geometry2D::NewtonIterationForLocCoord(), Nektar::SpatialDomains::Geometry3D::NewtonIterationForLocCoord(), v_ContainsPoint(), Nektar::SpatialDomains::Geometry2D::v_FindDistance(), Nektar::SpatialDomains::SegGeom::v_FindDistance(), Nektar::SpatialDomains::Geometry2D::v_GetLocCoords(), and Nektar::SpatialDomains::Geometry3D::v_GetLocCoords().

◆ ClearBoundingBox()

void Nektar::SpatialDomains::Geometry::ClearBoundingBox ( )

Definition at line 465 of file Geometry.cpp.

466{
467 m_boundingBox = {};
468}
Array< OneD, NekDouble > m_boundingBox
Array containing bounding box.
Definition: Geometry.h:208

References m_boundingBox.

◆ ContainsPoint() [1/3]

bool Nektar::SpatialDomains::Geometry::ContainsPoint ( const Array< OneD, const NekDouble > &  gloCoord,
Array< OneD, NekDouble > &  locCoord,
NekDouble  tol 
)
inline

Determine whether an element contains a particular Cartesian coordinate \((x,y,z)\).

See also
Geometry::ContainsPoint

Definition at line 481 of file Geometry.h.

484{
485 NekDouble dist;
486 return v_ContainsPoint(gloCoord, locCoord, tol, dist);
487}
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:229
double NekDouble

References v_ContainsPoint().

◆ ContainsPoint() [2/3]

bool Nektar::SpatialDomains::Geometry::ContainsPoint ( const Array< OneD, const NekDouble > &  gloCoord,
Array< OneD, NekDouble > &  locCoord,
NekDouble  tol,
NekDouble dist 
)
inline

Determine whether an element contains a particular Cartesian coordinate \(\vec{x} = (x,y,z)\).

For curvilinear and non-affine elements (i.e. where the Jacobian varies as a function of the standard element coordinates), this is a non-linear optimisation problem that requires the use of a Newton iteration. Note therefore that this can be an expensive operation.

The parameter tol which is by default 0, can be used to expand the coordinate range of the standard element from \([-1,1]^d\) to \([-1-\epsilon,1+\epsilon\) to handle challenging edge cases. The function also returns the local coordinates corresponding to gloCoord that can be used to speed up subsequent searches.

Parameters
gloCoordThe coordinate \( (x,y,z) \).
locCoordOn exit, this is the local coordinate \(\vec{\xi}\) such that \(\chi(\vec{\xi}) = \vec{x}\).
tolThe tolerance used to dictate the bounding box of the standard coordinates \(\vec{\xi}\).
distOn exit, returns the minimum distance between gloCoord and the quadrature points inside the element.
Returns
true if the coordinate gloCoord is contained in the element; false otherwise.
See also
Geometry::GetLocCoords.

Definition at line 517 of file Geometry.h.

520{
521 return v_ContainsPoint(gloCoord, locCoord, tol, dist);
522}

References v_ContainsPoint().

◆ ContainsPoint() [3/3]

bool Nektar::SpatialDomains::Geometry::ContainsPoint ( const Array< OneD, const NekDouble > &  gloCoord,
NekDouble  tol = 0.0 
)
inline

Determine whether an element contains a particular Cartesian coordinate \((x,y,z)\).

See also
Geometry::ContainsPoint

Definition at line 467 of file Geometry.h.

469{
470 Array<OneD, NekDouble> locCoord(GetCoordim(), 0.0);
471 NekDouble dist;
472 return v_ContainsPoint(gloCoord, locCoord, tol, dist);
473}
int GetCoordim() const
Return the coordinate dimension of this object (i.e. the dimension of the space in which this object ...
Definition: Geometry.h:283

References GetCoordim(), and v_ContainsPoint().

◆ FillGeom()

void Nektar::SpatialDomains::Geometry::FillGeom ( )
inline

Populate the coordinate mapping Geometry::m_coeffs information from any children geometry elements.

See also
v_FillGeom()

Definition at line 456 of file Geometry.h.

457{
458 v_FillGeom();
459}
virtual void v_FillGeom()
Populate the coordinate mapping Geometry::m_coeffs information from any children geometry elements.
Definition: Geometry.cpp:355

References v_FillGeom().

◆ FindDistance()

NekDouble Nektar::SpatialDomains::Geometry::FindDistance ( const Array< OneD, const NekDouble > &  xs,
Array< OneD, NekDouble > &  xi 
)
inline

Definition at line 560 of file Geometry.h.

562{
563 return v_FindDistance(xs, xi);
564}
virtual NekDouble v_FindDistance(const Array< OneD, const NekDouble > &xs, Array< OneD, NekDouble > &xi)
Definition: Geometry.cpp:264

References v_FindDistance().

◆ GenGeomFactors()

void Nektar::SpatialDomains::Geometry::GenGeomFactors ( )
inline

Handles generation of geometry factors.

Generate the geometric factors (i.e. derivatives of \(\chi\)) and related metrics.

See also
SpatialDomains::GeomFactors

Definition at line 692 of file Geometry.h.

693{
694 return v_GenGeomFactors();
695}

References v_GenGeomFactors().

Referenced by GetGeomFactors().

◆ GetBoundingBox()

std::array< NekDouble, 6 > Nektar::SpatialDomains::Geometry::GetBoundingBox ( )

Generates the bounding box for the element.

For regular elements, the vertices are sufficient to define the extent of the bounding box. For non-regular elements, the extremes of the quadrature point coordinates are used. A 10% margin is added around this computed region to account for convex hull elements where the true extent of the element may extend slightly beyond the quadrature points.

Definition at line 391 of file Geometry.cpp.

392{
393 if (m_boundingBox.size() == 6)
394 {
395 return {{m_boundingBox[0], m_boundingBox[1], m_boundingBox[2],
397 }
398 // NekDouble minx, miny, minz, maxx, maxy, maxz;
399 Array<OneD, NekDouble> min(3), max(3);
400
401 // Always get vertexes min/max
403 Array<OneD, NekDouble> x(3, 0.0);
404 p->GetCoords(x[0], x[1], x[2]);
405 for (int j = 0; j < 3; ++j)
406 {
407 min[j] = x[j];
408 max[j] = x[j];
409 }
410 for (int i = 1; i < GetNumVerts(); ++i)
411 {
412 p = GetVertex(i);
413 p->GetCoords(x[0], x[1], x[2]);
414 for (int j = 0; j < 3; ++j)
415 {
416 min[j] = (x[j] < min[j] ? x[j] : min[j]);
417 max[j] = (x[j] > max[j] ? x[j] : max[j]);
418 }
419 }
420 // If element is deformed loop over quadrature points
422 if (GetGeomFactors()->GetGtype() != eRegular)
423 {
424 marginFactor = 0.1;
425 const int nq = GetXmap()->GetTotPoints();
426 Array<OneD, Array<OneD, NekDouble>> xvec(3);
427 for (int j = 0; j < 3; ++j)
428 {
429 xvec[j] = Array<OneD, NekDouble>(nq, 0.0);
430 }
431 for (int j = 0; j < GetCoordim(); ++j)
432 {
433 GetXmap()->BwdTrans(m_coeffs[j], xvec[j]);
434 }
435 for (int j = 0; j < 3; ++j)
436 {
437 for (int i = 0; i < nq; ++i)
438 {
439 min[j] = (xvec[j][i] < min[j] ? xvec[j][i] : min[j]);
440 max[j] = (xvec[j][i] > max[j] ? xvec[j][i] : max[j]);
441 }
442 }
443 }
444 // Add margin to bounding box, in order to
445 // return the nearest element
446 for (int j = 0; j < 3; ++j)
447 {
448 NekDouble margin =
449 marginFactor * (max[j] - min[j]) + NekConstants::kFindDistanceMin;
450 min[j] -= margin;
451 max[j] += margin;
452 }
453
454 // save bounding box
455 m_boundingBox = Array<OneD, NekDouble>(6);
456 for (int j = 0; j < 3; ++j)
457 {
458 m_boundingBox[j] = min[j];
459 m_boundingBox[j + 3] = max[j];
460 }
461 // Return bounding box
462 return {{min[0], min[1], min[2], max[0], max[1], max[2]}};
463}
PointGeomSharedPtr GetVertex(int i) const
Returns vertex i of this object.
Definition: Geometry.h:357
Array< OneD, Array< OneD, NekDouble > > m_coeffs
Array containing expansion coefficients of m_xmap.
Definition: Geometry.h:206
StdRegions::StdExpansionSharedPtr GetXmap() const
Return the mapping object Geometry::m_xmap that represents the coordinate transformation from standar...
Definition: Geometry.h:435
int GetNumVerts() const
Get the number of vertices of this object.
Definition: Geometry.h:399
GeomFactorsSharedPtr GetGeomFactors()
Get the geometric factors for this object, generating them if required.
Definition: Geometry.h:301
static const NekDouble kGeomFactorsTol
static const NekDouble kFindDistanceMin
@ eRegular
Geometry is straight-sided with constant geometric factors.
std::shared_ptr< PointGeom > PointGeomSharedPtr
Definition: Geometry.h:57

References Nektar::SpatialDomains::eRegular, GetCoordim(), GetGeomFactors(), GetNumVerts(), GetVertex(), GetXmap(), Nektar::NekConstants::kFindDistanceMin, Nektar::NekConstants::kGeomFactorsTol, m_boundingBox, m_coeffs, and CellMLToNektar.cellml_metadata::p.

Referenced by MinMaxCheck().

◆ GetCoeffs()

const Array< OneD, const NekDouble > & Nektar::SpatialDomains::Geometry::GetCoeffs ( const int  i) const
inline

Return the coefficients of the transformation Geometry::m_xmap in coordinate direction i.

Definition at line 444 of file Geometry.h.

446{
447 return m_coeffs[i];
448}

References m_coeffs.

◆ GetCoord()

NekDouble Nektar::SpatialDomains::Geometry::GetCoord ( const int  i,
const Array< OneD, const NekDouble > &  Lcoord 
)
inline

Given local collapsed coordinate Lcoord, return the value of physical coordinate in direction i.

Definition at line 554 of file Geometry.h.

556{
557 return v_GetCoord(i, Lcoord);
558}
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:331

References v_GetCoord().

Referenced by Nektar::SpatialDomains::Geometry2D::v_FindDistance(), and Nektar::SpatialDomains::SegGeom::v_FindDistance().

◆ GetCoordim()

int Nektar::SpatialDomains::Geometry::GetCoordim ( ) const
inline

Return the coordinate dimension of this object (i.e. the dimension of the space in which this object is embedded).

Definition at line 283 of file Geometry.h.

284{
285 return m_coordim;
286}

References m_coordim.

Referenced by Nektar::SpatialDomains::MeshGraph::CheckRange(), ContainsPoint(), GetBoundingBox(), and Nektar::SpatialDomains::SegGeom::v_FindDistance().

◆ GetDir()

int Nektar::SpatialDomains::Geometry::GetDir ( const int  i,
const int  j = 0 
) const
inline

Returns the element coordinate direction corresponding to a given face coordinate direction.

Definition at line 657 of file Geometry.h.

658{
659 return v_GetDir(faceidx, facedir);
660}
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:320

References v_GetDir().

◆ GetEdge()

Geometry1DSharedPtr Nektar::SpatialDomains::Geometry::GetEdge ( int  i) const
inline

Returns edge i of this object.

Definition at line 365 of file Geometry.h.

366{
367 return v_GetEdge(i);
368}
virtual Geometry1DSharedPtr v_GetEdge(int i) const
Returns edge i of this object.
Definition: Geometry.cpp:128

References v_GetEdge().

Referenced by GetEid(), Nektar::SpatialDomains::HexGeom::SetUpLocalEdges(), Nektar::SpatialDomains::PrismGeom::SetUpLocalEdges(), Nektar::SpatialDomains::PyrGeom::SetUpLocalEdges(), and Nektar::SpatialDomains::TetGeom::SetUpLocalEdges().

◆ GetEdgeFaceMap()

int Nektar::SpatialDomains::Geometry::GetEdgeFaceMap ( int  i,
int  j 
) const
inline

Returns the standard element edge IDs that are connected to a given face.

For example, on a prism, edge 0 is connnected to faces 0 and 1; GetEdgeFaceMap(0,j) would therefore return the values 0 and 1 respectively. We assume that j runs between 0 and 1 inclusive, since every face is connected to precisely two faces for all 3D elements.

This function is used in the construction of the low-energy preconditioner.

Parameters
iThe edge to query connectivity for.
jThe local face index between 0 and 1 connected to this element.
See also
MultiRegions::PreconditionerLowEnergy

Definition at line 626 of file Geometry.h.

627{
628 return v_GetEdgeFaceMap(i, j);
629}
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:298

References v_GetEdgeFaceMap().

◆ GetEdgeNormalToFaceVert()

int Nektar::SpatialDomains::Geometry::GetEdgeNormalToFaceVert ( int  i,
int  j 
) const
inline

Returns the standard lement edge IDs that are normal to a given face vertex.

For example, on a hexahedron, on face 0 at vertices 0,1,2,3 the edges normal to that face are 4,5,6,7, ; so GetEdgeNormalToFaceVert(0,j) would therefore return the values 4, 5, 6 and 7 respectively. We assume that j runs between 0 and 3 inclusive on a quadrilateral face and between 0 and 2 inclusive on a triangular face.

This is used to help set up a length scale normal to an face

Parameters
iThe face to query for the normal edge
jThe local vertex index between 0 and nverts on this face

Definition at line 648 of file Geometry.h.

649{
650 return v_GetEdgeNormalToFaceVert(i, j);
651}
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:309

References v_GetEdgeNormalToFaceVert().

◆ GetEid()

int Nektar::SpatialDomains::Geometry::GetEid ( int  i) const

Get the ID of edge i of this object.

Definition at line 112 of file Geometry.cpp.

113{
114 return GetEdge(i)->GetGlobalID();
115}
Geometry1DSharedPtr GetEdge(int i) const
Returns edge i of this object.
Definition: Geometry.h:365

References GetEdge().

Referenced by GetTid(), Nektar::SpatialDomains::HexGeom::SetUpLocalEdges(), Nektar::SpatialDomains::PrismGeom::SetUpLocalEdges(), Nektar::SpatialDomains::PyrGeom::SetUpLocalEdges(), and Nektar::SpatialDomains::TetGeom::SetUpLocalEdges().

◆ GetEorient()

StdRegions::Orientation Nektar::SpatialDomains::Geometry::GetEorient ( const int  i) const
inline

Returns the orientation of edge i with respect to the ordering of edges in the standard element.

Definition at line 382 of file Geometry.h.

383{
384 return v_GetEorient(i);
385}
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:158

References v_GetEorient().

◆ GetFace()

Geometry2DSharedPtr Nektar::SpatialDomains::Geometry::GetFace ( int  i) const
inline

Returns face i of this object.

Definition at line 373 of file Geometry.h.

374{
375 return v_GetFace(i);
376}
virtual Geometry2DSharedPtr v_GetFace(int i) const
Returns face i of this object.
Definition: Geometry.cpp:138

References v_GetFace().

Referenced by GetFid().

◆ GetFid()

int Nektar::SpatialDomains::Geometry::GetFid ( int  i) const

Get the ID of face i of this object.

Definition at line 120 of file Geometry.cpp.

121{
122 return GetFace(i)->GetGlobalID();
123}
Geometry2DSharedPtr GetFace(int i) const
Returns face i of this object.
Definition: Geometry.h:373

References GetFace().

Referenced by GetTid().

◆ GetForient()

StdRegions::Orientation Nektar::SpatialDomains::Geometry::GetForient ( const int  i) const
inline

Returns the orientation of face i with respect to the ordering of faces in the standard element.

Definition at line 391 of file Geometry.h.

392{
393 return v_GetForient(i);
394}
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:169

References v_GetForient().

◆ GetGeomFactors()

GeomFactorsSharedPtr Nektar::SpatialDomains::Geometry::GetGeomFactors ( )
inline

Get the geometric factors for this object, generating them if required.

Definition at line 301 of file Geometry.h.

302{
305}
void GenGeomFactors()
Handles generation of geometry factors.
Definition: Geometry.h:692
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
GeomFactorsSharedPtr m_geomFactors
Geometric factors.
Definition: Geometry.h:190

References GenGeomFactors(), m_geomFactors, and ValidateRegGeomFactor().

Referenced by GetBoundingBox().

◆ GetGlobalID()

int Nektar::SpatialDomains::Geometry::GetGlobalID ( void  ) const
inline

◆ GetLocCoords()

NekDouble Nektar::SpatialDomains::Geometry::GetLocCoords ( const Array< OneD, const NekDouble > &  coords,
Array< OneD, NekDouble > &  Lcoords 
)
inline

Determine the local collapsed coordinates that correspond to a given Cartesian coordinate for this geometry object.

For curvilinear and non-affine elements (i.e. where the Jacobian varies as a function of the standard element coordinates), this is a non-linear optimisation problem that requires the use of a Newton iteration. Note therefore that this can be an expensive operation.

Note that, clearly, the provided Cartesian coordinate lie outside the element. The function therefore returns the minimum distance from some position in the element to . Lcoords will also be constrained to fit within the range \([-1,1]^d\) where \( d \) is the dimension of the element.

Parameters
coordsInput Cartesian global coordinates
LcoordsCorresponding local coordinates
Returns
Distance between obtained coordinates and provided ones.

Definition at line 544 of file Geometry.h.

546{
547 return v_GetLocCoords(coords, Lcoords);
548}
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:343

References v_GetLocCoords().

Referenced by v_ContainsPoint(), Nektar::SpatialDomains::Geometry2D::v_FindDistance(), and Nektar::SpatialDomains::SegGeom::v_FindDistance().

◆ GetMetricInfo()

GeomFactorsSharedPtr Nektar::SpatialDomains::Geometry::GetMetricInfo ( )
inline

◆ GetNumEdges()

int Nektar::SpatialDomains::Geometry::GetNumEdges ( ) const
inline

Get the number of edges of this object.

Definition at line 407 of file Geometry.h.

408{
409 return v_GetNumEdges();
410}
virtual int v_GetNumEdges() const
Get the number of edges of this object.
Definition: Geometry.cpp:180

References v_GetNumEdges().

◆ GetNumFaces()

int Nektar::SpatialDomains::Geometry::GetNumFaces ( ) const
inline

Get the number of faces of this object.

Definition at line 415 of file Geometry.h.

416{
417 return v_GetNumFaces();
418}
virtual int v_GetNumFaces() const
Get the number of faces of this object.
Definition: Geometry.cpp:188

References v_GetNumFaces().

◆ GetNumVerts()

int Nektar::SpatialDomains::Geometry::GetNumVerts ( ) const
inline

Get the number of vertices of this object.

Definition at line 399 of file Geometry.h.

400{
401 return v_GetNumVerts();
402}
virtual int v_GetNumVerts() const
Get the number of vertices of this object.
Definition: Geometry.cpp:148

References v_GetNumVerts().

Referenced by Nektar::SpatialDomains::MeshGraph::CheckRange(), and GetBoundingBox().

◆ GetRefGeomFactors()

GeomFactorsSharedPtr Nektar::SpatialDomains::Geometry::GetRefGeomFactors ( const Array< OneD, const LibUtilities::BasisSharedPtr > &  tbasis)

◆ GetShapeDim()

int Nektar::SpatialDomains::Geometry::GetShapeDim ( ) const
inline

Get the object's shape dimension.

For example, a segment is one dimensional and quadrilateral is two dimensional.

Definition at line 426 of file Geometry.h.

427{
428 return v_GetShapeDim();
429}
virtual int v_GetShapeDim() const
Get the object's shape dimension.
Definition: Geometry.cpp:196

References v_GetShapeDim().

Referenced by ClampLocCoords(), GetTid(), and v_ContainsPoint().

◆ GetShapeType()

LibUtilities::ShapeType Nektar::SpatialDomains::Geometry::GetShapeType ( void  )
inline

Get the geometric shape type of this object.

Definition at line 318 of file Geometry.h.

319{
320 return m_shapeType;
321}

References m_shapeType.

Referenced by Nektar::SpatialDomains::MeshGraph::CheckRange().

◆ GetTid()

int Nektar::SpatialDomains::Geometry::GetTid ( int  i) const
inline

Get the ID of trace i of this object.

The trace element is the facet one dimension lower than the object; for example, a quadrilateral has four trace segments forming its boundary.

Definition at line 345 of file Geometry.h.

346{
347 const int nDim = GetShapeDim();
348 return nDim == 1 ? GetVid(i)
349 : nDim == 2 ? GetEid(i)
350 : nDim == 3 ? GetFid(i)
351 : 0;
352}
int GetVid(int i) const
Get the ID of vertex i of this object.
Definition: Geometry.cpp:104
int GetFid(int i) const
Get the ID of face i of this object.
Definition: Geometry.cpp:120
int GetEid(int i) const
Get the ID of edge i of this object.
Definition: Geometry.cpp:112

References GetEid(), GetFid(), GetShapeDim(), and GetVid().

◆ GetVertex()

PointGeomSharedPtr Nektar::SpatialDomains::Geometry::GetVertex ( int  i) const
inline

◆ GetVertexEdgeMap()

int Nektar::SpatialDomains::Geometry::GetVertexEdgeMap ( int  i,
int  j 
) const
inline

Returns the standard element edge IDs that are connected to a given vertex.

For example, on a prism, vertex 0 is connnected to edges 0, 3, and 4; GetVertexEdgeMap(0,j) would therefore return the values 0, 1 and 4 respectively. We assume that j runs between 0 and 2 inclusive, which is true for every 3D element asides from the pyramid.

This function is used in the construction of the low-energy preconditioner.

Parameters
iThe vertex to query connectivity for.
jThe local edge index between 0 and 2 connected to this element.
Todo:
Expand to work with pyramid elements.
See also
MultiRegions::PreconditionerLowEnergy

Definition at line 583 of file Geometry.h.

584{
585 return v_GetVertexEdgeMap(i, j);
586}
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:276

References v_GetVertexEdgeMap().

◆ GetVertexFaceMap()

int Nektar::SpatialDomains::Geometry::GetVertexFaceMap ( int  i,
int  j 
) const
inline

Returns the standard element face IDs that are connected to a given vertex.

For example, on a hexahedron, vertex 0 is connnected to faces 0, 1, and 4; GetVertexFaceMap(0,j) would therefore return the values 0, 1 and 4 respectively. We assume that j runs between 0 and 2 inclusive, which is true for every 3D element asides from the pyramid.

This is used in the construction of the low-energy preconditioner.

Parameters
iThe vertex to query connectivity for.
jThe local face index between 0 and 2 connected to this element.
Todo:
Expand to work with pyramid elements.
See also
MultiRegions::PreconditionerLowEnergy

Definition at line 605 of file Geometry.h.

606{
607 return v_GetVertexFaceMap(i, j);
608}
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:287

References v_GetVertexFaceMap().

◆ GetVid()

int Nektar::SpatialDomains::Geometry::GetVid ( int  i) const

◆ GetXmap()

StdRegions::StdExpansionSharedPtr Nektar::SpatialDomains::Geometry::GetXmap ( ) const
inline

Return the mapping object Geometry::m_xmap that represents the coordinate transformation from standard element to physical element.

Definition at line 435 of file Geometry.h.

436{
437 return v_GetXmap();
438}
virtual StdRegions::StdExpansionSharedPtr v_GetXmap() const
Return the mapping object Geometry::m_xmap that represents the coordinate transformation from standar...
Definition: Geometry.cpp:218

References v_GetXmap().

Referenced by GetBoundingBox(), Nektar::SpatialDomains::HexGeom::SetUpXmap(), Nektar::SpatialDomains::PrismGeom::SetUpXmap(), Nektar::SpatialDomains::PyrGeom::SetUpXmap(), Nektar::SpatialDomains::QuadGeom::SetUpXmap(), Nektar::SpatialDomains::TetGeom::SetUpXmap(), Nektar::SpatialDomains::TriGeom::SetUpXmap(), Nektar::SpatialDomains::Geometry2D::v_AllLeftCheck(), Nektar::SpatialDomains::Geometry3D::v_FillGeom(), Nektar::SpatialDomains::QuadGeom::v_FillGeom(), and Nektar::SpatialDomains::TriGeom::v_FillGeom().

◆ MinMaxCheck()

bool Nektar::SpatialDomains::Geometry::MinMaxCheck ( const Array< OneD, const NekDouble > &  gloCoord)

Check if given global coord is within the BoundingBox of the element.

Parameters
coordsInput Cartesian global coordinates
Returns
True if within distance or False otherwise.

Definition at line 507 of file Geometry.cpp.

508{
509 // Validation checks
510 ASSERTL1(gloCoord.size() >= m_coordim,
511 "Expects number of global coordinates supplied to be greater than "
512 "or equal to the mesh dimension.");
513
514 std::array<NekDouble, 6> minMax = GetBoundingBox();
515 for (int i = 0; i < m_coordim; ++i)
516 {
517 if ((gloCoord[i] < minMax[i]) || (gloCoord[i] > minMax[i + 3]))
518 {
519 return false;
520 }
521 }
522 return true;
523}
std::array< NekDouble, 6 > GetBoundingBox()
Generates the bounding box for the element.
Definition: Geometry.cpp:391

References ASSERTL1, GetBoundingBox(), and m_coordim.

Referenced by PreliminaryCheck().

◆ PreliminaryCheck()

int Nektar::SpatialDomains::Geometry::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 elements, this check is unnecessary.

Parameters
coordsInput Cartesian global coordinates
Returns
1 is inside of the element. 0 maybe inside -1 outside of the element

Definition at line 481 of file Geometry.cpp.

482{
483 // bounding box check
484 if (!MinMaxCheck(gloCoord))
485 {
486 return -1;
487 }
488
489 // regular element check
490 if (GetMetricInfo()->GetGtype() == eRegular)
491 {
492 return 0;
493 }
494
495 // All left check for straight edges/plane surfaces
496 return v_AllLeftCheck(gloCoord);
497}
virtual int v_AllLeftCheck(const Array< OneD, const NekDouble > &gloCoord)
Definition: Geometry.cpp:203
GeomFactorsSharedPtr GetMetricInfo()
Get the geometric factors for this object.
Definition: Geometry.h:310
bool MinMaxCheck(const Array< OneD, const NekDouble > &gloCoord)
Check if given global coord is within the BoundingBox of the element.
Definition: Geometry.cpp:507

References Nektar::SpatialDomains::eRegular, GetMetricInfo(), MinMaxCheck(), and v_AllLeftCheck().

Referenced by v_ContainsPoint().

◆ Reset()

void Nektar::SpatialDomains::Geometry::Reset ( CurveMap curvedEdges,
CurveMap curvedFaces 
)
inline

Reset this geometry object: unset the current state, zero Geometry::m_coeffs and remove allocated GeomFactors.

Definition at line 666 of file Geometry.h.

667{
668 v_Reset(curvedEdges, curvedFaces);
669}
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:364

References v_Reset().

◆ ResetNonRecursive()

void Nektar::SpatialDomains::Geometry::ResetNonRecursive ( CurveMap curvedEdges,
CurveMap curvedFaces 
)
inline

Reset this geometry object non-recursively: unset the current state, zero Geometry::m_coeffs and remove allocated GeomFactors.

Definition at line 675 of file Geometry.h.

677{
678 Geometry::v_Reset(curvedEdges, curvedFaces);
679}

References v_Reset().

◆ SetCoordim()

void Nektar::SpatialDomains::Geometry::SetCoordim ( int  coordim)
inline

Sets the coordinate dimension of this object (i.e. the dimension of the space in which this object is embedded).

Definition at line 292 of file Geometry.h.

293{
294 m_coordim = dim;
295}

References m_coordim.

◆ SetGlobalID()

void Nektar::SpatialDomains::Geometry::SetGlobalID ( int  globalid)
inline

Set the ID of this object.

Definition at line 334 of file Geometry.h.

335{
336 m_globalID = globalid;
337}

References m_globalID.

◆ Setup()

void Nektar::SpatialDomains::Geometry::Setup ( )
inline

Definition at line 681 of file Geometry.h.

682{
683 v_Setup();
684}

References v_Setup().

◆ SetUpCoeffs()

void Nektar::SpatialDomains::Geometry::SetUpCoeffs ( const int  nCoeffs)
inlineprotected

◆ v_AllLeftCheck()

int Nektar::SpatialDomains::Geometry::v_AllLeftCheck ( const Array< OneD, const NekDouble > &  gloCoord)
protectedvirtual

Reimplemented in Nektar::SpatialDomains::Geometry2D, and Nektar::SpatialDomains::Geometry3D.

Definition at line 203 of file Geometry.cpp.

205{
206 return 0;
207}

Referenced by PreliminaryCheck().

◆ v_CalculateInverseIsoParam()

void Nektar::SpatialDomains::Geometry::v_CalculateInverseIsoParam ( )
protectedvirtual

Reimplemented in Nektar::SpatialDomains::Geometry2D, and Nektar::SpatialDomains::Geometry3D.

Definition at line 209 of file Geometry.cpp.

210{
212 "This function is only valid for shape type geometries");
213}
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:202

References Nektar::ErrorUtil::efatal, and NEKERROR.

◆ v_ContainsPoint()

bool Nektar::SpatialDomains::Geometry::v_ContainsPoint ( const Array< OneD, const NekDouble > &  gloCoord,
Array< OneD, NekDouble > &  locCoord,
NekDouble  tol,
NekDouble dist 
)
protectedvirtual

Determine whether an element contains a particular Cartesian coordinate \(\vec{x} = (x,y,z)\).

For curvilinear and non-affine elements (i.e. where the Jacobian varies as a function of the standard element coordinates), this is a non-linear optimisation problem that requires the use of a Newton iteration. Note therefore that this can be an expensive operation.

The parameter tol which is by default 0, can be used to expand the coordinate range of the standard element from \([-1,1]^d\) to \([-1-\epsilon,1+\epsilon\) to handle challenging edge cases. The function also returns the local coordinates corresponding to gloCoord that can be used to speed up subsequent searches.

Parameters
gloCoordThe coordinate \( (x,y,z) \).
locCoordOn exit, this is the local coordinate \(\vec{\xi}\) such that \(\chi(\vec{\xi}) = \vec{x}\).
tolThe tolerance used to dictate the bounding box of the standard coordinates \(\vec{\xi}\).
distOn exit, returns the minimum distance between gloCoord and the quadrature points inside the element.
Returns
true if the coordinate gloCoord is contained in the element; false otherwise.
See also
Geometry::GetLocCoords. dist is assigned value for curved elements

Reimplemented in Nektar::SpatialDomains::Geometry0D.

Definition at line 229 of file Geometry.cpp.

232{
233 int inside = PreliminaryCheck(gloCoord);
234 if (inside == -1)
235 {
236 dist = std::numeric_limits<double>::max();
237 return false;
238 }
239 dist = GetLocCoords(gloCoord, locCoord);
240 if (inside == 1)
241 {
242 dist = 0.;
243 return true;
244 }
245 else
246 {
247 Array<OneD, NekDouble> eta(GetShapeDim(), 0.);
248 m_xmap->LocCoordToLocCollapsed(locCoord, eta);
249 if (ClampLocCoords(eta, tol))
250 {
251 if (GetMetricInfo()->GetGtype() == eRegular)
252 {
253 dist = std::numeric_limits<double>::max();
254 }
255 return false;
256 }
257 return 3 != m_coordim ||
260 dist <= tol;
261 }
262}
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:544
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:481
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:530
StdRegions::StdExpansionSharedPtr m_xmap
mapping containing isoparametric transformation.
Definition: Geometry.h:194

References ClampLocCoords(), Nektar::LibUtilities::eQuadrilateral, Nektar::SpatialDomains::eRegular, Nektar::LibUtilities::eTriangle, GetLocCoords(), GetMetricInfo(), GetShapeDim(), m_coordim, m_shapeType, m_xmap, and PreliminaryCheck().

Referenced by ContainsPoint().

◆ v_FillGeom()

void Nektar::SpatialDomains::Geometry::v_FillGeom ( )
protectedvirtual

Populate the coordinate mapping Geometry::m_coeffs information from any children geometry elements.

See also
v_FillGeom()

Reimplemented in Nektar::SpatialDomains::Geometry3D, Nektar::SpatialDomains::QuadGeom, Nektar::SpatialDomains::SegGeom, and Nektar::SpatialDomains::TriGeom.

Definition at line 355 of file Geometry.cpp.

356{
358 "This function is only valid for expansion type geometries");
359}

References Nektar::ErrorUtil::efatal, and NEKERROR.

Referenced by FillGeom(), Nektar::SpatialDomains::Geometry1D::v_GetLocCoords(), and Nektar::SpatialDomains::Geometry2D::v_GetLocCoords().

◆ v_FindDistance()

NekDouble Nektar::SpatialDomains::Geometry::v_FindDistance ( const Array< OneD, const NekDouble > &  xs,
Array< OneD, NekDouble > &  xi 
)
protectedvirtual

Reimplemented in Nektar::SpatialDomains::Geometry2D, and Nektar::SpatialDomains::SegGeom.

Definition at line 264 of file Geometry.cpp.

267{
269 "This function has not been defined for this geometry");
270 return false;
271}

References Nektar::ErrorUtil::efatal, and NEKERROR.

Referenced by FindDistance().

◆ v_GenGeomFactors()

virtual void Nektar::SpatialDomains::Geometry::v_GenGeomFactors ( )
protectedpure virtual

◆ v_GetCoord()

NekDouble Nektar::SpatialDomains::Geometry::v_GetCoord ( const int  i,
const Array< OneD, const NekDouble > &  Lcoord 
)
protectedvirtual

Given local collapsed coordinate Lcoord, return the value of physical coordinate in direction i.

Reimplemented in Nektar::SpatialDomains::Geometry3D, Nektar::SpatialDomains::QuadGeom, Nektar::SpatialDomains::SegGeom, and Nektar::SpatialDomains::TriGeom.

Definition at line 331 of file Geometry.cpp.

334{
336 "This function is only valid for expansion type geometries");
337 return 0.0;
338}

References Nektar::ErrorUtil::efatal, and NEKERROR.

Referenced by GetCoord().

◆ v_GetDir()

int Nektar::SpatialDomains::Geometry::v_GetDir ( const int  i,
const int  j 
) const
protectedvirtual

Returns the element coordinate direction corresponding to a given face coordinate direction.

Reimplemented in Nektar::SpatialDomains::HexGeom, Nektar::SpatialDomains::PrismGeom, Nektar::SpatialDomains::PyrGeom, Nektar::SpatialDomains::QuadGeom, Nektar::SpatialDomains::TetGeom, and Nektar::SpatialDomains::TriGeom.

Definition at line 320 of file Geometry.cpp.

322{
324 "This function has not been defined for this geometry");
325 return 0;
326}

References Nektar::ErrorUtil::efatal, and NEKERROR.

Referenced by GetDir().

◆ v_GetEdge()

Geometry1DSharedPtr Nektar::SpatialDomains::Geometry::v_GetEdge ( int  i) const
protectedvirtual

Returns edge i of this object.

Reimplemented in Nektar::SpatialDomains::Geometry2D, and Nektar::SpatialDomains::Geometry3D.

Definition at line 128 of file Geometry.cpp.

129{
131 "This function is only valid for shape type geometries");
132 return Geometry1DSharedPtr();
133}
std::shared_ptr< Geometry1D > Geometry1DSharedPtr
Definition: Geometry.h:61

References Nektar::ErrorUtil::efatal, and NEKERROR.

Referenced by GetEdge().

◆ v_GetEdgeFaceMap()

int Nektar::SpatialDomains::Geometry::v_GetEdgeFaceMap ( int  i,
int  j 
) const
protectedvirtual

Returns the standard element edge IDs that are connected to a given face.

For example, on a prism, edge 0 is connnected to faces 0 and 1; GetEdgeFaceMap(0,j) would therefore return the values 0 and 1 respectively. We assume that j runs between 0 and 1 inclusive, since every face is connected to precisely two faces for all 3D elements.

This function is used in the construction of the low-energy preconditioner.

Parameters
iThe edge to query connectivity for.
jThe local face index between 0 and 1 connected to this element.
See also
MultiRegions::PreconditionerLowEnergy

Reimplemented in Nektar::SpatialDomains::HexGeom, Nektar::SpatialDomains::PrismGeom, and Nektar::SpatialDomains::TetGeom.

Definition at line 298 of file Geometry.cpp.

300{
302 "This function has not been defined for this geometry");
303 return 0;
304}

References Nektar::ErrorUtil::efatal, and NEKERROR.

Referenced by GetEdgeFaceMap().

◆ v_GetEdgeNormalToFaceVert()

int Nektar::SpatialDomains::Geometry::v_GetEdgeNormalToFaceVert ( const int  i,
const int  j 
) const
protectedvirtual

Returns the standard lement edge IDs that are normal to a given face vertex.

For example, on a hexahedron, on face 0 at vertices 0,1,2,3 the edges normal to that face are 4,5,6,7, ; so GetEdgeNormalToFaceVert(0,j) would therefore return the values 4, 5, 6 and 7 respectively. We assume that j runs between 0 and 3 inclusive on a quadrilateral face and between 0 and 2 inclusive on a triangular face.

This is used to help set up a length scale normal to an face

Parameters
iThe face to query for the normal edge
jThe local vertex index between 0 and nverts on this face

Reimplemented in Nektar::SpatialDomains::HexGeom, Nektar::SpatialDomains::PrismGeom, Nektar::SpatialDomains::PyrGeom, and Nektar::SpatialDomains::TetGeom.

Definition at line 309 of file Geometry.cpp.

311{
313 "This function has not been defined for this geometry");
314 return 0;
315}

References Nektar::ErrorUtil::efatal, and NEKERROR.

Referenced by GetEdgeNormalToFaceVert().

◆ v_GetEorient()

StdRegions::Orientation Nektar::SpatialDomains::Geometry::v_GetEorient ( const int  i) const
protectedvirtual

Returns the orientation of edge i with respect to the ordering of edges in the standard element.

Reimplemented in Nektar::SpatialDomains::Geometry2D, and Nektar::SpatialDomains::Geometry3D.

Definition at line 158 of file Geometry.cpp.

160{
162 "This function is not valid for this geometry.");
164}

References Nektar::ErrorUtil::efatal, Nektar::StdRegions::eForwards, and NEKERROR.

Referenced by GetEorient().

◆ v_GetFace()

Geometry2DSharedPtr Nektar::SpatialDomains::Geometry::v_GetFace ( int  i) const
protectedvirtual

Returns face i of this object.

Reimplemented in Nektar::SpatialDomains::Geometry3D.

Definition at line 138 of file Geometry.cpp.

139{
141 "This function is only valid for shape type geometries");
142 return Geometry2DSharedPtr();
143}
std::shared_ptr< Geometry2D > Geometry2DSharedPtr
Definition: Geometry.h:62

References Nektar::ErrorUtil::efatal, and NEKERROR.

Referenced by GetFace().

◆ v_GetForient()

StdRegions::Orientation Nektar::SpatialDomains::Geometry::v_GetForient ( const int  i) const
protectedvirtual

Returns the orientation of face i with respect to the ordering of faces in the standard element.

Reimplemented in Nektar::SpatialDomains::Geometry3D.

Definition at line 169 of file Geometry.cpp.

171{
173 "This function is not valid for this geometry.");
174 return StdRegions::eFwd;
175}

References Nektar::ErrorUtil::efatal, Nektar::StdRegions::eFwd, and NEKERROR.

Referenced by GetForient().

◆ v_GetLocCoords()

NekDouble Nektar::SpatialDomains::Geometry::v_GetLocCoords ( const Array< OneD, const NekDouble > &  coords,
Array< OneD, NekDouble > &  Lcoords 
)
protectedvirtual

Determine the local collapsed coordinates that correspond to a given Cartesian coordinate for this geometry object.

For curvilinear and non-affine elements (i.e. where the Jacobian varies as a function of the standard element coordinates), this is a non-linear optimisation problem that requires the use of a Newton iteration. Note therefore that this can be an expensive operation.

Note that, clearly, the provided Cartesian coordinate lie outside the element. The function therefore returns the minimum distance from some position in the element to . Lcoords will also be constrained to fit within the range \([-1,1]^d\) where \( d \) is the dimension of the element.

Parameters
coordsInput Cartesian global coordinates
LcoordsCorresponding local coordinates
Returns
Distance between obtained coordinates and provided ones.

Reimplemented in Nektar::SpatialDomains::Geometry1D, Nektar::SpatialDomains::Geometry2D, and Nektar::SpatialDomains::Geometry3D.

Definition at line 343 of file Geometry.cpp.

346{
348 "This function is only valid for expansion type geometries");
349 return 0.0;
350}

References Nektar::ErrorUtil::efatal, and NEKERROR.

Referenced by GetLocCoords().

◆ v_GetNumEdges()

int Nektar::SpatialDomains::Geometry::v_GetNumEdges ( ) const
protectedvirtual

Get the number of edges of this object.

Reimplemented in Nektar::SpatialDomains::Geometry2D, and Nektar::SpatialDomains::Geometry3D.

Definition at line 180 of file Geometry.cpp.

181{
182 return 0;
183}

Referenced by GetNumEdges().

◆ v_GetNumFaces()

int Nektar::SpatialDomains::Geometry::v_GetNumFaces ( ) const
protectedvirtual

Get the number of faces of this object.

Reimplemented in Nektar::SpatialDomains::Geometry3D.

Definition at line 188 of file Geometry.cpp.

189{
190 return 0;
191}

Referenced by GetNumFaces().

◆ v_GetNumVerts()

int Nektar::SpatialDomains::Geometry::v_GetNumVerts ( ) const
protectedvirtual

Get the number of vertices of this object.

Reimplemented in Nektar::SpatialDomains::Geometry2D, Nektar::SpatialDomains::Geometry3D, and Nektar::SpatialDomains::SegGeom.

Definition at line 148 of file Geometry.cpp.

149{
151 "This function is only valid for shape type geometries");
152 return 0;
153}

References Nektar::ErrorUtil::efatal, and NEKERROR.

Referenced by GetNumVerts().

◆ v_GetShapeDim()

int Nektar::SpatialDomains::Geometry::v_GetShapeDim ( ) const
protectedvirtual

Get the object's shape dimension.

For example, a segment is one dimensional and quadrilateral is two dimensional.

Reimplemented in Nektar::SpatialDomains::Geometry0D, Nektar::SpatialDomains::Geometry1D, Nektar::SpatialDomains::Geometry2D, and Nektar::SpatialDomains::Geometry3D.

Definition at line 196 of file Geometry.cpp.

197{
199 "This function is only valid for shape type geometries");
200 return 0;
201}

References Nektar::ErrorUtil::efatal, and NEKERROR.

Referenced by GetShapeDim().

◆ v_GetVertex()

virtual PointGeomSharedPtr Nektar::SpatialDomains::Geometry::v_GetVertex ( int  i) const
protectedpure virtual

◆ v_GetVertexEdgeMap()

int Nektar::SpatialDomains::Geometry::v_GetVertexEdgeMap ( int  i,
int  j 
) const
protectedvirtual

Returns the standard element edge IDs that are connected to a given vertex.

For example, on a prism, vertex 0 is connnected to edges 0, 3, and 4; GetVertexEdgeMap(0,j) would therefore return the values 0, 1 and 4 respectively. We assume that j runs between 0 and 2 inclusive, which is true for every 3D element asides from the pyramid.

This function is used in the construction of the low-energy preconditioner.

Parameters
iThe vertex to query connectivity for.
jThe local edge index between 0 and 2 connected to this element.
Todo:
Expand to work with pyramid elements.
See also
MultiRegions::PreconditionerLowEnergy

Reimplemented in Nektar::SpatialDomains::HexGeom, Nektar::SpatialDomains::PrismGeom, and Nektar::SpatialDomains::TetGeom.

Definition at line 276 of file Geometry.cpp.

278{
280 "This function has not been defined for this geometry");
281 return 0;
282}

References Nektar::ErrorUtil::efatal, and NEKERROR.

Referenced by GetVertexEdgeMap().

◆ v_GetVertexFaceMap()

int Nektar::SpatialDomains::Geometry::v_GetVertexFaceMap ( int  i,
int  j 
) const
protectedvirtual

Returns the standard element face IDs that are connected to a given vertex.

For example, on a hexahedron, vertex 0 is connnected to faces 0, 1, and 4; GetVertexFaceMap(0,j) would therefore return the values 0, 1 and 4 respectively. We assume that j runs between 0 and 2 inclusive, which is true for every 3D element asides from the pyramid.

This is used in the construction of the low-energy preconditioner.

Parameters
iThe vertex to query connectivity for.
jThe local face index between 0 and 2 connected to this element.
Todo:
Expand to work with pyramid elements.
See also
MultiRegions::PreconditionerLowEnergy

Reimplemented in Nektar::SpatialDomains::HexGeom, Nektar::SpatialDomains::PrismGeom, and Nektar::SpatialDomains::TetGeom.

Definition at line 287 of file Geometry.cpp.

289{
291 "This function has not been defined for this geometry");
292 return 0;
293}

References Nektar::ErrorUtil::efatal, and NEKERROR.

Referenced by GetVertexFaceMap().

◆ v_GetXmap()

StdRegions::StdExpansionSharedPtr Nektar::SpatialDomains::Geometry::v_GetXmap ( ) const
protectedvirtual

Return the mapping object Geometry::m_xmap that represents the coordinate transformation from standard element to physical element.

Definition at line 218 of file Geometry.cpp.

219{
220 return m_xmap;
221}

References m_xmap.

Referenced by GetXmap().

◆ v_Reset()

void Nektar::SpatialDomains::Geometry::v_Reset ( CurveMap curvedEdges,
CurveMap curvedFaces 
)
protectedvirtual

◆ v_Setup()

void Nektar::SpatialDomains::Geometry::v_Setup ( )
protectedvirtual

◆ ValidateRegGeomFactor()

GeomFactorsSharedPtr Nektar::SpatialDomains::Geometry::ValidateRegGeomFactor ( GeomFactorsSharedPtr  geomFactor)
staticprotected

Check to see if a geometric factor has already been created that contains the same regular information.

The principle behind this is that many regular (i.e. constant Jacobian) elements have identicial geometric factors. Memory may therefore be reduced by storing only the unique factors.

Parameters
geomFactorThe GeomFactor to check.
Returns
Either the cached GeomFactor or geomFactor.
Todo:
Currently this method is disabled since the lookup is very expensive.

Definition at line 81 of file Geometry.cpp.

83{
84 GeomFactorsSharedPtr returnval = geomFactor;
85
86 return returnval;
87}

Referenced by GetGeomFactors().

Member Data Documentation

◆ m_boundingBox

Array<OneD, NekDouble> Nektar::SpatialDomains::Geometry::m_boundingBox
protected

Array containing bounding box.

Definition at line 208 of file Geometry.h.

Referenced by ClearBoundingBox(), and GetBoundingBox().

◆ m_coeffs

Array<OneD, Array<OneD, NekDouble> > Nektar::SpatialDomains::Geometry::m_coeffs
protected

◆ m_coordim

int Nektar::SpatialDomains::Geometry::m_coordim
protected

Coordinate dimension of this geometry object.

Definition at line 188 of file Geometry.h.

Referenced by Nektar::SpatialDomains::PointGeom::Add(), Nektar::SpatialDomains::Geometry2D::Geometry2D(), Nektar::SpatialDomains::Geometry3D::Geometry3D(), GetCoordim(), Nektar::SpatialDomains::PointGeom::GetCoords(), MinMaxCheck(), Nektar::SpatialDomains::PointGeom::Mult(), Nektar::SpatialDomains::PointGeom::PointGeom(), Nektar::SpatialDomains::QuadGeom::QuadGeom(), Nektar::SpatialDomains::SegGeom::SegGeom(), SetCoordim(), SetUpCoeffs(), Nektar::SpatialDomains::HexGeom::SetUpFaceOrientation(), Nektar::SpatialDomains::PrismGeom::SetUpFaceOrientation(), Nektar::SpatialDomains::PyrGeom::SetUpFaceOrientation(), Nektar::SpatialDomains::TetGeom::SetUpFaceOrientation(), Nektar::SpatialDomains::PointGeom::Sub(), Nektar::SpatialDomains::TriGeom::TriGeom(), Nektar::SpatialDomains::Geometry2D::v_AllLeftCheck(), v_ContainsPoint(), Nektar::SpatialDomains::Geometry3D::v_FillGeom(), Nektar::SpatialDomains::QuadGeom::v_FillGeom(), Nektar::SpatialDomains::SegGeom::v_FillGeom(), Nektar::SpatialDomains::TriGeom::v_FillGeom(), Nektar::SpatialDomains::SegGeom::v_FindDistance(), Nektar::SpatialDomains::HexGeom::v_GenGeomFactors(), Nektar::SpatialDomains::PrismGeom::v_GenGeomFactors(), Nektar::SpatialDomains::PyrGeom::v_GenGeomFactors(), Nektar::SpatialDomains::QuadGeom::v_GenGeomFactors(), Nektar::SpatialDomains::SegGeom::v_GenGeomFactors(), Nektar::SpatialDomains::TetGeom::v_GenGeomFactors(), Nektar::SpatialDomains::TriGeom::v_GenGeomFactors(), Nektar::SpatialDomains::Geometry1D::v_GetLocCoords(), and Nektar::SpatialDomains::Geometry2D::v_GetLocCoords().

◆ m_geomFactors

GeomFactorsSharedPtr Nektar::SpatialDomains::Geometry::m_geomFactors
protected

◆ m_geomFactorsState

GeomState Nektar::SpatialDomains::Geometry::m_geomFactorsState
protected

◆ m_geomType

GeomType Nektar::SpatialDomains::Geometry::m_geomType
protected

Type of geometry.

Definition at line 200 of file Geometry.h.

◆ m_globalID

int Nektar::SpatialDomains::Geometry::m_globalID
protected

◆ m_invIsoParam

Array<OneD, Array<OneD, NekDouble> > Nektar::SpatialDomains::Geometry::m_invIsoParam
protected

◆ m_isoParameter

Array<OneD, Array<OneD, NekDouble> > Nektar::SpatialDomains::Geometry::m_isoParameter
protected

◆ m_regGeomFactorsManager

GeomFactorsVector Nektar::SpatialDomains::Geometry::m_regGeomFactorsManager
staticprotected

Definition at line 185 of file Geometry.h.

◆ m_setupState

bool Nektar::SpatialDomains::Geometry::m_setupState
protected

◆ m_shapeType

LibUtilities::ShapeType Nektar::SpatialDomains::Geometry::m_shapeType
protected

◆ m_state

GeomState Nektar::SpatialDomains::Geometry::m_state
protected

◆ m_straightEdge

int Nektar::SpatialDomains::Geometry::m_straightEdge
protected

◆ m_xmap

StdRegions::StdExpansionSharedPtr Nektar::SpatialDomains::Geometry::m_xmap
protected

\(\chi\) mapping containing isoparametric transformation.

Definition at line 194 of file Geometry.h.

Referenced by Nektar::SpatialDomains::SegGeom::GenerateOneSpaceDimGeom(), Nektar::SpatialDomains::Geometry3D::NewtonIterationForLocCoord(), Nektar::SpatialDomains::Geometry2D::NewtonIterationForLocCoord(), Nektar::SpatialDomains::SegGeom::SegGeom(), Nektar::SpatialDomains::HexGeom::SetUpXmap(), Nektar::SpatialDomains::PrismGeom::SetUpXmap(), Nektar::SpatialDomains::PyrGeom::SetUpXmap(), Nektar::SpatialDomains::QuadGeom::SetUpXmap(), Nektar::SpatialDomains::SegGeom::SetUpXmap(), Nektar::SpatialDomains::TetGeom::SetUpXmap(), Nektar::SpatialDomains::TriGeom::SetUpXmap(), v_ContainsPoint(), Nektar::SpatialDomains::Geometry3D::v_FillGeom(), Nektar::SpatialDomains::QuadGeom::v_FillGeom(), Nektar::SpatialDomains::SegGeom::v_FillGeom(), Nektar::SpatialDomains::TriGeom::v_FillGeom(), Nektar::SpatialDomains::Geometry2D::v_FindDistance(), Nektar::SpatialDomains::SegGeom::v_FindDistance(), Nektar::SpatialDomains::HexGeom::v_GenGeomFactors(), Nektar::SpatialDomains::PrismGeom::v_GenGeomFactors(), Nektar::SpatialDomains::PyrGeom::v_GenGeomFactors(), Nektar::SpatialDomains::QuadGeom::v_GenGeomFactors(), Nektar::SpatialDomains::SegGeom::v_GenGeomFactors(), Nektar::SpatialDomains::TetGeom::v_GenGeomFactors(), Nektar::SpatialDomains::TriGeom::v_GenGeomFactors(), Nektar::SpatialDomains::Geometry3D::v_GetCoord(), Nektar::SpatialDomains::QuadGeom::v_GetCoord(), Nektar::SpatialDomains::SegGeom::v_GetCoord(), Nektar::SpatialDomains::TriGeom::v_GetCoord(), Nektar::SpatialDomains::Geometry1D::v_GetLocCoords(), Nektar::SpatialDomains::Geometry2D::v_GetLocCoords(), Nektar::SpatialDomains::Geometry3D::v_GetLocCoords(), v_GetXmap(), Nektar::SpatialDomains::HexGeom::v_Reset(), Nektar::SpatialDomains::PyrGeom::v_Reset(), Nektar::SpatialDomains::QuadGeom::v_Reset(), Nektar::SpatialDomains::SegGeom::v_Reset(), Nektar::SpatialDomains::TriGeom::v_Reset(), Nektar::SpatialDomains::HexGeom::v_Setup(), Nektar::SpatialDomains::PrismGeom::v_Setup(), Nektar::SpatialDomains::PyrGeom::v_Setup(), Nektar::SpatialDomains::QuadGeom::v_Setup(), Nektar::SpatialDomains::SegGeom::v_Setup(), Nektar::SpatialDomains::TetGeom::v_Setup(), and Nektar::SpatialDomains::TriGeom::v_Setup().