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 destructor. More...
 
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...
 
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...
 
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)
 Clamp local coords to be within standard regions [-1, 1]^dim. More...
 
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 Setup ()
 

Protected Member Functions

void GenGeomFactors ()
 Handles generation of geometry factors. More...
 
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)
 
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 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...
 

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...
 

Static Protected Attributes

static GeomFactorsVector m_regGeomFactorsManager
 

Detailed Description

Base class for shape geometry information.

Definition at line 81 of file Geometry.h.

Constructor & Destructor Documentation

◆ Geometry() [1/2]

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

Default constructor.

Definition at line 54 of file Geometry.cpp.

57  m_globalID(-1)
58 {
59 }
bool m_setupState
Wether or not the setup routines have been run.
Definition: Geometry.h:190
GeomState m_state
Enumeration to dictate whether coefficients are filled.
Definition: Geometry.h:188
LibUtilities::ShapeType m_shapeType
Type of shape.
Definition: Geometry.h:194
GeomState m_geomFactorsState
State of the geometric factors.
Definition: Geometry.h:184
int m_coordim
Coordinate dimension of this geometry object.
Definition: Geometry.h:180
@ 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 64 of file Geometry.cpp.

◆ ~Geometry()

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

Default destructor.

Definition at line 74 of file Geometry.cpp.

75 {
76 }

Member Function Documentation

◆ ClampLocCoords()

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

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

Parameters
LcoordsCorresponding local coordinates

Definition at line 496 of file Geometry.cpp.

497 {
498  // Validation checks
499  ASSERTL1(locCoord.size() >= GetShapeDim(),
500  "Expects local coordinates to be same or "
501  "larger than shape dimension.");
502 
503  // If out of range clamp locCoord to be within [-1,1]^dim
504  // since any larger value will be very oscillatory if
505  // called by 'returnNearestElmt' option in
506  // ExpList::GetExpIndex
507  bool clamp = false;
508  for (int i = 0; i < GetShapeDim(); ++i)
509  {
510  if (!std::isfinite(locCoord[i]))
511  {
512  locCoord[i] = 0.;
513  clamp = true;
514  }
515  else if (locCoord[i] < -(1. + tol))
516  {
517  locCoord[i] = -(1. + tol);
518  clamp = true;
519  }
520  else if (locCoord[i] > (1. + tol))
521  {
522  locCoord[i] = 1. + tol;
523  clamp = true;
524  }
525  }
526  return clamp;
527 }
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:249
int GetShapeDim() const
Get the object's shape dimension.
Definition: Geometry.h:414

References ASSERTL1, and GetShapeDim().

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

◆ 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 469 of file Geometry.h.

472 {
473  NekDouble dist;
474  return v_ContainsPoint(gloCoord, locCoord, tol, dist);
475 }
virtual bool v_ContainsPoint(const Array< OneD, const NekDouble > &gloCoord, Array< OneD, NekDouble > &locCoord, NekDouble tol, NekDouble &dist)
Definition: Geometry.cpp:253
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 505 of file Geometry.h.

508 {
509  return v_ContainsPoint(gloCoord, locCoord, tol, dist);
510 }

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 455 of file Geometry.h.

