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

3D geometry information More...

#include <Geometry3D.h>

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

Public Member Functions

 Geometry3D ()
 
 Geometry3D (const int coordim)
 
virtual ~Geometry3D ()
 
- 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 Attributes

static const int kDim = 3
 

Protected Member Functions

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, const Array< OneD, const NekDouble > &ptsz, Array< OneD, NekDouble > &Lcoords, NekDouble &dist)
 
void NewtonIterationForLocCoord (const Array< OneD, const NekDouble > &coords, Array< OneD, NekDouble > &Lcoords)
 
virtual void v_FillGeom () override
 Put all quadrature information into face/edge structure and backward transform. More...
 
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_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 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 int v_GetNumFaces () const override
 Get the number of faces of this object. 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 Geometry2DSharedPtr v_GetFace (int i) const override
 Returns face i 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 StdRegions::Orientation v_GetForient (const int i) const override
 Returns the orientation of face i with respect to the ordering of faces in the standard element. More...
 
- 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 ()
 

Protected Attributes

PointGeomVector m_verts
 
SegGeomVector m_edges
 
Geometry2DVector m_faces
 
std::vector< StdRegions::Orientationm_eorient
 
std::vector< StdRegions::Orientationm_forient
 
int m_eid
 
bool m_ownverts
 
- 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
 

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...
 
- Static Protected Attributes inherited from Nektar::SpatialDomains::Geometry
static GeomFactorsVector m_regGeomFactorsManager
 

Detailed Description

3D geometry information

Definition at line 67 of file Geometry3D.h.

Constructor & Destructor Documentation

◆ Geometry3D() [1/2]

Nektar::SpatialDomains::Geometry3D::Geometry3D ( )

Definition at line 51 of file Geometry3D.cpp.

52{
53}

◆ Geometry3D() [2/2]

Nektar::SpatialDomains::Geometry3D::Geometry3D ( const int  coordim)

Definition at line 55 of file Geometry3D.cpp.

55 : Geometry(coordim)
56{
58 "Coordinate dimension should be at least 3 for a 3D geometry.");
59}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
Geometry()
Default constructor.
Definition: Geometry.cpp:54
int m_coordim
Coordinate dimension of this geometry object.
Definition: Geometry.h:186

References ASSERTL0, and Nektar::SpatialDomains::Geometry::m_coordim.

◆ ~Geometry3D()

Nektar::SpatialDomains::Geometry3D::~Geometry3D ( )
virtual

Definition at line 61 of file Geometry3D.cpp.

62{
63}

Member Function Documentation

◆ NewtonIterationForLocCoord() [1/2]

void Nektar::SpatialDomains::Geometry3D::NewtonIterationForLocCoord ( const Array< OneD, const NekDouble > &  coords,
Array< OneD, NekDouble > &  Lcoords 
)
protected

Definition at line 72 of file Geometry3D.cpp.

