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

3D geometry information More...

#include <Geometry3D.h>

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

Public Member Functions

 Geometry3D ()
 
 Geometry3D (const int coordim)
 
virtual ~Geometry3D ()
 
int GetDir (const int faceidx, const int facedir) const
 Returns the element coordinate direction corresponding to a given face coordinate direction. More...
 
- 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...
 
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 &resid)
 Determine whether an element contains a particular Cartesian coordinate \(\vec{x} = (x,y,z)\). More...
 
NekDouble GetLocCoords (const Array< OneD, const NekDouble > &coords, Array< OneD, NekDouble > &Lcoords)
 Determine the local collapsed coordinates that correspond to a given Cartesian coordinate for this geometry object. More...
 
NekDouble GetCoord (const int i, const Array< OneD, const NekDouble > &Lcoord)
 Given local collapsed coordinate Lcoord, return the value of physical coordinate in direction i. More...
 
bool MinMaxCheck (const Array< OneD, const NekDouble > &gloCoord)
 Check if given global coord is within twice the min/max distance of the element. More...
 
void ClampLocCoords (Array< OneD, NekDouble > &locCoord, NekDouble tol)
 Clamp local coords to be within standard regions [-1, 1]^dim. More...
 
int GetVertexEdgeMap (int i, int j) const
 Returns the standard element edge IDs that are connected to a given vertex. More...
 
int GetVertexFaceMap (int i, int j) const
 Returns the standard element face IDs that are connected to a given vertex. More...
 
int GetEdgeFaceMap (int i, int j) const
 Returns the standard element edge IDs that are connected to a given face. More...
 
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

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 &resid)
 
virtual void v_FillGeom ()
 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)
 Given local collapsed coordinate Lcoord return the value of physical coordinate in direction i. More...
 
virtual int v_GetShapeDim () const
 Get the object's shape dimension. 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 PointGeomSharedPtr v_GetVertex (int i) const
 
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_GetDir (const int faceidx, const int facedir) const =0
 
- Protected Member Functions inherited from Nektar::SpatialDomains::Geometry
void GenGeomFactors ()
 Handles generation of geometry factors. 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 bool v_ContainsPoint (const Array< OneD, const NekDouble > &gloCoord, Array< OneD, NekDouble > &locCoord, NekDouble tol, NekDouble &resid)
 
virtual NekDouble v_GetLocCoords (const Array< OneD, const NekDouble > &coords, Array< OneD, NekDouble > &Lcoords)
 Determine the local collapsed coordinates that correspond to a given Cartesian coordinate for this geometry object. More...
 
virtual int v_GetVertexEdgeMap (int i, int j) const
 Returns the standard element edge IDs that are connected to a given vertex. More...
 
virtual int v_GetVertexFaceMap (int i, int j) const
 Returns the standard element face IDs that are connected to a given vertex. More...
 
virtual int v_GetEdgeFaceMap (int i, int j) const
 Returns the standard element edge IDs that are connected to a given face. More...
 
virtual 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...
 

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

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.

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

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

◆ ~Geometry3D()

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

Definition at line 61 of file Geometry3D.cpp.

62 {
63 }

Member Function Documentation

◆ GetDir()

int Nektar::SpatialDomains::Geometry3D::GetDir ( const int  faceidx,
const int  facedir 
) const

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

Definition at line 73 of file Geometry3D.cpp.

References v_GetDir().

Referenced by Nektar::SpatialDomains::MeshGraph::GetFaceBasisKey().

74 {
75  return v_GetDir(faceidx, facedir);
76 }
virtual int v_GetDir(const int faceidx, const int facedir) const =0

◆ NewtonIterationForLocCoord()

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 resid 
)
protected

Definition at line 81 of file Geometry3D.cpp.

References Nektar::SpatialDomains::Geometry::m_geomFactors, Nektar::SpatialDomains::Geometry::m_xmap, Vmath::Vsum(), and WARNINGL1.

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