457 {
458  Array<OneD, NekDouble> locCoord(GetCoordim(), 0.0);
459  NekDouble dist;
460  return v_ContainsPoint(gloCoord, locCoord, tol, dist);
461 }
int GetCoordim() const
Return the coordinate dimension of this object (i.e. the dimension of the space in which this object ...
Definition: Geometry.h:271

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 444 of file Geometry.h.

445 {
446  v_FillGeom();
447 }
virtual void v_FillGeom()
Populate the coordinate mapping Geometry::m_coeffs information from any children geometry elements.
Definition: Geometry.cpp:358

References v_FillGeom().

Referenced by export_Geometry().

◆ GenGeomFactors()

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

Handles generation of geometry factors.

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

See also
SpatialDomains::GeomFactors

Definition at line 663 of file Geometry.h.

664 {
665  return v_GenGeomFactors();
666 }

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 394 of file Geometry.cpp.

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

References Nektar::SpatialDomains::eRegular, GetCoordim(), GetGeomFactors(), GetNumVerts(), GetVertex(), GetXmap(), 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 432 of file Geometry.h.

434 {
435  return m_coeffs[i];
436 }

References m_coeffs.

Referenced by export_Geometry().

◆ 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 542 of file Geometry.h.

544 {
545  return v_GetCoord(i, Lcoord);
546 }
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:334

References v_GetCoord().

◆ 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 271 of file Geometry.h.

272 {
273  return m_coordim;
274 }

References m_coordim.

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

◆ 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 639 of file Geometry.h.

640 {
641  return v_GetDir(faceidx, facedir);
642 }
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:323

References v_GetDir().

◆ GetEdge()

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

Returns edge i of this object.

Definition at line 353 of file Geometry.h.

354 {
355  return v_GetEdge(i);
356 }
virtual Geometry1DSharedPtr v_GetEdge(int i) const
Returns edge i of this object.
Definition: Geometry.cpp:163

References v_GetEdge().

Referenced by export_Geometry(), 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 608 of file Geometry.h.

609 {
610  return v_GetEdgeFaceMap(i, j);
611 }
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:301

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 630 of file Geometry.h.

631 {
632  return v_GetEdgeNormalToFaceVert(i, j);
633 }
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:312

References v_GetEdgeNormalToFaceVert().

◆ GetEid()

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

Get the ID of edge i of this object.

Definition at line 147 of file Geometry.cpp.

148 {
149  return GetEdge(i)->GetGlobalID();
150 }
Geometry1DSharedPtr GetEdge(int i) const
Returns edge i of this object.
Definition: Geometry.h:353

References GetEdge().

Referenced by export_Geometry(), 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 370 of file Geometry.h.

371 {
372  return v_GetEorient(i);
373 }
virtual StdRegions::Orientation v_GetEorient(const int i) const
Returns the orientation of edge i with respect to the ordering of edges in the standard element.
Definition: Geometry.cpp:195

References v_GetEorient().

Referenced by export_Geometry().

◆ GetFace()

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

Returns face i of this object.

Definition at line 361 of file Geometry.h.

362 {
363  return v_GetFace(i);
364 }
virtual Geometry2DSharedPtr v_GetFace(int i) const
Returns face i of this object.
Definition: Geometry.cpp:174

References v_GetFace().

Referenced by export_Geometry(), and GetFid().

◆ GetFid()

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

Get the ID of face i of this object.

Definition at line 155 of file Geometry.cpp.

156 {
157  return GetFace(i)->GetGlobalID();
158 }
Geometry2DSharedPtr GetFace(int i) const
Returns face i of this object.
Definition: Geometry.h:361

References GetFace().

Referenced by export_Geometry(), and 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 379 of file Geometry.h.

380 {
381  return v_GetForient(i);
382 }
virtual StdRegions::Orientation v_GetForient(const int i) const
Returns the orientation of face i with respect to the ordering of faces in the standard element.
Definition: Geometry.cpp:206

References v_GetForient().

Referenced by export_Geometry().

◆ GetGeomFactors()

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

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

Definition at line 289 of file Geometry.h.

290 {
291  GenGeomFactors();
293 }
void GenGeomFactors()
Handles generation of geometry factors.
Definition: Geometry.h:663
static GeomFactorsSharedPtr ValidateRegGeomFactor(GeomFactorsSharedPtr geomFactor)
Check to see if a geometric factor has already been created that contains the same regular informatio...
Definition: Geometry.cpp:92
GeomFactorsSharedPtr m_geomFactors
Geometric factors.
Definition: Geometry.h:182

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 532 of file Geometry.h.

534 {
535  return v_GetLocCoords(coords, Lcoords);
536 }
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:346

References v_GetLocCoords().

Referenced by v_ContainsPoint().

◆ GetMetricInfo()

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

Get the geometric factors for this object.

Definition at line 298 of file Geometry.h.

299 {
300  return m_geomFactors;
301 }

References m_geomFactors.

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

◆ GetNumEdges()

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

Get the number of edges of this object.

Definition at line 395 of file Geometry.h.

396 {
397  return v_GetNumEdges();
398 }
virtual int v_GetNumEdges() const
Get the number of edges of this object.
Definition: Geometry.cpp:217

References v_GetNumEdges().

Referenced by export_Geometry().

◆ GetNumFaces()

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

Get the number of faces of this object.

Definition at line 403 of file Geometry.h.

404 {
405  return v_GetNumFaces();
406 }
virtual int v_GetNumFaces() const
Get the number of faces of this object.
Definition: Geometry.cpp:225

References v_GetNumFaces().

Referenced by export_Geometry().

◆ GetNumVerts()

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

Get the number of vertices of this object.

Definition at line 387 of file Geometry.h.

388 {
389  return v_GetNumVerts();
390 }
virtual int v_GetNumVerts() const
Get the number of vertices of this object.
Definition: Geometry.cpp:185

References v_GetNumVerts().

Referenced by Nektar::SpatialDomains::MeshGraph::CheckRange(), export_Geometry(), 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 414 of file Geometry.h.

415 {
416  return v_GetShapeDim();
417 }
virtual int v_GetShapeDim() const
Get the object's shape dimension.
Definition: Geometry.cpp:233

References v_GetShapeDim().

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

◆ GetShapeType()

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

Get the geometric shape type of this object.

Definition at line 306 of file Geometry.h.

307 {
308  return m_shapeType;
309 }

References m_shapeType.

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

◆ 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 333 of file Geometry.h.

334 {
335  const int nDim = GetShapeDim();
336  return nDim == 1 ? GetVid(i)
337  : nDim == 2 ? GetEid(i)
338  : nDim == 3 ? GetFid(i)
339  : 0;
340 }
int GetVid(int i) const
Get the ID of vertex i of this object.
Definition: Geometry.cpp:139
int GetFid(int i) const
Get the ID of face i of this object.
Definition: Geometry.cpp:155
int GetEid(int i) const
Get the ID of edge i of this object.
Definition: Geometry.cpp:147

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

Referenced by export_Geometry().

◆ 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 565 of file Geometry.h.

566 {
567  return v_GetVertexEdgeMap(i, j);
568 }
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:279

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 587 of file Geometry.h.

588 {
589  return v_GetVertexFaceMap(i, j);
590 }
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:290

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 423 of file Geometry.h.

424 {
425  return v_GetXmap();
426 }
virtual StdRegions::StdExpansionSharedPtr v_GetXmap() const
Return the mapping object Geometry::m_xmap that represents the coordinate transformation from standar...
Definition: Geometry.cpp:243

References v_GetXmap().

Referenced by export_Geometry(), 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::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 473 of file Geometry.cpp.

474 {
475  // Validation checks
476  ASSERTL1(gloCoord.size() >= m_coordim,
477  "Expects number of global coordinates supplied to be greater than "
478  "or equal to the mesh dimension.");
479 
480  std::array<NekDouble, 6> minMax = GetBoundingBox();
481  for (int i = 0; i < m_coordim; ++i)
482  {
483  if ((gloCoord[i] < minMax[i]) || (gloCoord[i] > minMax[i + 3]))
484  {
485  return false;
486  }
487  }
488  return true;
489 }
std::array< NekDouble, 6 > GetBoundingBox()
Generates the bounding box for the element.
Definition: Geometry.cpp:394

References ASSERTL1, GetBoundingBox(), and m_coordim.

◆ 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 648 of file Geometry.h.

649 {
650  v_Reset(curvedEdges, curvedFaces);
651 }
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:367

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 280 of file Geometry.h.

281 {
282  m_coordim = dim;
283 }

References m_coordim.

◆ SetGlobalID()

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

Set the ID of this object.

Definition at line 322 of file Geometry.h.

323 {
324  m_globalID = globalid;
325 }

References m_globalID.

◆ Setup()

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

Definition at line 652 of file Geometry.h.

653 {
654  v_Setup();
655 }

References v_Setup().

Referenced by export_Geometry().

◆ SetUpCoeffs()

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

◆ v_ContainsPoint()

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

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.

Reimplemented in Nektar::SpatialDomains::Geometry0D.

Definition at line 253 of file Geometry.cpp.

256 {
257  // Convert to the local (xi) coordinates.
258  dist = GetLocCoords(gloCoord, locCoord);
259  if (dist <= tol + NekConstants::kNekMachineEpsilon)
260  {
261  return true;
262  }
263  Array<OneD, NekDouble> eta(GetShapeDim(), 0.);
264  m_xmap->LocCoordToLocCollapsed(locCoord, eta);
265  if (ClampLocCoords(eta, tol))
266  {
267  m_xmap->LocCollapsedToLocCoord(eta, locCoord);
268  return false;
269  }
270  else
271  {
272  return true;
273  }
274 }
bool ClampLocCoords(Array< OneD, NekDouble > &locCoord, NekDouble tol)
Clamp local coords to be within standard regions [-1, 1]^dim.
Definition: Geometry.cpp:496
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:532
StdRegions::StdExpansionSharedPtr m_xmap
mapping containing isoparametric transformation.
Definition: Geometry.h:186
static const NekDouble kNekMachineEpsilon

References ClampLocCoords(), GetLocCoords(), GetShapeDim(), Nektar::NekConstants::kNekMachineEpsilon, and m_xmap.

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::TriGeom, Nektar::SpatialDomains::SegGeom, Nektar::SpatialDomains::QuadGeom, and Nektar::SpatialDomains::Geometry3D.

Definition at line 358 of file Geometry.cpp.

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

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

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

◆ 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::TriGeom, Nektar::SpatialDomains::SegGeom, Nektar::SpatialDomains::QuadGeom, and Nektar::SpatialDomains::Geometry3D.

Definition at line 334 of file Geometry.cpp.

336 {
337  boost::ignore_unused(i, Lcoord);
339  "This function is only valid for expansion type geometries");
340  return 0.0;
341 }

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::TriGeom, Nektar::SpatialDomains::TetGeom, Nektar::SpatialDomains::QuadGeom, Nektar::SpatialDomains::PyrGeom, Nektar::SpatialDomains::PrismGeom, and Nektar::SpatialDomains::HexGeom.

