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 () override
 
- Public Member Functions inherited from Nektar::SpatialDomains::Geometry2D
 Geometry2D ()
 
 Geometry2D (const int coordim, CurveSharedPtr curve)
 
 ~Geometry2D () override
 
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 ResetNonRecursive (CurveMap &curvedEdges, CurveMap &curvedFaces)
 Reset this geometry object non-recursively: unset the current state, zero Geometry::m_coeffs and remove allocated GeomFactors. More...
 
void Setup ()
 
void GenGeomFactors ()
 Handles generation of geometry factors. More...
 

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

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...
 
void v_GenGeomFactors () override
 
void v_FillGeom () override
 
int v_GetDir (const int faceidx, const int facedir) const override
 Returns the element coordinate direction corresponding to a given face coordinate direction. More...
 
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...
 
void v_Setup () override
 
void PreSolveStraightEdge ()
 
- Protected Member Functions inherited from Nektar::SpatialDomains::Geometry2D
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 SolveStraightEdgeQuad (const Array< OneD, const NekDouble > &coords, Array< OneD, NekDouble > &Lcoords)
 
void v_CalculateInverseIsoParam () override
 
int v_AllLeftCheck (const Array< OneD, const NekDouble > &gloCoord) override
 
int v_GetShapeDim () const override
 Get the object's shape dimension. More...
 
PointGeomSharedPtr v_GetVertex (int i) const override
 
Geometry1DSharedPtr v_GetEdge (int i) const override
 Returns edge i of this object. More...
 
int v_GetNumVerts () const override
 Get the number of vertices of this object. More...
 
int v_GetNumEdges () const override
 Get the number of edges of this object. More...
 
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...
 
NekDouble v_FindDistance (const Array< OneD, const NekDouble > &xs, Array< OneD, NekDouble > &xi) override
 
- Protected Member Functions inherited from Nektar::SpatialDomains::Geometry
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
 
int m_straightEdge
 
- Static Protected Attributes inherited from Nektar::SpatialDomains::Geometry
static GeomFactorsVector m_regGeomFactorsManager
 

Detailed Description

Definition at line 54 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 86 of file QuadGeom.cpp.

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

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 54 of file QuadGeom.cpp.

56 : Geometry2D(edges[0]->GetVertex(0)->GetCoordim(), curve)
57{
58 int j;
59
61 m_globalID = id;
62
63 /// Copy the edge shared pointers.
64 m_edges.insert(m_edges.begin(), edges, edges + QuadGeom::kNedges);
65 m_eorient.resize(kNedges);
66
67 for (j = 0; j < kNedges; ++j)
68 {
69 m_eorient[j] =
70 SegGeom::GetEdgeOrientation(*edges[j], *edges[(j + 1) % kNedges]);
71 m_verts.push_back(
72 edges[j]->GetVertex(m_eorient[j] == StdRegions::eForwards ? 0 : 1));
73 }
74
75 for (j = 2; j < kNedges; ++j)
76 {
80 }
81
82 m_coordim = edges[0]->GetVertex(0)->GetCoordim();
83 ASSERTL0(m_coordim > 1, "Cannot call function with dim == 1");
84}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
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:189
static StdRegions::Orientation GetEdgeOrientation(const SegGeom &edge1, const SegGeom &edge2)
Get the orientation of edge1.
Definition: SegGeom.cpp:194

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 ( )
override

Definition at line 101 of file QuadGeom.cpp.

102{
103}

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 149 of file QuadGeom.cpp.

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

140{
141 return GetFaceOrientation(face1.m_verts, face2.m_verts, doRot, dir, angle,
142 tol);
143}
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:135

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

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

◆ PreSolveStraightEdge()

void Nektar::SpatialDomains::QuadGeom::PreSolveStraightEdge ( )
protected

Definition at line 464 of file QuadGeom.cpp.