88 {
89  // maximum iterations for convergence
90  const int MaxIterations = 51;
91  // |x-xp|^2 < EPSILON error tolerance
92  const NekDouble Tol = 1.e-8;
93  // |r,s| > LcoordDIV stop the search
94  const NekDouble LcoordDiv = 15.0;
95 
96  Array<OneD, const NekDouble> Jac =
97  m_geomFactors->GetJac(m_xmap->GetPointsKeys());
98 
99  NekDouble ScaledTol = Vmath::Vsum(Jac.num_elements(), Jac, 1) /
100  ((NekDouble)Jac.num_elements());
101  ScaledTol *= Tol;
102 
103  NekDouble xmap, ymap, zmap, F1, F2, F3;
104 
105  NekDouble derx_1, derx_2, derx_3, dery_1, dery_2, dery_3, derz_1, derz_2,
106  derz_3, jac;
107 
108  // save intiial guess for later reference if required.
109  NekDouble init0 = Lcoords[0], init1 = Lcoords[1], init2 = Lcoords[2];
110 
111  Array<OneD, NekDouble> DxD1(ptsx.num_elements());
112  Array<OneD, NekDouble> DxD2(ptsx.num_elements());
113  Array<OneD, NekDouble> DxD3(ptsx.num_elements());
114  Array<OneD, NekDouble> DyD1(ptsx.num_elements());
115  Array<OneD, NekDouble> DyD2(ptsx.num_elements());
116  Array<OneD, NekDouble> DyD3(ptsx.num_elements());
117  Array<OneD, NekDouble> DzD1(ptsx.num_elements());
118  Array<OneD, NekDouble> DzD2(ptsx.num_elements());
119  Array<OneD, NekDouble> DzD3(ptsx.num_elements());
120 
121  // Ideally this will be stored in m_geomfactors
122  m_xmap->PhysDeriv(ptsx, DxD1, DxD2, DxD3);
123  m_xmap->PhysDeriv(ptsy, DyD1, DyD2, DyD3);
124  m_xmap->PhysDeriv(ptsz, DzD1, DzD2, DzD3);
125 
126  int cnt = 0;
127  Array<OneD, DNekMatSharedPtr> I(3);
128  Array<OneD, NekDouble> eta(3);
129 
130  F1 = F2 = F3 = 2000; // Starting value of Function
131 
132  while (cnt++ < MaxIterations)
133  {
134  // evaluate lagrange interpolant at Lcoords
135  m_xmap->LocCoordToLocCollapsed(Lcoords, eta);
136  I[0] = m_xmap->GetBasis(0)->GetI(eta);
137  I[1] = m_xmap->GetBasis(1)->GetI(eta + 1);
138  I[2] = m_xmap->GetBasis(2)->GetI(eta + 2);
139 
140  // calculate the global point `corresponding to Lcoords
141  xmap = m_xmap->PhysEvaluate(I, ptsx);
142  ymap = m_xmap->PhysEvaluate(I, ptsy);
143  zmap = m_xmap->PhysEvaluate(I, ptsz);
144 
145  F1 = coords[0] - xmap;
146  F2 = coords[1] - ymap;
147  F3 = coords[2] - zmap;
148 
149  if (F1 * F1 + F2 * F2 + F3 * F3 < ScaledTol)
150  {
151  resid = sqrt(F1 * F1 + F2 * F2 + F3 * F3);
152  break;
153  }
154 
155  // Interpolate derivative metric at Lcoords
156  derx_1 = m_xmap->PhysEvaluate(I, DxD1);
157  derx_2 = m_xmap->PhysEvaluate(I, DxD2);
158  derx_3 = m_xmap->PhysEvaluate(I, DxD3);
159  dery_1 = m_xmap->PhysEvaluate(I, DyD1);
160  dery_2 = m_xmap->PhysEvaluate(I, DyD2);
161  dery_3 = m_xmap->PhysEvaluate(I, DyD3);
162  derz_1 = m_xmap->PhysEvaluate(I, DzD1);
163  derz_2 = m_xmap->PhysEvaluate(I, DzD2);
164  derz_3 = m_xmap->PhysEvaluate(I, DzD3);
165 
166  jac = derx_1 * (dery_2 * derz_3 - dery_3 * derz_2) -
167  derx_2 * (dery_1 * derz_3 - dery_3 * derz_1) +
168  derx_3 * (dery_1 * derz_2 - dery_2 * derz_1);
169 
170  // use analytical inverse of derivitives which are also similar to
171  // those of metric factors.
172  Lcoords[0] =
173  Lcoords[0] +
174  ((dery_2 * derz_3 - dery_3 * derz_2) * (coords[0] - xmap) -
175  (derx_2 * derz_3 - derx_3 * derz_2) * (coords[1] - ymap) +
176  (derx_2 * dery_3 - derx_3 * dery_2) * (coords[2] - zmap)) /
177  jac;
178 
179  Lcoords[1] =
180  Lcoords[1] -
181  ((dery_1 * derz_3 - dery_3 * derz_1) * (coords[0] - xmap) -
182  (derx_1 * derz_3 - derx_3 * derz_1) * (coords[1] - ymap) +
183  (derx_1 * dery_3 - derx_3 * dery_1) * (coords[2] - zmap)) /
184  jac;
185 
186  Lcoords[2] =
187  Lcoords[2] +
188  ((dery_1 * derz_2 - dery_2 * derz_1) * (coords[0] - xmap) -
189  (derx_1 * derz_2 - derx_2 * derz_1) * (coords[1] - ymap) +
190  (derx_1 * dery_2 - derx_2 * dery_1) * (coords[2] - zmap)) /
191  jac;
192 
193  if (fabs(Lcoords[0]) > LcoordDiv || fabs(Lcoords[1]) > LcoordDiv ||
194  fabs(Lcoords[0]) > LcoordDiv)
195  {
196  break; // lcoords have diverged so stop iteration
197  }
198  }
199 
200  resid = sqrt(F1 * F1 + F2 * F2 + F3 * F3);
201 
202  if (cnt >= MaxIterations)
203  {
204  Array<OneD, NekDouble> collCoords(3);
205  m_xmap->LocCoordToLocCollapsed(Lcoords, collCoords);
206 
207  // if coordinate is inside element dump error!
208  if ((collCoords[0] >= -1.0 && collCoords[0] <= 1.0) &&
209  (collCoords[1] >= -1.0 && collCoords[1] <= 1.0) &&
210  (collCoords[2] >= -1.0 && collCoords[2] <= 1.0))
211  {
212  std::ostringstream ss;
213 
214  ss << "Reached MaxIterations (" << MaxIterations
215  << ") in Newton iteration ";
216  ss << "Init value (" << setprecision(4) << init0 << "," << init1
217  << "," << init2 << ") ";
218  ss << "Fin value (" << Lcoords[0] << "," << Lcoords[1] << ","
219  << Lcoords[2] << ") ";
220  ss << "Resid = " << resid << " Tolerance = " << sqrt(ScaledTol);
221 
222  WARNINGL1(cnt < MaxIterations, ss.str());
223  }
224  }
225 }
StdRegions::StdExpansionSharedPtr m_xmap
mapping containing isoparametric transformation.
Definition: Geometry.h:189
GeomFactorsSharedPtr m_geomFactors
Geometric factors.
Definition: Geometry.h:185
double NekDouble
#define WARNINGL1(condition, msg)
Definition: ErrorUtil.hpp:251
T Vsum(int n, const T *x, const int incx)
Subtract return sum(x)
Definition: Vmath.cpp:740