Definition at line 323 of file Geometry.cpp.

324 {
325  boost::ignore_unused(i, j);
327  "This function has not been defined for this geometry");
328  return 0;
329 }

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::Geometry3D, and Nektar::SpatialDomains::Geometry2D.

Definition at line 163 of file Geometry.cpp.

164 {
165  boost::ignore_unused(i);
167  "This function is only valid for shape type geometries");
168  return Geometry1DSharedPtr();
169 }
std::shared_ptr< Geometry1D > Geometry1DSharedPtr
Definition: Geometry.h:63

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::TetGeom, Nektar::SpatialDomains::PrismGeom, and Nektar::SpatialDomains::HexGeom.

Definition at line 301 of file Geometry.cpp.

302 {
303  boost::ignore_unused(i, j);
305  "This function has not been defined for this geometry");
306  return 0;
307 }

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::TetGeom, Nektar::SpatialDomains::PyrGeom, Nektar::SpatialDomains::PrismGeom, and Nektar::SpatialDomains::HexGeom.

Definition at line 312 of file Geometry.cpp.

313 {
314  boost::ignore_unused(i, j);
316  "This function has not been defined for this geometry");
317  return 0;
318 }

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::Geometry3D, and Nektar::SpatialDomains::Geometry2D.

Definition at line 195 of file Geometry.cpp.