74{
75 // maximum iterations for convergence
76 const int MaxIterations = 51;
77 // |x-xp|^2 < EPSILON error tolerance
78 const NekDouble Tol = 1.e-8;
79 // |r,s| > LcoordDIV stop the search
80 const NekDouble LcoordDiv = 15.0;
81
82 Array<OneD, const NekDouble> Jac =
83 m_geomFactors->GetJac(m_xmap->GetPointsKeys());
84
85 NekDouble ScaledTol =
86 Vmath::Vsum(Jac.size(), Jac, 1) / ((NekDouble)Jac.size());
87 ScaledTol *= Tol * Tol;
88
89 NekDouble xmap, ymap, zmap, res;
90 int cnt = 0;
91 Array<OneD, NekDouble> var(8, 1.);
92 Array<OneD, NekDouble> deriv(9);
93 Array<OneD, NekDouble> tmp(3);
94 while (cnt++ < MaxIterations)
95 {
96 var[1] = Lcoords[0];
97 var[2] = Lcoords[1];
98 var[3] = Lcoords[2];
99 var[4] = Lcoords[0] * Lcoords[1];
100 var[5] = Lcoords[1] * Lcoords[2];
101 var[6] = Lcoords[0] * Lcoords[2];
102 var[7] = var[4] * Lcoords[2];
103 // calculate the global point corresponding to Lcoords
104 xmap =
105 Vmath::Dot(m_isoParameter[0].size(), var, 1, m_isoParameter[0], 1);
106 ymap =
107 Vmath::Dot(m_isoParameter[0].size(), var, 1, m_isoParameter[1], 1);
108 zmap =
109 Vmath::Dot(m_isoParameter[0].size(), var, 1, m_isoParameter[2], 1);
110
111 tmp[0] = coords[0] - xmap;
112 tmp[1] = coords[1] - ymap;
113 tmp[2] = coords[2] - zmap;
114
115 res = tmp[0] * tmp[0] + tmp[1] * tmp[1] + tmp[2] * tmp[2];
116 if (res < ScaledTol)
117 {
118 break;
119 }
120
121 // Interpolate derivative metric at Lcoords (ddx1, ddx2, ddx3)
122 deriv[0] = m_isoParameter[0][1] + m_isoParameter[0][4] * Lcoords[1];
123 deriv[1] = m_isoParameter[0][2] + m_isoParameter[0][4] * Lcoords[0];
124 deriv[2] = m_isoParameter[0][3];
125 deriv[3] = m_isoParameter[1][1] + m_isoParameter[1][4] * Lcoords[1];
126 deriv[4] = m_isoParameter[1][2] + m_isoParameter[1][4] * Lcoords[0];
127 deriv[5] = m_isoParameter[1][3];
128 deriv[6] = m_isoParameter[2][1] + m_isoParameter[2][4] * Lcoords[1];
129 deriv[7] = m_isoParameter[2][2] + m_isoParameter[2][4] * Lcoords[0];
130 deriv[8] = m_isoParameter[2][3];
131 if (m_isoParameter[0].size() >= 6)
132 {
133 deriv[1] += m_isoParameter[0][5] * Lcoords[2];
134 deriv[2] += m_isoParameter[0][5] * Lcoords[1];
135 deriv[4] += m_isoParameter[1][5] * Lcoords[2];
136 deriv[5] += m_isoParameter[1][5] * Lcoords[1];
137 deriv[7] += m_isoParameter[2][5] * Lcoords[2];
138 deriv[8] += m_isoParameter[2][5] * Lcoords[1];
139 }
140 if (m_isoParameter[0].size() >= 8)
141 {
142 deriv[0] += m_isoParameter[0][6] * Lcoords[2] +
143 m_isoParameter[0][7] * var[5];
144 deriv[1] += m_isoParameter[0][7] * var[6];
145 deriv[2] += m_isoParameter[0][6] * Lcoords[0] +
146 m_isoParameter[0][7] * var[4];
147 deriv[3] += m_isoParameter[1][6] * Lcoords[2] +
148 m_isoParameter[1][7] * var[5];
149 deriv[4] += m_isoParameter[1][7] * var[6];
150 deriv[5] += m_isoParameter[1][6] * Lcoords[0] +
151 m_isoParameter[1][7] * var[4];
152 deriv[6] += m_isoParameter[2][6] * Lcoords[2] +
153 m_isoParameter[2][7] * var[5];
154 deriv[7] += m_isoParameter[2][7] * var[6];
155 deriv[8] += m_isoParameter[2][6] * Lcoords[0] +
156 m_isoParameter[2][7] * var[4];
157 }
158 DNekMatSharedPtr mat =
160 Vmath::Vcopy(9, deriv, 1, mat->GetPtr(), 1);
161 mat->Invert();
162 Lcoords[0] += Vmath::Dot(3, mat->GetPtr(), 1, tmp, 1);
163 Lcoords[1] += Vmath::Dot(3, mat->GetPtr() + 3, 1, tmp, 1);
164 Lcoords[2] += Vmath::Dot(3, mat->GetPtr() + 6, 1, tmp, 1);
165
166 if (!(std::isfinite(Lcoords[0]) && std::isfinite(Lcoords[1]) &&
167 std::isfinite(Lcoords[2])) ||
168 fabs(Lcoords[0]) > LcoordDiv || fabs(Lcoords[1]) > LcoordDiv ||
169 fabs(Lcoords[2]) > LcoordDiv)
170 {
171 std::ostringstream ss;
172 ss << "Iteration has diverged in NewtonIterationForLocCoord in "
173 "element "
174 << GetGlobalID();
175 WARNINGL1(false, ss.str());
176 return;
177 }
178 }
179
180 if (cnt >= MaxIterations)
181 {
182 std::ostringstream ss;
183
184 ss << "Reached MaxIterations (" << MaxIterations
185 << ") in Newton iteration ";
186
187 WARNINGL1(cnt < MaxIterations, ss.str());
188 }
189}
#define WARNINGL1(condition, msg)
Definition: ErrorUtil.hpp:250
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
Array< OneD, Array< OneD, NekDouble > > m_isoParameter
Definition: Geometry.h:207
int GetGlobalID(void) const
Get the ID of this object.
Definition: Geometry.h:327
StdRegions::StdExpansionSharedPtr m_xmap
mapping containing isoparametric transformation.
Definition: Geometry.h:192
GeomFactorsSharedPtr m_geomFactors
Geometric factors.
Definition: Geometry.h:188
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:75
double NekDouble
T Vsum(int n, const T *x, const int incx)
Subtract return sum(x)
Definition: Vmath.cpp:890
T Dot(int n, const T *w, const T *x)
dot (vector times vector): z = w*x
Definition: Vmath.cpp:1095
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1191

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Vmath::Dot(), Nektar::eFULL, Nektar::SpatialDomains::Geometry::GetGlobalID(), Nektar::SpatialDomains::Geometry::m_geomFactors, Nektar::SpatialDomains::Geometry::m_isoParameter, Nektar::SpatialDomains::Geometry::m_xmap, Vmath::Vcopy(), Vmath::Vsum(), and WARNINGL1.

