37#ifndef NEKTAR_SPATIALDOMAINS_GEOMETRY_H
38#define NEKTAR_SPATIALDOMAINS_GEOMETRY_H
45#include <unordered_map>
66typedef std::unordered_map<int, CurveSharedPtr>
CurveMap;
72 const std::shared_ptr<Geometry> &lhs,
const std::shared_ptr<Geometry> &rhs);
75 const std::shared_ptr<Geometry> &lhs,
const std::shared_ptr<Geometry> &rhs);
158 NekDouble tol = std::numeric_limits<NekDouble>::epsilon());
171 const int j = 0)
const;
246 virtual int v_GetDir(
const int faceidx,
const int facedir)
const;
266 int nVert =
p->GetNumVerts();
267 std::vector<unsigned int> ids(nVert);
269 for (i = 0; i < nVert; ++i)
271 ids[i] =
p->GetVid(i);
273 std::sort(ids.begin(), ids.end());
349 return nDim == 1 ?
GetVid(i)
669 v_Reset(curvedEdges, curvedFaces);
#define SPATIAL_DOMAINS_EXPORT
Base class for shape geometry information.
virtual int v_GetNumEdges() const
Get the number of edges of this object.
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.
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.
LibUtilities::ShapeType GetShapeType(void)
Get the geometric shape type of this object.
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.
Array< OneD, NekDouble > m_boundingBox
Array containing bounding box.
GeomFactorsSharedPtr GetRefGeomFactors(const Array< OneD, const LibUtilities::BasisSharedPtr > &tbasis)
bool m_setupState
Wether or not the setup routines have been run.
PointGeomSharedPtr GetVertex(int i) const
Returns vertex i of this object.
GeomState m_state
Enumeration to dictate whether coefficients are filled.
virtual void v_CalculateInverseIsoParam()
void SetUpCoeffs(const int nCoeffs)
Initialise the Geometry::m_coeffs array.
virtual int v_AllLeftCheck(const Array< OneD, const NekDouble > &gloCoord)
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...
GeomType m_geomType
Type of geometry.
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...
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.
Geometry1DSharedPtr GetEdge(int i) const
Returns edge i of this object.
void SetGlobalID(int globalid)
Set the ID of this object.
int GetNumFaces() const
Get the number of faces of this object.
int GetShapeDim() const
Get the object's shape dimension.
virtual void v_FillGeom()
Populate the coordinate mapping Geometry::m_coeffs information from any children geometry elements.
Geometry()
Default constructor.
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...
virtual Geometry2DSharedPtr v_GetFace(int i) const
Returns face i of this object.
int GetVid(int i) const
Get the ID of vertex i of this object.
virtual int v_GetNumFaces() const
Get the number of faces of this object.
virtual StdRegions::StdExpansionSharedPtr v_GetXmap() const
Return the mapping object Geometry::m_xmap that represents the coordinate transformation from standar...
Array< OneD, Array< OneD, NekDouble > > m_isoParameter
Array< OneD, Array< OneD, NekDouble > > m_invIsoParam
virtual void v_Reset(CurveMap &curvedEdges, CurveMap &curvedFaces)
Reset this geometry object: unset the current state, zero Geometry::m_coeffs and remove allocated Geo...
bool ClampLocCoords(Array< OneD, NekDouble > &locCoord, NekDouble tol=std::numeric_limits< NekDouble >::epsilon())
Clamp local coords to be within standard regions [-1, 1]^dim.
const Array< OneD, const NekDouble > & GetCoeffs(const int i) const
Return the coefficients of the transformation Geometry::m_xmap in coordinate direction i.
int GetGlobalID(void) const
Get the ID of this object.
int GetCoordim() const
Return the coordinate dimension of this object (i.e. the dimension of the space in which this object ...
int GetFid(int i) const
Get the ID of face i of this object.
virtual int v_GetEdgeFaceMap(int i, int j) const
Returns the standard element edge IDs that are connected to a given face.
LibUtilities::ShapeType m_shapeType
Type of shape.
virtual void v_GenGeomFactors()=0
int GetDir(const int i, const int j=0) const
Returns the element coordinate direction corresponding to a given face coordinate direction.
void FillGeom()
Populate the coordinate mapping Geometry::m_coeffs information from any children geometry elements.
Array< OneD, Array< OneD, NekDouble > > m_coeffs
Array containing expansion coefficients of m_xmap.
virtual ~Geometry()
Default destructor.
int GetVertexEdgeMap(int i, int j) const
Returns the standard element edge IDs that are connected to a given vertex.
GeomState m_geomFactorsState
State of the geometric factors.
void ResetNonRecursive(CurveMap &curvedEdges, CurveMap &curvedFaces)
Reset this geometry object non-recursively: unset the current state, zero Geometry::m_coeffs and remo...
virtual PointGeomSharedPtr v_GetVertex(int i) const =0
GeomFactorsSharedPtr GetMetricInfo()
Get the geometric factors for this object.
StdRegions::StdExpansionSharedPtr m_xmap
mapping containing isoparametric transformation.
virtual Geometry1DSharedPtr v_GetEdge(int i) const
Returns edge i of this object.
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.
StdRegions::StdExpansionSharedPtr GetXmap() const
Return the mapping object Geometry::m_xmap that represents the coordinate transformation from standar...
NekDouble FindDistance(const Array< OneD, const NekDouble > &xs, Array< OneD, NekDouble > &xi)
void GenGeomFactors()
Handles generation of geometry factors.
static GeomFactorsVector m_regGeomFactorsManager
int GetNumEdges() const
Get the number of edges of this object.
std::array< NekDouble, 6 > GetBoundingBox()
Generates the bounding box for the element.
int GetNumVerts() const
Get the number of vertices of this object.
Geometry2DSharedPtr GetFace(int i) const
Returns face i of this object.
int GetVertexFaceMap(int i, int j) const
Returns the standard element face IDs that are connected to a given vertex.
bool ContainsPoint(const Array< OneD, const NekDouble > &gloCoord, NekDouble tol=0.0)
Determine whether an element contains a particular Cartesian coordinate .
int GetEdgeFaceMap(int i, int j) const
Returns the standard element edge IDs that are connected to a given face.
bool MinMaxCheck(const Array< OneD, const NekDouble > &gloCoord)
Check if given global coord is within the BoundingBox of the element.
virtual int v_GetVertexEdgeMap(int i, int j) const
Returns the standard element edge IDs that are connected to a given vertex.
int GetEdgeNormalToFaceVert(int i, int j) const
Returns the standard lement edge IDs that are normal to a given face vertex.
virtual NekDouble v_FindDistance(const Array< OneD, const NekDouble > &xs, Array< OneD, NekDouble > &xi)
GeomFactorsSharedPtr GetGeomFactors()
Get the geometric factors for this object, generating them if required.
static GeomFactorsSharedPtr ValidateRegGeomFactor(GeomFactorsSharedPtr geomFactor)
Check to see if a geometric factor has already been created that contains the same regular informatio...
StdRegions::Orientation GetEorient(const int i) const
Returns the orientation of edge i with respect to the ordering of edges in the standard element.
virtual int v_GetNumVerts() const
Get the number of vertices of this object.
GeomFactorsSharedPtr m_geomFactors
Geometric factors.
int m_coordim
Coordinate dimension of this geometry object.
int GetTid(int i) const
Get the ID of trace i of this object.
virtual int v_GetDir(const int faceidx, const int facedir) const
Returns the element coordinate direction corresponding to a given face coordinate direction.
void SetCoordim(int coordim)
Sets the coordinate dimension of this object (i.e. the dimension of the space in which this object is...
virtual int v_GetShapeDim() const
Get the object's shape dimension.
void Reset(CurveMap &curvedEdges, CurveMap &curvedFaces)
Reset this geometry object: unset the current state, zero Geometry::m_coeffs and remove allocated Geo...
virtual int v_GetVertexFaceMap(int i, int j) const
Returns the standard element face IDs that are connected to a given vertex.
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 .
int GetEid(int i) const
Get the ID of edge i of this object.
StdRegions::Orientation GetForient(const int i) const
Returns the orientation of face i with respect to the ordering of faces in the standard element.
std::vector< GeomFactorsSharedPtr > GeomFactorsVector
A vector of GeomFactor pointers.
bool GlobalIdEquality(const std::shared_ptr< Geometry > &lhs, const std::shared_ptr< Geometry > &rhs)
std::shared_ptr< GeomFactors > GeomFactorsSharedPtr
Pointer to a GeomFactors object.
std::shared_ptr< GeometryVector > GeometryVectorSharedPtr
bool SortByGlobalId(const std::shared_ptr< Geometry > &lhs, const std::shared_ptr< Geometry > &rhs)
Less than operator to sort Geometry objects by global id when sorting STL containers.
std::unordered_set< GeometrySharedPtr > GeometrySet
std::shared_ptr< Curve > CurveSharedPtr
std::unordered_map< int, CurveSharedPtr > CurveMap
GeomType
Indicates the type of element geometry.
static CurveMap NullCurveMap
std::shared_ptr< PointGeom > PointGeomSharedPtr
std::shared_ptr< Geometry2D > Geometry2DSharedPtr
GeomState
Indicates if the geometric information for an element has been populated.
std::shared_ptr< Geometry > GeometrySharedPtr
std::shared_ptr< Geometry1D > Geometry1DSharedPtr
std::vector< GeometrySharedPtr > GeometryVector
std::shared_ptr< StdExpansion > StdExpansionSharedPtr
std::size_t hash_range(Iter first, Iter last)
Unary function that constructs a hash of a Geometry object, based on the vertex IDs.
std::size_t operator()(GeometrySharedPtr const &p) const