Nektar++
Public Member Functions | Static Public Member Functions | Static Public Attributes | Protected Member Functions | Private Member Functions | List of all members
Nektar::SpatialDomains::QuadGeom Class Reference

#include <QuadGeom.h>

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

Public Member Functions

 QuadGeom ()
 
 QuadGeom (const QuadGeom &in)
 
 QuadGeom (const int id, const SegGeomSharedPtr edges[], const CurveSharedPtr curve=CurveSharedPtr())
 
 ~QuadGeom ()
 
- Public Member Functions inherited from Nektar::SpatialDomains::Geometry2D
 Geometry2D ()
 
 Geometry2D (const int coordim, CurveSharedPtr curve)
 
virtual ~Geometry2D ()
 
CurveSharedPtr GetCurve ()
 
- Public Member Functions inherited from Nektar::SpatialDomains::Geometry
 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...
 
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 Setup ()
 

Static Public Member Functions

static StdRegions::Orientation GetFaceOrientation (const QuadGeom &face1, const QuadGeom &face2, bool doRot=false, int dir=0, NekDouble angle=0.0, NekDouble tol=1e-8)
 Get the orientation of face1. More...
 
static StdRegions::Orientation GetFaceOrientation (const PointGeomVector &face1, const PointGeomVector &face2, bool doRot=false, int dir=0, NekDouble angle=0.0, NekDouble tol=1e-8)
 

Static Public Attributes

static const int kNverts = 4
 
static const int kNedges = 4
 
static const std::string XMLElementType
 
- Static Public Attributes inherited from Nektar::SpatialDomains::Geometry2D
static const int kDim = 2
 

Protected Member Functions

virtual NekDouble v_GetCoord (const int i, const Array< OneD, const NekDouble > &Lcoord) override
 Given local collapsed coordinate Lcoord, return the value of physical coordinate in direction i. More...
 
virtual void v_GenGeomFactors () override
 
virtual void v_FillGeom () override
 
virtual int v_GetDir (const int faceidx, const int facedir) const override
 Returns the element coordinate direction corresponding to a given face coordinate direction. More...
 
virtual void v_Reset (CurveMap &curvedEdges, CurveMap &curvedFaces) override
 Reset this geometry object: unset the current state, zero Geometry::m_coeffs and remove allocated GeomFactors. More...
 
virtual void v_Setup () override
 
- Protected Member Functions inherited from Nektar::SpatialDomains::Geometry2D
virtual NekDouble v_GetLocCoords (const Array< OneD, const NekDouble > &coords, Array< OneD, NekDouble > &Lcoords) override
 Determine the local collapsed coordinates that correspond to a given Cartesian coordinate for this geometry object. More...
 
void NewtonIterationForLocCoord (const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &ptsx, const Array< OneD, const NekDouble > &ptsy, Array< OneD, NekDouble > &Lcoords, NekDouble &dist)
 
void NewtonIterationForLocCoord (const Array< OneD, const NekDouble > &coords, Array< OneD, NekDouble > &Lcoords)
 
virtual void v_CalculateInverseIsoParam () override
 
virtual int v_AllLeftCheck (const Array< OneD, const NekDouble > &gloCoord) override
 
virtual int v_GetShapeDim () const override
 Get the object's shape dimension. More...
 
virtual PointGeomSharedPtr v_GetVertex (int i) const override
 
virtual Geometry1DSharedPtr v_GetEdge (int i) const override
 Returns edge i of this object. More...
 
virtual int v_GetNumVerts () const override
 Get the number of vertices of this object. More...
 
virtual int v_GetNumEdges () const override
 Get the number of edges of this object. More...
 
virtual StdRegions::Orientation v_GetEorient (const int i) const override
 Returns the orientation of edge i with respect to the ordering of edges in the standard element. More...
 
virtual NekDouble v_FindDistance (const Array< OneD, const NekDouble > &xs, Array< OneD, NekDouble > &xi) override
 
- Protected Member Functions inherited from Nektar::SpatialDomains::Geometry
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)
 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 ()
 

Private Member Functions

void SetUpXmap ()
 

Additional Inherited Members

- Static Protected Member Functions inherited from Nektar::SpatialDomains::Geometry
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 inherited from Nektar::SpatialDomains::Geometry2D
PointGeomVector m_verts
 