196 {
197  boost::ignore_unused(i);
199  "This function is not valid for this geometry.");
200  return StdRegions::eForwards;
201 }

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 174 of file Geometry.cpp.

175 {
176  boost::ignore_unused(i);
178  "This function is only valid for shape type geometries");
179  return Geometry2DSharedPtr();
180 }
std::shared_ptr< Geometry2D > Geometry2DSharedPtr
Definition: Geometry.h:65

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 206 of file Geometry.cpp.

207 {
208  boost::ignore_unused(i);
210  "This function is not valid for this geometry.");
211  return StdRegions::eFwd;
212 }

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::Geometry3D, Nektar::SpatialDomains::Geometry2D, and Nektar::SpatialDomains::Geometry1D.

Definition at line 346 of file Geometry.cpp.

348 {
349  boost::ignore_unused(coords, Lcoords);
351  "This function is only valid for expansion type geometries");
352  return 0.0;
353 }

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::Geometry3D, and Nektar::SpatialDomains::Geometry2D.

Definition at line 217 of file Geometry.cpp.

218 {
219  return 0;
220 }

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 225 of file Geometry.cpp.

226 {
227  return 0;
228 }

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::SegGeom, Nektar::SpatialDomains::Geometry3D, and Nektar::SpatialDomains::Geometry2D.