◆ v_FillGeom()

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

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 235 of file Geometry3D.cpp.

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

236 {
237  if (m_state == ePtsFilled)
238  {
239  return;
240  }
241 
242  int i, j, k;
243 
244  for (i = 0; i < m_forient.size(); i++)
245  {
246  m_faces[i]->FillGeom();
247 
248  int nFaceCoeffs = m_faces[i]->GetXmap()->GetNcoeffs();
249 
250  Array<OneD, unsigned int> mapArray(nFaceCoeffs);
251  Array<OneD, int> signArray(nFaceCoeffs);
252 
253  if (m_forient[i] < 9)
254  {
255  m_xmap->GetFaceToElementMap(
256  i,
257  m_forient[i],
258  mapArray,
259  signArray,
260  m_faces[i]->GetXmap()->GetEdgeNcoeffs(0),
261  m_faces[i]->GetXmap()->GetEdgeNcoeffs(1));
262  }
263  else
264  {
265  m_xmap->GetFaceToElementMap(
266  i,
267  m_forient[i],
268  mapArray,
269  signArray,
270  m_faces[i]->GetXmap()->GetEdgeNcoeffs(1),
271  m_faces[i]->GetXmap()->GetEdgeNcoeffs(0));
272  }
273 
274  for (j = 0; j < m_coordim; j++)
275  {
276  const Array<OneD, const NekDouble> &coeffs =
277  m_faces[i]->GetCoeffs(j);
278 
279  for (k = 0; k < nFaceCoeffs; k++)
280  {
281  NekDouble v = signArray[k] * coeffs[k];
282  m_coeffs[j][mapArray[k]] = v;
283  }
284  }
285  }
286 
288 }
StdRegions::StdExpansionSharedPtr m_xmap
mapping containing isoparametric transformation.
Definition: Geometry.h:189
StdRegions::StdExpansionSharedPtr GetXmap() const
Return the mapping object Geometry::m_xmap that represents the coordinate transformation from standar...
Definition: Geometry.h:421
std::vector< StdRegions::Orientation > m_forient
Definition: Geometry3D.h:87
double NekDouble
Array< OneD, Array< OneD, NekDouble > > m_coeffs
Array containing expansion coefficients of m_xmap.
Definition: Geometry.h:201
Geometric information has been generated.
GeomState m_state
Enumeration to dictate whether coefficients are filled.
Definition: Geometry.h:191
int m_coordim
Coordinate dimension of this geometry object.
Definition: Geometry.h:183