SegGeomVector m_edges
 
std::vector< StdRegions::Orientationm_eorient
 
CurveSharedPtr m_curve
 
Array< OneD, int > m_manifold
 
Array< OneD, Array< OneD, NekDouble > > m_edgeNormal
 
- Protected Attributes inherited from Nektar::SpatialDomains::Geometry
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
 
bool m_straightEdge
 
- Static Protected Attributes inherited from Nektar::SpatialDomains::Geometry
static GeomFactorsVector m_regGeomFactorsManager
 

Detailed Description

Definition at line 56 of file QuadGeom.h.

Constructor & Destructor Documentation

◆ QuadGeom() [1/3]

Nektar::SpatialDomains::QuadGeom::QuadGeom ( )

◆ QuadGeom() [2/3]

Nektar::SpatialDomains::QuadGeom::QuadGeom ( const QuadGeom in)

Definition at line 88 of file QuadGeom.cpp.

88 : Geometry2D(in)
89{
90 // From Geometry
91 m_shapeType = in.m_shapeType;
92 m_globalID = in.m_globalID;
93
94 // From QuadGeom
95 m_verts = in.m_verts;
96 m_edges = in.m_edges;
97 for (int i = 0; i < kNedges; i++)
98 {
99 m_eorient[i] = in.m_eorient[i];
100 }
101}
std::vector< StdRegions::Orientation > m_eorient
Definition: Geometry2D.h:84
static const int kNedges
Definition: QuadGeom.h:76

References kNedges, Nektar::SpatialDomains::Geometry2D::m_edges, Nektar::SpatialDomains::Geometry2D::m_eorient, Nektar::SpatialDomains::Geometry::m_globalID, Nektar::SpatialDomains::Geometry::m_shapeType, and Nektar::SpatialDomains::Geometry2D::m_verts.

◆ QuadGeom() [3/3]

Nektar::SpatialDomains::QuadGeom::QuadGeom ( const int  id,
const SegGeomSharedPtr  edges[],
const CurveSharedPtr  curve = CurveSharedPtr() 
)

Copy the edge shared pointers.

Definition at line 56 of file QuadGeom.cpp.