◆ NewtonIterationForLocCoord() [2/2]

void Nektar::SpatialDomains::Geometry3D::NewtonIterationForLocCoord ( const Array< OneD, const NekDouble > &  coords,
const Array< OneD, const NekDouble > &  ptsx,
const Array< OneD, const NekDouble > &  ptsy,
const Array< OneD, const NekDouble > &  ptsz,
Array< OneD, NekDouble > &  Lcoords,
NekDouble dist 
)
protected

Definition at line 198 of file Geometry3D.cpp.

204{
205 // maximum iterations for convergence
206 const int MaxIterations = NekConstants::kNewtonIterations;
207 // |x-xp|^2 < EPSILON error tolerance
208 const NekDouble Tol = 1.e-8;
209 // |r,s| > LcoordDIV stop the search
210 const NekDouble LcoordDiv = 15.0;
211
212 Array<OneD, const NekDouble> Jac =
213 m_geomFactors->GetJac(m_xmap->GetPointsKeys());
214
215 NekDouble ScaledTol =
216 Vmath::Vsum(Jac.size(), Jac, 1) / ((NekDouble)Jac.size());
217 ScaledTol *= Tol;
218
219 NekDouble xmap, ymap, zmap, F1, F2, F3;
220
221 NekDouble derx_1, derx_2, derx_3, dery_1, dery_2, dery_3, derz_1, derz_2,
222 derz_3, jac;
223
224 // save intiial guess for later reference if required.
225 NekDouble init0 = Lcoords[0], init1 = Lcoords[1], init2 = Lcoords[2];
226
227 Array<OneD, NekDouble> DxD1(ptsx.size());
228 Array<OneD, NekDouble> DxD2(ptsx.size());
229 Array<OneD, NekDouble> DxD3(ptsx.size());
230 Array<OneD, NekDouble> DyD1(ptsx.size());
231 Array<OneD, NekDouble> DyD2(ptsx.size());
232 Array<OneD, NekDouble> DyD3(ptsx.size());
233 Array<OneD, NekDouble> DzD1(ptsx.size());
234 Array<OneD, NekDouble> DzD2(ptsx.size());
235 Array<OneD, NekDouble> DzD3(ptsx.size());
236
237 // Ideally this will be stored in m_geomfactors
238 m_xmap->PhysDeriv(ptsx, DxD1, DxD2, DxD3);
239 m_xmap->PhysDeriv(ptsy, DyD1, DyD2, DyD3);
240 m_xmap->PhysDeriv(ptsz, DzD1, DzD2, DzD3);
241
242 int cnt = 0;
243 Array<OneD, DNekMatSharedPtr> I(3);
244 Array<OneD, NekDouble> eta(3);
245
246 F1 = F2 = F3 = 2000; // Starting value of Function
247 NekDouble resid = sqrt(F1 * F1 + F2 * F2 + F3 * F3);
248 while (cnt++ < MaxIterations)
249 {
250 // evaluate lagrange interpolant at Lcoords
251 m_xmap->LocCoordToLocCollapsed(Lcoords, eta);
252 I[0] = m_xmap->GetBasis(0)->GetI(eta);
253 I[1] = m_xmap->GetBasis(1)->GetI(eta + 1);
254 I[2] = m_xmap->GetBasis(2)->GetI(eta + 2);
255
256 // calculate the global point `corresponding to Lcoords
257 xmap = m_xmap->PhysEvaluate(I, ptsx);
258 ymap = m_xmap->PhysEvaluate(I, ptsy);
259 zmap = m_xmap->PhysEvaluate(I, ptsz);
260
261 F1 = coords[0] - xmap;
262 F2 = coords[1] - ymap;
263 F3 = coords[2] - zmap;
264
265 if (F1 * F1 + F2 * F2 + F3 * F3 < ScaledTol)
266 {
267 resid = sqrt(F1 * F1 + F2 * F2 + F3 * F3);
268 break;
269 }
270
271 // Interpolate derivative metric at Lcoords
272 derx_1 = m_xmap->PhysEvaluate(I, DxD1);
273 derx_2 = m_xmap->PhysEvaluate(I, DxD2);
274 derx_3 = m_xmap->PhysEvaluate(I, DxD3);
275 dery_1 = m_xmap->PhysEvaluate(I, DyD1);
276 dery_2 = m_xmap->PhysEvaluate(I, DyD2);
277 dery_3 = m_xmap->PhysEvaluate(I, DyD3);
278 derz_1 = m_xmap->PhysEvaluate(I, DzD1);
279 derz_2 = m_xmap->PhysEvaluate(I, DzD2);
280 derz_3 = m_xmap->PhysEvaluate(I, DzD3);
281
282 jac = derx_1 * (dery_2 * derz_3 - dery_3 * derz_2) -
283 derx_2 * (dery_1 * derz_3 - dery_3 * derz_1) +
284 derx_3 * (dery_1 * derz_2 - dery_2 * derz_1);
285
286 // use analytical inverse of derivitives which are also similar to
287 // those of metric factors.
288 Lcoords[0] =
289 Lcoords[0] +
290 ((dery_2 * derz_3 - dery_3 * derz_2) * (coords[0] - xmap) -
291 (derx_2 * derz_3 - derx_3 * derz_2) * (coords[1] - ymap) +
292 (derx_2 * dery_3 - derx_3 * dery_2) * (coords[2] - zmap)) /
293 jac;
294
295 Lcoords[1] =
296 Lcoords[1] -
297 ((dery_1 * derz_3 - dery_3 * derz_1) * (coords[0] - xmap) -
298 (derx_1 * derz_3 - derx_3 * derz_1) * (coords[1] - ymap) +
299 (derx_1 * dery_3 - derx_3 * dery_1) * (coords[2] - zmap)) /
300 jac;
301
302 Lcoords[2] =
303 Lcoords[2] +
304 ((dery_1 * derz_2 - dery_2 * derz_1) * (coords[0] - xmap) -
305 (derx_1 * derz_2 - derx_2 * derz_1) * (coords[1] - ymap) +
306 (derx_1 * dery_2 - derx_2 * dery_1) * (coords[2] - zmap)) /
307 jac;
308
309 if (!(std::isfinite(Lcoords[0]) && std::isfinite(Lcoords[1]) &&
310 std::isfinite(Lcoords[2])))
311 {
312 dist = 1e16;
313 std::ostringstream ss;
314 ss << "nan or inf found in NewtonIterationForLocCoord in element "
315 << GetGlobalID();
316 WARNINGL1(false, ss.str());
317 return;
318 }
319 if (fabs(Lcoords[0]) > LcoordDiv || fabs(Lcoords[1]) > LcoordDiv ||
320 fabs(Lcoords[2]) > LcoordDiv)
321 {
322 break; // lcoords have diverged so stop iteration
323 }
324 }
325
326 m_xmap->LocCoordToLocCollapsed(Lcoords, eta);
327 if (ClampLocCoords(eta, 0.))
328 {
329 I[0] = m_xmap->GetBasis(0)->GetI(eta);
330 I[1] = m_xmap->GetBasis(1)->GetI(eta + 1);
331 I[2] = m_xmap->GetBasis(2)->GetI(eta + 2);
332 // calculate the global point corresponding to Lcoords
333 xmap = m_xmap->PhysEvaluate(I, ptsx);
334 ymap = m_xmap->PhysEvaluate(I, ptsy);
335 zmap = m_xmap->PhysEvaluate(I, ptsz);
336 F1 = coords[0] - xmap;
337 F2 = coords[1] - ymap;
338 F3 = coords[2] - zmap;
339 dist = sqrt(F1 * F1 + F2 * F2 + F3 * F3);
340 }
341 else
342 {
343 dist = 0.;
344 }
345
346 if (cnt >= MaxIterations)
347 {
348 Array<OneD, NekDouble> collCoords(3);
349 m_xmap->LocCoordToLocCollapsed(Lcoords, collCoords);
350
351 // if coordinate is inside element dump error!
352 if ((collCoords[0] >= -1.0 && collCoords[0] <= 1.0) &&
353 (collCoords[1] >= -1.0 && collCoords[1] <= 1.0) &&
354 (collCoords[2] >= -1.0 && collCoords[2] <= 1.0))
355 {
356 std::ostringstream ss;
357
358 ss << "Reached MaxIterations (" << MaxIterations
359 << ") in Newton iteration ";
360 ss << "Init value (" << setprecision(4) << init0 << "," << init1
361 << "," << init2 << ") ";
362 ss << "Fin value (" << Lcoords[0] << "," << Lcoords[1] << ","
363 << Lcoords[2] << ") ";
364 ss << "Resid = " << resid << " Tolerance = " << sqrt(ScaledTol);
365
366 WARNINGL1(cnt < MaxIterations, ss.str());
367 }
368 }
369}
bool ClampLocCoords(Array< OneD, NekDouble > &locCoord, NekDouble tol=std::numeric_limits< NekDouble >::epsilon())
Clamp local coords to be within standard regions [-1, 1]^dim.
Definition: Geometry.cpp:556
static const unsigned int kNewtonIterations
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294