Definition at line 185 of file Geometry.cpp.

186 {
188  "This function is only valid for shape type geometries");
189  return 0;
190 }

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::Geometry3D, Nektar::SpatialDomains::Geometry2D, Nektar::SpatialDomains::Geometry1D, and Nektar::SpatialDomains::Geometry0D.

Definition at line 233 of file Geometry.cpp.

234 {
236  "This function is only valid for shape type geometries");
237  return 0;
238 }

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::TetGeom, Nektar::SpatialDomains::PrismGeom, and Nektar::SpatialDomains::HexGeom.

Definition at line 279 of file Geometry.cpp.

280 {
281  boost::ignore_unused(i, j);
283  "This function has not been defined for this geometry");
284  return 0;
285 }

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::TetGeom, Nektar::SpatialDomains::PrismGeom, and Nektar::SpatialDomains::HexGeom.

Definition at line 290 of file Geometry.cpp.

291 {
292  boost::ignore_unused(i, j);
294  "This function has not been defined for this geometry");
295  return 0;
296 }

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 243 of file Geometry.cpp.

244 {
245  return m_xmap;
246 }

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.
Todo:
should this '#if 0' statement be removed?

Definition at line 92 of file Geometry.cpp.

94 {
95  GeomFactorsSharedPtr returnval = geomFactor;
96 
97 /// \todo should this '#if 0' statement be removed?
98 #if 0
99  bool found = false;
100  if (geomFactor->GetGtype() == eRegular)
101  {
102  for (GeomFactorsVectorIter iter = m_regGeomFactorsManager.begin();
103  iter != m_regGeomFactorsManager.end();
104  ++iter)
105  {
106  if (**iter == *geomFactor)
107  {
108  returnval = *iter;
109  found = true;
110  break;
111  }
112  }
113 
114  if (!found)
115  {
116  m_regGeomFactorsManager.push_back(geomFactor);
117  returnval = geomFactor;
118  }
119  }
120 #endif
121  return returnval;
122 }
static GeomFactorsVector m_regGeomFactorsManager
Definition: Geometry.h:177

References Nektar::SpatialDomains::eRegular, and m_regGeomFactorsManager.

Referenced by GetGeomFactors().

Member Data Documentation

◆ m_boundingBox

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

Array containing bounding box.

Definition at line 200 of file Geometry.h.

Referenced by 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 180 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::Geometry3D::v_FillGeom(), Nektar::SpatialDomains::QuadGeom::v_FillGeom(), Nektar::SpatialDomains::SegGeom::v_FillGeom(), Nektar::SpatialDomains::TriGeom::v_FillGeom(), 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(), Nektar::SpatialDomains::Geometry2D::v_GetLocCoords(), and Nektar::SpatialDomains::Geometry3D::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 192 of file Geometry.h.

◆ m_globalID

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

◆ m_regGeomFactorsManager

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

Definition at line 177 of file Geometry.h.

Referenced by ValidateRegGeomFactor().

◆ 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_xmap

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

\(\chi\) mapping containing isoparametric transformation.

Definition at line 186 of file Geometry.h.

Referenced by Nektar::SpatialDomains::SegGeom::GenerateOneSpaceDimGeom(), Nektar::SpatialDomains::Geometry2D::NewtonIterationForLocCoord(), Nektar::SpatialDomains::Geometry3D::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::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().