58 : Geometry2D(edges[0]->GetVertex(0)->GetCoordim(), curve)
59{
60 int j;
61
63 m_globalID = id;
64
65 /// Copy the edge shared pointers.
66 m_edges.insert(m_edges.begin(), edges, edges + QuadGeom::kNedges);
67 m_eorient.resize(kNedges);
68
69 for (j = 0; j < kNedges; ++j)
70 {
71 m_eorient[j] =
72 SegGeom::GetEdgeOrientation(*edges[j], *edges[(j + 1) % kNedges]);
73 m_verts.push_back(
74 edges[j]->GetVertex(m_eorient[j] == StdRegions::eForwards ? 0 : 1));
75 }
76
77 for (j = 2; j < kNedges; ++j)
78 {
82 }
83
84 m_coordim = edges[0]->GetVertex(0)->GetCoordim();
85 ASSERTL0(m_coordim > 1, "Cannot call function with dim == 1");
86}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
PointGeomSharedPtr GetVertex(int i) const
Returns vertex i of this object.
Definition: Geometry.h:358
int GetCoordim() const
Return the coordinate dimension of this object (i.e. the dimension of the space in which this object ...
Definition: Geometry.h:284
int m_coordim
Coordinate dimension of this geometry object.
Definition: Geometry.h:186
static StdRegions::Orientation GetEdgeOrientation(const SegGeom &edge1, const SegGeom &edge2)
Get the orientation of edge1.
Definition: SegGeom.cpp:196

References ASSERTL0, Nektar::StdRegions::eBackwards, Nektar::StdRegions::eForwards, Nektar::LibUtilities::eQuadrilateral, Nektar::SpatialDomains::SegGeom::GetEdgeOrientation(), Nektar::SpatialDomains::Geometry::GetVertex(), kNedges, Nektar::SpatialDomains::Geometry::m_coordim, Nektar::SpatialDomains::Geometry2D::m_edges, Nektar::SpatialDomains::Geometry2D::m_eorient, Nektar::SpatialDomains::Geometry::m_globalID, Nektar::SpatialDomains::Geometry::m_shapeType, and Nektar::SpatialDomains::Geometry2D::m_verts.

◆ ~QuadGeom()

Nektar::SpatialDomains::QuadGeom::~QuadGeom ( )

Definition at line 103 of file QuadGeom.cpp.

104{
105}

Member Function Documentation

◆ GetFaceOrientation() [1/2]

StdRegions::Orientation Nektar::SpatialDomains::QuadGeom::GetFaceOrientation ( const PointGeomVector face1,
const PointGeomVector face2,
bool  doRot = false,
int  dir = 0,
NekDouble  angle = 0.0,
NekDouble  tol = 1e-8 
)
static

Calculate the orientation of face2 to face1 (note this is not face1 to face2!).

Definition at line 151 of file QuadGeom.cpp.

154{
155 int i, j, vmap[4] = {-1, -1, -1, -1};
156
157 if (doRot)
158 {
159 PointGeom rotPt;
160
161 for (i = 0; i < 4; ++i)
162 {
163 rotPt.Rotate((*face1[i]), dir, angle);
164 for (j = 0; j < 4; ++j)
165 {
166 if (rotPt.dist(*face2[j]) < tol)
167 {
168 vmap[j] = i;
169 break;
170 }
171 }
172 }
173 }
174 else
175 {
176
177 NekDouble x, y, z, x1, y1, z1, cx = 0.0, cy = 0.0, cz = 0.0;
178
179 // For periodic faces, we calculate the vector between the centre
180 // points of the two faces. (For connected faces this will be
181 // zero). We can then use this to determine alignment later in the
182 // algorithm.
183 for (i = 0; i < 4; ++i)
184 {
185 cx += (*face2[i])(0) - (*face1[i])(0);
186 cy += (*face2[i])(1) - (*face1[i])(1);
187 cz += (*face2[i])(2) - (*face1[i])(2);
188 }
189 cx /= 4;
190 cy /= 4;
191 cz /= 4;
192
193 // Now construct a mapping which takes us from the vertices of one
194 // face to the other. That is, vertex j of face2 corresponds to
195 // vertex vmap[j] of face1.
196 for (i = 0; i < 4; ++i)
197 {
198 x = (*face1[i])(0);
199 y = (*face1[i])(1);
200 z = (*face1[i])(2);
201 for (j = 0; j < 4; ++j)
202 {
203 x1 = (*face2[j])(0) - cx;
204 y1 = (*face2[j])(1) - cy;
205 z1 = (*face2[j])(2) - cz;
206 if (sqrt((x1 - x) * (x1 - x) + (y1 - y) * (y1 - y) +
207 (z1 - z) * (z1 - z)) < 1e-8)
208 {
209 vmap[j] = i;
210 break;
211 }
212 }
213 }
214 }
215
216 // Use the mapping to determine the eight alignment options between
217 // faces.
218 if (vmap[1] == (vmap[0] + 1) % 4)
219 {
220 switch (vmap[0])
221 {
222 case 0:
224 break;
225 case 1:
227 break;
228 case 2:
230 break;
231 case 3:
233 break;
234 }
235 }
236 else
237 {
238 switch (vmap[0])
239 {
240 case 0:
242 break;
243 case 1:
245 break;
246 case 2:
248 break;
249 case 3:
251 break;
252 }
253 }
254
255 ASSERTL0(false, "unable to determine face orientation");
257}
std::vector< double > z(NPUPPER)
double NekDouble
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294

References ASSERTL0, Nektar::SpatialDomains::PointGeom::dist(), Nektar::StdRegions::eDir1BwdDir1_Dir2BwdDir2, Nektar::StdRegions::eDir1BwdDir1_Dir2FwdDir2, Nektar::StdRegions::eDir1BwdDir2_Dir2BwdDir1, Nektar::StdRegions::eDir1BwdDir2_Dir2FwdDir1, Nektar::StdRegions::eDir1FwdDir1_Dir2BwdDir2, Nektar::StdRegions::eDir1FwdDir1_Dir2FwdDir2, Nektar::StdRegions::eDir1FwdDir2_Dir2BwdDir1, Nektar::StdRegions::eDir1FwdDir2_Dir2FwdDir1, Nektar::SpatialDomains::PointGeom::Rotate(), tinysimd::sqrt(), and Nektar::UnitTests::z().

◆ GetFaceOrientation() [2/2]

StdRegions::Orientation Nektar::SpatialDomains::QuadGeom::GetFaceOrientation ( const QuadGeom face1,
const QuadGeom face2,
bool  doRot = false,
int  dir = 0,
NekDouble  angle = 0.0,
NekDouble  tol = 1e-8 
)
static

Get the orientation of face1.

Definition at line 137 of file QuadGeom.cpp.

142{
143 return GetFaceOrientation(face1.m_verts, face2.m_verts, doRot, dir, angle,
144 tol);
145}
static StdRegions::Orientation GetFaceOrientation(const QuadGeom &face1, const QuadGeom &face2, bool doRot=false, int dir=0, NekDouble angle=0.0, NekDouble tol=1e-8)
Get the orientation of face1.
Definition: QuadGeom.cpp:137

References GetFaceOrientation(), and Nektar::SpatialDomains::Geometry2D::m_verts.

Referenced by Nektar::MultiRegions::DisContField::FindPeriodicTraces(), and GetFaceOrientation().

◆ SetUpXmap()

void Nektar::SpatialDomains::QuadGeom::SetUpXmap ( )
private

Definition at line 107 of file QuadGeom.cpp.

108{
109 int order0 = max(m_edges[0]->GetXmap()->GetBasis(0)->GetNumModes(),
110 m_edges[2]->GetXmap()->GetBasis(0)->GetNumModes());
111 int order1 = max(m_edges[1]->GetXmap()->GetBasis(0)->GetNumModes(),
112 m_edges[3]->GetXmap()->GetBasis(0)->GetNumModes());
113
114 const LibUtilities::BasisKey B0(
116 LibUtilities::PointsKey(order0 + 1,
118 const LibUtilities::BasisKey B1(
120 LibUtilities::PointsKey(order1 + 1,
122
124}
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
StdRegions::StdExpansionSharedPtr m_xmap
mapping containing isoparametric transformation.
Definition: Geometry.h:192
StdRegions::StdExpansionSharedPtr GetXmap() const
Return the mapping object Geometry::m_xmap that represents the coordinate transformation from standar...
Definition: Geometry.h:436
@ eGaussLobattoLegendre
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:53
@ eModified_A
Principle Modified Functions .
Definition: BasisType.h:50

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::LibUtilities::eGaussLobattoLegendre, Nektar::LibUtilities::eModified_A, Nektar::SpatialDomains::Geometry::GetXmap(), Nektar::SpatialDomains::Geometry2D::m_edges, and Nektar::SpatialDomains::Geometry::m_xmap.

Referenced by v_Reset(), and v_Setup().

◆ v_FillGeom()

void Nektar::SpatialDomains::QuadGeom::v_FillGeom ( )
overrideprotectedvirtual

Note verts and edges are listed according to anticlockwise convention but points in _coeffs have to be in array format from left to right.

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 386 of file QuadGeom.cpp.

387{
388 // check to see if geometry structure is already filled
389 if (m_state != ePtsFilled)
390 {
391 int i, j, k;
392 int nEdgeCoeffs;
393
394 if (m_curve)
395 {
396 int npts = m_curve->m_points.size();
397 int nEdgePts = (int)sqrt(static_cast<NekDouble>(npts));
398 Array<OneD, NekDouble> tmp(npts);
399 Array<OneD, NekDouble> tmp2(m_xmap->GetTotPoints());
400 LibUtilities::PointsKey curveKey(nEdgePts, m_curve->m_ptype);
401
402 // Sanity checks:
403 // - Curved faces should have square number of points;
404 // - Each edge should have sqrt(npts) points.
405 ASSERTL0(nEdgePts * nEdgePts == npts,
406 "NUMPOINTS should be a square number in"
407 " quadrilteral " +
408 boost::lexical_cast<string>(m_globalID));
409
410 for (i = 0; i < kNedges; ++i)
411 {
412 ASSERTL0(m_edges[i]->GetXmap()->GetNcoeffs() == nEdgePts,
413 "Number of edge points does not correspond to "
414 "number of face points in quadrilateral " +
415 boost::lexical_cast<string>(m_globalID));
416 }
417
418 for (i = 0; i < m_coordim; ++i)
419 {
420 for (j = 0; j < npts; ++j)
421 {
422 tmp[j] = (m_curve->m_points[j]->GetPtr())[i];
423 }
424
425 // Interpolate m_curve points to GLL points
426 LibUtilities::Interp2D(curveKey, curveKey, tmp,
427 m_xmap->GetBasis(0)->GetPointsKey(),
428 m_xmap->GetBasis(1)->GetPointsKey(),
429 tmp2);
430
431 // Forwards transform to get coefficient space.
432 m_xmap->FwdTrans(tmp2, m_coeffs[i]);
433 }
434 }
435
436 // Now fill in edges.
437 Array<OneD, unsigned int> mapArray;
438 Array<OneD, int> signArray;
439
440 for (i = 0; i < kNedges; i++)
441 {
442 m_edges[i]->FillGeom();
443 m_xmap->GetTraceToElementMap(i, mapArray, signArray, m_eorient[i]);
444
445 nEdgeCoeffs = m_edges[i]->GetXmap()->GetNcoeffs();
446
447 for (j = 0; j < m_coordim; j++)
448 {
449 for (k = 0; k < nEdgeCoeffs; k++)
450 {
451 m_coeffs[j][mapArray[k]] =
452 signArray[k] * (m_edges[i]->GetCoeffs(j))[k];
453 }
454 }
455 }
456
458 }
459}
GeomState m_state
Enumeration to dictate whether coefficients are filled.
Definition: Geometry.h:194
Array< OneD, Array< OneD, NekDouble > > m_coeffs
Array containing expansion coefficients of m_xmap.
Definition: Geometry.h:204
void Interp2D(const BasisKey &fbasis0, const BasisKey &fbasis1, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, const BasisKey &tbasis1, Array< OneD, NekDouble > &to)
this function interpolates a 2D function evaluated at the quadrature points of the 2D basis,...
Definition: Interp.cpp:103
@ ePtsFilled
Geometric information has been generated.

References ASSERTL0, Nektar::SpatialDomains::ePtsFilled, Nektar::SpatialDomains::Geometry::GetXmap(), Nektar::LibUtilities::Interp2D(), kNedges, Nektar::SpatialDomains::Geometry::m_coeffs, Nektar::SpatialDomains::Geometry::m_coordim, Nektar::SpatialDomains::Geometry2D::m_curve, Nektar::SpatialDomains::Geometry2D::m_edges, Nektar::SpatialDomains::Geometry2D::m_eorient, Nektar::SpatialDomains::Geometry::m_globalID, Nektar::SpatialDomains::Geometry::m_state, Nektar::SpatialDomains::Geometry::m_xmap, and tinysimd::sqrt().

Referenced by v_GenGeomFactors().

◆ v_GenGeomFactors()

void Nektar::SpatialDomains::QuadGeom::v_GenGeomFactors ( )
overrideprotectedvirtual

Set up GeoFac for this geometry using Coord quadrature distribution

Implements Nektar::SpatialDomains::Geometry.

Definition at line 262 of file QuadGeom.cpp.

263{
264 if (!m_setupState)
265 {
267 }
268
270 {
271 GeomType Gtype = eRegular;
272
274
275 // We will first check whether we have a regular or deformed
276 // geometry. We will define regular as those cases where the
277 // Jacobian and the metric terms of the derivative are constants
278 // (i.e. not coordinate dependent)
279
280 // Check to see if expansions are linear
281 // If not linear => deformed geometry
282 m_straightEdge = true;
283 if ((m_xmap->GetBasisNumModes(0) != 2) ||
284 (m_xmap->GetBasisNumModes(1) != 2))
285 {
286 Gtype = eDeformed;
287 m_straightEdge = false;
288 }
289
290 // For linear expansions, the mapping from standard to local
291 // element is given by the relation:
292 // x_i = 0.25 * [ ( x_i^A + x_i^B + x_i^C + x_i^D) +
293 // (-x_i^A + x_i^B + x_i^C - x_i^D)*xi_1 +
294 // (-x_i^A - x_i^B + x_i^C + x_i^D)*xi_2 +
295 // ( x_i^A - x_i^B + x_i^C - x_i^D)*xi_1*xi_2 ]
296 //
297 // The jacobian of the transformation and the metric terms
298 // dxi_i/dx_j, involve only terms of the form dx_i/dxi_j (both
299 // for coordim == 2 or 3). Inspecting the formula above, it can
300 // be appreciated that the derivatives dx_i/dxi_j will be
301 // constant, if the coefficient of the non-linear term is zero.
302 //
303 // That is why for regular geometry, we require
304 //
305 // x_i^A - x_i^B + x_i^C - x_i^D = 0
306 //
307 // or equivalently
308 //
309 // x_i^A - x_i^B = x_i^D - x_i^C
310 //
311 // This corresponds to quadrilaterals which are paralellograms.
312 m_manifold = Array<OneD, int>(2);
313 m_manifold[0] = 0;
314 m_manifold[1] = 1;
315 if (m_coordim == 3)
316 {
317 PointGeom e01, e21, norm;
318 e01.Sub(*m_verts[0], *m_verts[1]);
319 e21.Sub(*m_verts[3], *m_verts[1]);
320 norm.Mult(e01, e21);
321 int tmpi = 0;
322 double tmp = std::fabs(norm[0]);
323 if (tmp < fabs(norm[1]))
324 {
325 tmp = fabs(norm[1]);
326 tmpi = 1;
327 }
328 if (tmp < fabs(norm[2]))
329 {
330 tmpi = 2;
331 }
332 m_manifold[0] = (tmpi + 1) % 3;
333 m_manifold[1] = (tmpi + 2) % 3;
334 }
335
336 if (Gtype == eRegular)
337 {
338 Array<OneD, Array<OneD, NekDouble>> verts(m_verts.size());
339 for (int i = 0; i < m_verts.size(); ++i)
340 {
341 verts[i] = Array<OneD, NekDouble>(3);
342 m_verts[i]->GetCoords(verts[i]);
343 }
344 // a00 + a01 xi1 + a02 xi2 + a03 xi1 xi2
345 // a10 + a11 xi1 + a12 xi2 + a03 xi1 xi2
346 m_isoParameter = Array<OneD, Array<OneD, NekDouble>>(2);
347 for (int i = 0; i < 2; i++)
348 {
349 unsigned int d = m_manifold[i];
350 m_isoParameter[i] = Array<OneD, NekDouble>(4, 0.);
351 // Karniadakis, Sherwin 2005, Appendix D
352 NekDouble A = verts[0][d];
353 NekDouble B = verts[1][d];
354 NekDouble D = verts[2][d];
355 NekDouble C = verts[3][d];
356 m_isoParameter[i][0] = 0.25 * (A + B + C + D); // 1
357 m_isoParameter[i][1] = 0.25 * (-A + B - C + D); // xi1
358 m_isoParameter[i][2] = 0.25 * (-A - B + C + D); // xi2
359 m_isoParameter[i][3] = 0.25 * (A - B - C + D); // xi1*xi2
360 NekDouble tmp =
361 fabs(m_isoParameter[i][1]) + fabs(m_isoParameter[i][2]);
362 if (fabs(m_isoParameter[i][3]) >
364 {
365 Gtype = eDeformed;
366 }
367 }
368 }
369
370 if (Gtype == eRegular)
371 {
373 }
374
376 Gtype, m_coordim, m_xmap, m_coeffs);
378 }
379}
Array< OneD, int > m_manifold
Definition: Geometry2D.h:86
virtual void v_CalculateInverseIsoParam() override
Definition: Geometry2D.cpp:625
bool m_setupState
Wether or not the setup routines have been run.
Definition: Geometry.h:196
Array< OneD, Array< OneD, NekDouble > > m_isoParameter
Definition: Geometry.h:207
GeomState m_geomFactorsState
State of the geometric factors.
Definition: Geometry.h:190
GeomFactorsSharedPtr m_geomFactors
Geometric factors.
Definition: Geometry.h:188
virtual void v_Setup() override
Definition: QuadGeom.cpp:487
virtual void v_FillGeom() override
Definition: QuadGeom.cpp:386
static const NekDouble kNekZeroTol
GeomType
Indicates the type of element geometry.
@ eRegular
Geometry is straight-sided with constant geometric factors.
@ eDeformed
Geometry is curved or has non-constant factors.
std::vector< double > d(NPUPPER *NPUPPER)

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::UnitTests::d(), Nektar::SpatialDomains::eDeformed, Nektar::SpatialDomains::ePtsFilled, Nektar::SpatialDomains::eRegular, Nektar::NekConstants::kNekZeroTol, Nektar::SpatialDomains::Geometry::m_coeffs, Nektar::SpatialDomains::Geometry::m_coordim, Nektar::SpatialDomains::Geometry::m_geomFactors, Nektar::SpatialDomains::Geometry::m_geomFactorsState, Nektar::SpatialDomains::Geometry::m_isoParameter, Nektar::SpatialDomains::Geometry2D::m_manifold, Nektar::SpatialDomains::Geometry::m_setupState, Nektar::SpatialDomains::Geometry::m_straightEdge, Nektar::SpatialDomains::Geometry2D::m_verts, Nektar::SpatialDomains::Geometry::m_xmap, Nektar::SpatialDomains::PointGeom::Mult(), Nektar::SpatialDomains::PointGeom::Sub(), Nektar::SpatialDomains::Geometry2D::v_CalculateInverseIsoParam(), v_FillGeom(), and v_Setup().

◆ v_GetCoord()

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

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

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 126 of file QuadGeom.cpp.

128{
129 ASSERTL1(m_state == ePtsFilled, "Geometry is not in physical space");
130
131 Array<OneD, NekDouble> tmp(m_xmap->GetTotPoints());
132 m_xmap->BwdTrans(m_coeffs[i], tmp);
133
134 return m_xmap->PhysEvaluate(Lcoord, tmp);
135}
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:249

References ASSERTL1, Nektar::SpatialDomains::ePtsFilled, Nektar::SpatialDomains::Geometry::m_coeffs, Nektar::SpatialDomains::Geometry::m_state, and Nektar::SpatialDomains::Geometry::m_xmap.

◆ v_GetDir()

int Nektar::SpatialDomains::QuadGeom::v_GetDir ( const int  i,
const int  j 
) const
overrideprotectedvirtual

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

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 461 of file QuadGeom.cpp.

462{
463 boost::ignore_unused(j); // required in 3D shapes
464
465 return i % 2;
466}

◆ v_Reset()

void Nektar::SpatialDomains::QuadGeom::v_Reset ( CurveMap curvedEdges,
CurveMap curvedFaces 
)
overrideprotectedvirtual

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

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 468 of file QuadGeom.cpp.

469{
470 Geometry::v_Reset(curvedEdges, curvedFaces);
471 CurveMap::iterator it = curvedFaces.find(m_globalID);
472
473 if (it != curvedFaces.end())
474 {
475 m_curve = it->second;
476 }
477
478 for (int i = 0; i < 4; ++i)
479 {
480 m_edges[i]->Reset(curvedEdges, curvedFaces);
481 }
482
483 SetUpXmap();
484 SetUpCoeffs(m_xmap->GetNcoeffs());
485}
void SetUpCoeffs(const int nCoeffs)
Initialise the Geometry::m_coeffs array.
Definition: Geometry.h:690
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:390

References Nektar::SpatialDomains::Geometry2D::m_curve, Nektar::SpatialDomains::Geometry2D::m_edges, Nektar::SpatialDomains::Geometry::m_globalID, Nektar::SpatialDomains::Geometry::m_xmap, Nektar::SpatialDomains::Geometry::SetUpCoeffs(), SetUpXmap(), and Nektar::SpatialDomains::Geometry::v_Reset().

◆ v_Setup()

void Nektar::SpatialDomains::QuadGeom::v_Setup ( )
overrideprotectedvirtual

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 487 of file QuadGeom.cpp.

488{
489 if (!m_setupState)
490 {
491 for (int i = 0; i < 4; ++i)
492 {
493 m_edges[i]->Setup();
494 }
495 SetUpXmap();
496 SetUpCoeffs(m_xmap->GetNcoeffs());
497 m_setupState = true;
498 }
499}

References Nektar::SpatialDomains::Geometry2D::m_edges, Nektar::SpatialDomains::Geometry::m_setupState, Nektar::SpatialDomains::Geometry::m_xmap, Nektar::SpatialDomains::Geometry::SetUpCoeffs(), and SetUpXmap().

Referenced by v_GenGeomFactors().

Member Data Documentation

◆ kNedges

const int Nektar::SpatialDomains::QuadGeom::kNedges = 4
static

◆ kNverts

const int Nektar::SpatialDomains::QuadGeom::kNverts = 4
static

◆ XMLElementType

const std::string Nektar::SpatialDomains::QuadGeom::XMLElementType
static

Definition at line 77 of file QuadGeom.h.