References Nektar::SpatialDomains::Geometry::ClampLocCoords(), Nektar::SpatialDomains::Geometry::GetGlobalID(), Nektar::NekConstants::kNewtonIterations, Nektar::SpatialDomains::Geometry::m_geomFactors, Nektar::SpatialDomains::Geometry::m_xmap, tinysimd::sqrt(), Vmath::Vsum(), and WARNINGL1.

Referenced by v_GetLocCoords().

◆ v_AllLeftCheck()

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

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 191 of file Geometry3D.cpp.

192{
193 boost::ignore_unused(gloCoord);
194 // todo: check only plane surfaces
195 return 0;
196}

◆ v_CalculateInverseIsoParam()

void Nektar::SpatialDomains::Geometry3D::v_CalculateInverseIsoParam ( )
overrideprotectedvirtual

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 572 of file Geometry3D.cpp.

573{
574 DNekMatSharedPtr mat =
576 for (int i = 0; i < 3; ++i)
577 {
578 for (int j = 0; j < 3; ++j)
579 {
580
581 mat->SetValue(i, j, m_isoParameter[i][j + 1]);
582 }
583 }
584 mat->Invert();
585 m_invIsoParam = Array<OneD, Array<OneD, NekDouble>>(3);
586 for (int i = 0; i < 3; ++i)
587 {
588 m_invIsoParam[i] = Array<OneD, NekDouble>(3);
589 for (int j = 0; j < 3; ++j)
590 {
591 m_invIsoParam[i][j] = mat->GetValue(i, j);
592 }
593 }
594}
Array< OneD, Array< OneD, NekDouble > > m_invIsoParam
Definition: Geometry.h:208

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::eFULL, Nektar::SpatialDomains::Geometry::m_invIsoParam, and Nektar::SpatialDomains::Geometry::m_isoParameter.

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