◆ v_GetCoord()

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

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

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 294 of file Geometry3D.cpp.

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

296 {
297  ASSERTL1(m_state == ePtsFilled, "Geometry is not in physical space");
298 
299  Array<OneD, NekDouble> tmp(m_xmap->GetTotPoints());
300  m_xmap->BwdTrans(m_coeffs[i], tmp);
301 
302  return m_xmap->PhysEvaluate(Lcoord, tmp);
303 }
StdRegions::StdExpansionSharedPtr m_xmap
mapping containing isoparametric transformation.
Definition: Geometry.h:189
Array< OneD, Array< OneD, NekDouble > > m_coeffs
Array containing expansion coefficients of m_xmap.
Definition: Geometry.h:201
Geometric information has been generated.
GeomState m_state
Enumeration to dictate whether coefficients are filled.
Definition: Geometry.h:191
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:250

◆ v_GetDir()

virtual int Nektar::SpatialDomains::Geometry3D::v_GetDir ( const int  faceidx,
const int  facedir 
) const
protectedpure virtual

◆ v_GetEdge()

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

Returns edge i of this object.

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 334 of file Geometry3D.cpp.

References ASSERTL2, and m_edges.

335 {
336  ASSERTL2(i >= 0 && i <= m_edges.size() - 1,
337  "Edge ID must be between 0 and " +
338  boost::lexical_cast<string>(m_edges.size() - 1));
339  return m_edges[i];
340 }
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed t...
Definition: ErrorUtil.hpp:274

◆ v_GetEorient()

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

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 348 of file Geometry3D.cpp.

References ASSERTL2, m_edges, and m_eorient.

349 {
350  ASSERTL2(i >= 0 && i <= m_edges.size() - 1,
351  "Edge ID must be between 0 and " +
352  boost::lexical_cast<string>(m_edges.size() - 1));
353  return m_eorient[i];
354 }
std::vector< StdRegions::Orientation > m_eorient
Definition: Geometry3D.h:86
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed t...
Definition: ErrorUtil.hpp:274

◆ v_GetFace()

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

Returns face i of this object.

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 342 of file Geometry3D.cpp.

References ASSERTL2, and m_faces.

343 {
344  ASSERTL2((i >= 0) && (i <= 5), "Edge id must be between 0 and 4");
345  return m_faces[i];
346 }
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed t...
Definition: ErrorUtil.hpp:274

◆ v_GetForient()

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

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

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 356 of file Geometry3D.cpp.

References ASSERTL2, m_faces, and m_forient.

357 {
358  ASSERTL2(i >= 0 && i <= m_faces.size() - 1,
359  "Face ID must be between 0 and " +
360  boost::lexical_cast<string>(m_faces.size() - 1));
361  return m_forient[i];
362 }
std::vector< StdRegions::Orientation > m_forient
Definition: Geometry3D.h:87
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed t...
Definition: ErrorUtil.hpp:274

◆ v_GetNumEdges()

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

Get the number of edges of this object.

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 319 of file Geometry3D.cpp.

References m_edges.

320 {
321  return m_edges.size();
322 }

◆ v_GetNumFaces()

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

Get the number of faces of this object.

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 324 of file Geometry3D.cpp.

References m_faces.

325 {
326  return m_faces.size();
327 }

◆ v_GetNumVerts()

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

Get the number of vertices of this object.

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 314 of file Geometry3D.cpp.

References m_verts.

315 {
316  return m_verts.size();
317 }

◆ v_GetShapeDim()

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

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 309 of file Geometry3D.cpp.

310 {
311  return 3;
312 }

◆ v_GetVertex()

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

Implements Nektar::SpatialDomains::Geometry.

Definition at line 329 of file Geometry3D.cpp.

References m_verts.

330 {
331  return m_verts[i];
332 }

Member Data Documentation

◆ kDim

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

Definition at line 80 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 88 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 89 of file Geometry3D.h.

◆ m_verts

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