465{
466 int i0, i1, j1, j2;
467 if (fabs(m_isoParameter[0][3]) >= fabs(m_isoParameter[1][3]))
468 {
469 i0 = 0;
470 i1 = 1;
471 }
472 else
473 {
474 i1 = 0;
475 i0 = 1;
476 m_straightEdge |= 2;
477 }
478 NekDouble gamma = m_isoParameter[i1][3] / m_isoParameter[i0][3];
479 std::vector<NekDouble> c(3);
480 for (int i = 0; i < 3; ++i)
481 {
482 c[i] = m_isoParameter[i1][i] - gamma * m_isoParameter[i0][i];
483 }
484 if (fabs(c[1]) >= fabs(c[2]))
485 {
486 j1 = 1;
487 j2 = 2;
488 }
489 else
490 {
491 j1 = 2;
492 j2 = 1;
493 m_straightEdge |= 4;
494 }
495 NekDouble beta = c[j2] / c[j1];
496 if (i0 == 1)
497 {
499 }
500 if (j1 == 2)
501 {
502 NekDouble temp = m_isoParameter[0][j1];
503 m_isoParameter[0][j1] = m_isoParameter[0][j2];
504 m_isoParameter[0][j2] = temp;
505 }
506 m_isoParameter[0][2] -= m_isoParameter[0][1] * beta;
507 m_isoParameter[1][0] = c[0];
508 m_isoParameter[1][1] = 1. / c[j1];
509 m_isoParameter[1][2] = beta;
510 m_isoParameter[1][3] = gamma;
511}
Array< OneD, Array< OneD, NekDouble > > m_isoParameter
Definition: Geometry.h:210
@ beta
Gauss Radau pinned at x=-1,.
Definition: PointsType.h:59
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825

References Nektar::LibUtilities::beta, Nektar::SpatialDomains::Geometry::m_isoParameter, Nektar::SpatialDomains::Geometry::m_straightEdge, and Vmath::Vcopy().

Referenced by v_GenGeomFactors().

◆ SetUpXmap()

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

Definition at line 105 of file QuadGeom.cpp.

106{
107 int order0 = max(m_edges[0]->GetXmap()->GetBasis(0)->GetNumModes(),
108 m_edges[2]->GetXmap()->GetBasis(0)->GetNumModes());
109 int order1 = max(m_edges[1]->GetXmap()->GetBasis(0)->GetNumModes(),
110 m_edges[3]->GetXmap()->GetBasis(0)->GetNumModes());
111
112 const LibUtilities::BasisKey B0(
114 LibUtilities::PointsKey(order0 + 1,
116 const LibUtilities::BasisKey B1(
118 LibUtilities::PointsKey(order1 + 1,
120
122}
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:195
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:51
@ eModified_A
Principle Modified Functions .
Definition: BasisType.h:48

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 389 of file QuadGeom.cpp.

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

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

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

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 513 of file QuadGeom.cpp.

514{
515 return i % 2;
516}

◆ 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 518 of file QuadGeom.cpp.

519{
520 Geometry::v_Reset(curvedEdges, curvedFaces);
521 CurveMap::iterator it = curvedFaces.find(m_globalID);
522
523 if (it != curvedFaces.end())
524 {
525 m_curve = it->second;
526 }
527
528 for (int i = 0; i < 4; ++i)
529 {
530 m_edges[i]->Reset(curvedEdges, curvedFaces);
531 }
532
533 SetUpXmap();
534 SetUpCoeffs(m_xmap->GetNcoeffs());
535}
void SetUpCoeffs(const int nCoeffs)
Initialise the Geometry::m_coeffs array.
Definition: Geometry.h:701
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:395

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 537 of file QuadGeom.cpp.

538{
539 if (!m_setupState)
540 {
541 for (int i = 0; i < 4; ++i)
542 {
543 m_edges[i]->Setup();
544 }
545 SetUpXmap();
546 SetUpCoeffs(m_xmap->GetNcoeffs());
547 m_setupState = true;
548 }
549}

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 75 of file QuadGeom.h.