◆ v_FillGeom()

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

Put all quadrature information into face/edge structure and backward transform.

Note verts, edges, and faces 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 449 of file Geometry3D.cpp.

450{
451 if (m_state == ePtsFilled)
452 {
453 return;
454 }
455
456 int i, j, k;
457
458 for (i = 0; i < m_forient.size(); i++)
459 {
460 m_faces[i]->FillGeom();
461
462 int nFaceCoeffs = m_faces[i]->GetXmap()->GetNcoeffs();
463
464 Array<OneD, unsigned int> mapArray(nFaceCoeffs);
465 Array<OneD, int> signArray(nFaceCoeffs);
466
467 if (m_forient[i] < 9)
468 {
469 m_xmap->GetTraceToElementMap(
470 i, mapArray, signArray, m_forient[i],
471 m_faces[i]->GetXmap()->GetTraceNcoeffs(0),
472 m_faces[i]->GetXmap()->GetTraceNcoeffs(1));
473 }
474 else
475 {
476 m_xmap->GetTraceToElementMap(
477 i, mapArray, signArray, m_forient[i],
478 m_faces[i]->GetXmap()->GetTraceNcoeffs(1),
479 m_faces[i]->GetXmap()->GetTraceNcoeffs(0));
480 }
481
482 for (j = 0; j < m_coordim; j++)
483 {
484 const Array<OneD, const NekDouble> &coeffs =
485 m_faces[i]->GetCoeffs(j);
486
487 for (k = 0; k < nFaceCoeffs; k++)
488 {
489 NekDouble v = signArray[k] * coeffs[k];
490 m_coeffs[j][mapArray[k]] = v;
491 }
492 }
493 }
494
496}
std::vector< StdRegions::Orientation > m_forient
Definition: Geometry3D.h:84
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
StdRegions::StdExpansionSharedPtr GetXmap() const
Return the mapping object Geometry::m_xmap that represents the coordinate transformation from standar...
Definition: Geometry.h:436
@ ePtsFilled
Geometric information has been generated.

References Nektar::SpatialDomains::ePtsFilled, Nektar::SpatialDomains::Geometry::GetXmap(), Nektar::SpatialDomains::Geometry::m_coeffs, Nektar::SpatialDomains::Geometry::m_coordim, m_faces, m_forient, Nektar::SpatialDomains::Geometry::m_state, and Nektar::SpatialDomains::Geometry::m_xmap.

Referenced by Nektar::SpatialDomains::PyrGeom::v_GenGeomFactors(), Nektar::SpatialDomains::HexGeom::v_GenGeomFactors(), Nektar::SpatialDomains::PrismGeom::v_GenGeomFactors(), Nektar::SpatialDomains::TetGeom::v_GenGeomFactors(), and v_GetLocCoords().

◆ v_GetCoord()

NekDouble Nektar::SpatialDomains::Geometry3D::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 502 of file Geometry3D.cpp.

504{
505 ASSERTL1(m_state == ePtsFilled, "Geometry is not in physical space");
506
507 Array<OneD, NekDouble> tmp(m_xmap->GetTotPoints());
508 m_xmap->BwdTrans(m_coeffs[i], tmp);
509
510 return m_xmap->PhysEvaluate(Lcoord, tmp);
511}
#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_GetEdge()

Geometry1DSharedPtr Nektar::SpatialDomains::Geometry3D::v_GetEdge ( int  i) const
overrideprotectedvirtual

Returns edge i of this object.

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 542 of file Geometry3D.cpp.

543{
544 ASSERTL2(i >= 0 && i <= m_edges.size() - 1,
545 "Edge ID must be between 0 and " +
546 boost::lexical_cast<string>(m_edges.size() - 1));
547 return m_edges[i];
548}
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed to...
Definition: ErrorUtil.hpp:272

References ASSERTL2, and m_edges.

◆ v_GetEorient()

StdRegions::Orientation Nektar::SpatialDomains::Geometry3D::v_GetEorient ( const int  i) const
inlineoverrideprotectedvirtual

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

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 556 of file Geometry3D.cpp.

557{
558 ASSERTL2(i >= 0 && i <= m_edges.size() - 1,
559 "Edge ID must be between 0 and " +
560 boost::lexical_cast<string>(m_edges.size() - 1));
561 return m_eorient[i];
562}
std::vector< StdRegions::Orientation > m_eorient
Definition: Geometry3D.h:83

References ASSERTL2, m_edges, and m_eorient.

◆ v_GetFace()

Geometry2DSharedPtr Nektar::SpatialDomains::Geometry3D::v_GetFace ( int  i) const
overrideprotectedvirtual

Returns face i of this object.

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 550 of file Geometry3D.cpp.

551{
552 ASSERTL2((i >= 0) && (i <= 5), "Edge id must be between 0 and 4");
553 return m_faces[i];
554}

References ASSERTL2, and m_faces.

◆ v_GetForient()

StdRegions::Orientation Nektar::SpatialDomains::Geometry3D::v_GetForient ( const int  i) const
overrideprotectedvirtual

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

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 564 of file Geometry3D.cpp.

565{
566 ASSERTL2(i >= 0 && i <= m_faces.size() - 1,
567 "Face ID must be between 0 and " +
568 boost::lexical_cast<string>(m_faces.size() - 1));
569 return m_forient[i];
570}

References ASSERTL2, m_faces, and m_forient.

◆ v_GetLocCoords()

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

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 from Nektar::SpatialDomains::Geometry.

Definition at line 371 of file Geometry3D.cpp.

373{
374 NekDouble dist = std::numeric_limits<double>::max();
375 Array<OneD, NekDouble> tmpcoords(3);
376 tmpcoords[0] = coords[0];
377 tmpcoords[1] = coords[1];
378 tmpcoords[2] = coords[2];
379 if (GetMetricInfo()->GetGtype() == eRegular)
380 {
381 tmpcoords[0] -= m_isoParameter[0][0];
382 tmpcoords[1] -= m_isoParameter[1][0];
383 tmpcoords[2] -= m_isoParameter[2][0];
384 Lcoords[0] = m_invIsoParam[0][0] * tmpcoords[0] +
385 m_invIsoParam[0][1] * tmpcoords[1] +
386 m_invIsoParam[0][2] * tmpcoords[2];
387 Lcoords[1] = m_invIsoParam[1][0] * tmpcoords[0] +
388 m_invIsoParam[1][1] * tmpcoords[1] +
389 m_invIsoParam[1][2] * tmpcoords[2];
390 Lcoords[2] = m_invIsoParam[2][0] * tmpcoords[0] +
391 m_invIsoParam[2][1] * tmpcoords[1] +
392 m_invIsoParam[2][2] * tmpcoords[2];
393 }
394 else if (m_straightEdge)
395 {
396 ClampLocCoords(Lcoords, 0.);
397 NewtonIterationForLocCoord(tmpcoords, Lcoords);
398 }
399 else
400 {
401 v_FillGeom();
402 // Determine nearest point of coords to values in m_xmap
403 int npts = m_xmap->GetTotPoints();
404 Array<OneD, NekDouble> ptsx(npts), ptsy(npts), ptsz(npts);
405 Array<OneD, NekDouble> tmp1(npts), tmp2(npts);
406
407 m_xmap->BwdTrans(m_coeffs[0], ptsx);
408 m_xmap->BwdTrans(m_coeffs[1], ptsy);
409 m_xmap->BwdTrans(m_coeffs[2], ptsz);
410
411 const Array<OneD, const NekDouble> za = m_xmap->GetPoints(0);
412 const Array<OneD, const NekDouble> zb = m_xmap->GetPoints(1);
413 const Array<OneD, const NekDouble> zc = m_xmap->GetPoints(2);
414
415 // guess the first local coords based on nearest point
416 Vmath::Sadd(npts, -coords[0], ptsx, 1, tmp1, 1);
417 Vmath::Vmul(npts, tmp1, 1, tmp1, 1, tmp1, 1);
418 Vmath::Sadd(npts, -coords[1], ptsy, 1, tmp2, 1);
419 Vmath::Vvtvp(npts, tmp2, 1, tmp2, 1, tmp1, 1, tmp1, 1);
420 Vmath::Sadd(npts, -coords[2], ptsz, 1, tmp2, 1);
421 Vmath::Vvtvp(npts, tmp2, 1, tmp2, 1, tmp1, 1, tmp1, 1);
422
423 int min_i = Vmath::Imin(npts, tmp1, 1);
424
425 // Get Local coordinates
426 int qa = za.size(), qb = zb.size();
427
428 Array<OneD, NekDouble> eta(3, 0.);
429 eta[2] = zc[min_i / (qa * qb)];
430 min_i = min_i % (qa * qb);
431 eta[1] = zb[min_i / qa];
432 eta[0] = za[min_i % qa];
433 m_xmap->LocCollapsedToLocCoord(eta, Lcoords);
434
435 // Perform newton iteration to find local coordinates
436 NewtonIterationForLocCoord(coords, ptsx, ptsy, ptsz, Lcoords, dist);
437 }
438 return dist;
439}
virtual void v_FillGeom() override
Put all quadrature information into face/edge structure and backward transform.
Definition: Geometry3D.cpp:449
void NewtonIterationForLocCoord(const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &ptsx, const Array< OneD, const NekDouble > &ptsy, const Array< OneD, const NekDouble > &ptsz, Array< OneD, NekDouble > &Lcoords, NekDouble &dist)
Definition: Geometry3D.cpp:198
GeomFactorsSharedPtr GetMetricInfo()
Get the geometric factors for this object.
Definition: Geometry.h:311
@ eRegular
Geometry is straight-sided with constant geometric factors.
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:207
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.cpp:569
int Imin(int n, const T *x, const int incx)
Return the index of the minimum element in x.
Definition: Vmath.cpp:1018
void Sadd(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Add scalar y = alpha + x.
Definition: Vmath.cpp:379

References Nektar::SpatialDomains::Geometry::ClampLocCoords(), Nektar::SpatialDomains::eRegular, Nektar::SpatialDomains::Geometry::GetMetricInfo(), Vmath::Imin(), Nektar::SpatialDomains::Geometry::m_coeffs, Nektar::SpatialDomains::Geometry::m_invIsoParam, Nektar::SpatialDomains::Geometry::m_isoParameter, Nektar::SpatialDomains::Geometry::m_straightEdge, Nektar::SpatialDomains::Geometry::m_xmap, NewtonIterationForLocCoord(), Vmath::Sadd(), v_FillGeom(), Vmath::Vmul(), and Vmath::Vvtvp().

◆ v_GetNumEdges()

int Nektar::SpatialDomains::Geometry3D::v_GetNumEdges ( ) const
overrideprotectedvirtual

Get the number of edges of this object.

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 527 of file Geometry3D.cpp.

528{
529 return m_edges.size();
530}

References m_edges.

◆ v_GetNumFaces()

int Nektar::SpatialDomains::Geometry3D::v_GetNumFaces ( ) const
overrideprotectedvirtual

Get the number of faces of this object.

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 532 of file Geometry3D.cpp.

533{
534 return m_faces.size();
535}

References m_faces.

◆ v_GetNumVerts()

int Nektar::SpatialDomains::Geometry3D::v_GetNumVerts ( ) const
overrideprotectedvirtual

Get the number of vertices of this object.

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 522 of file Geometry3D.cpp.

523{
524 return m_verts.size();
525}

References m_verts.

◆ v_GetShapeDim()

int Nektar::SpatialDomains::Geometry3D::v_GetShapeDim ( ) const
overrideprotectedvirtual

Get the object's shape dimension.

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

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 517 of file Geometry3D.cpp.

518{
519 return 3;
520}

◆ v_GetVertex()

PointGeomSharedPtr Nektar::SpatialDomains::Geometry3D::v_GetVertex ( int  i) const
overrideprotectedvirtual

Implements Nektar::SpatialDomains::Geometry.

Definition at line 537 of file Geometry3D.cpp.

538{
539 return m_verts[i];
540}

References m_verts.

Member Data Documentation

◆ kDim

const int Nektar::SpatialDomains::Geometry3D::kDim = 3
static

Definition at line 77 of file Geometry3D.h.

◆ m_edges

SegGeomVector Nektar::SpatialDomains::Geometry3D::m_edges
protected

◆ m_eid

int Nektar::SpatialDomains::Geometry3D::m_eid
protected

Definition at line 85 of file Geometry3D.h.

◆ m_eorient

std::vector<StdRegions::Orientation> Nektar::SpatialDomains::Geometry3D::m_eorient
protected

◆ m_faces

Geometry2DVector Nektar::SpatialDomains::Geometry3D::m_faces
protected

◆ m_forient

std::vector<StdRegions::Orientation> Nektar::SpatialDomains::Geometry3D::m_forient
protected

◆ m_ownverts

bool Nektar::SpatialDomains::Geometry3D::m_ownverts
protected

Definition at line 86 of file Geometry3D.h.

◆ m_verts

PointGeomVector Nektar::SpatialDomains::Geometry3D::m_verts
protected