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

#include <PyrGeom.h>

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

Public Member Functions

 PyrGeom ()
 
 PyrGeom (int id, const Geometry2DSharedPtr faces[])
 
 ~PyrGeom ()
 
- Public Member Functions inherited from Nektar::SpatialDomains::Geometry3D
 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 kNverts = 5
 
static const int kNedges = 8
 
static const int kNqfaces = 1
 
static const int kNtfaces = 4
 
static const int kNfaces = kNqfaces + kNtfaces
 
static const std::string XMLElementType
 
- Static Public Attributes inherited from Nektar::SpatialDomains::Geometry3D
static const int kDim = 3
 

Protected Member Functions

virtual void v_GenGeomFactors ()
 
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 ()
 
- Protected Member Functions inherited from Nektar::SpatialDomains::Geometry3D
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 ()
 

Private Member Functions

void SetUpLocalEdges ()
 
void SetUpLocalVertices ()
 
void SetUpEdgeOrientation ()
 
void SetUpFaceOrientation ()
 
void SetUpXmap ()
 Set up the m_xmap object by determining the order of each direction from derived faces. More...
 

Static Private Attributes

static const unsigned int EdgeNormalToFaceVert [5][4]
 

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

Detailed Description

Definition at line 46 of file PyrGeom.h.

Constructor & Destructor Documentation

◆ PyrGeom() [1/2]

Nektar::SpatialDomains::PyrGeom::PyrGeom ( )

Definition at line 52 of file PyrGeom.cpp.

53{
55}
LibUtilities::ShapeType m_shapeType
Type of shape.
Definition: Geometry.h:200

References Nektar::LibUtilities::ePyramid, and Nektar::SpatialDomains::Geometry::m_shapeType.

◆ PyrGeom() [2/2]

Nektar::SpatialDomains::PyrGeom::PyrGeom ( int  id,
const Geometry2DSharedPtr  faces[] 
)

Copy the face shared pointers

Set up orientation vectors with correct amount of elements.

Definition at line 57 of file PyrGeom.cpp.

58 : Geometry3D(faces[0]->GetEdge(0)->GetVertex(0)->GetCoordim())
59{
61 m_globalID = id;
62
63 /// Copy the face shared pointers
64 m_faces.insert(m_faces.begin(), faces, faces + PyrGeom::kNfaces);
65
66 /// Set up orientation vectors with correct amount of elements.
67 m_eorient.resize(kNedges);
68 m_forient.resize(kNfaces);
69
74}
std::vector< StdRegions::Orientation > m_forient
Definition: Geometry3D.h:84
std::vector< StdRegions::Orientation > m_eorient
Definition: Geometry3D.h:83
PointGeomSharedPtr GetVertex(int i) const
Returns vertex i of this object.
Definition: Geometry.h:358
Geometry1DSharedPtr GetEdge(int i) const
Returns edge i of this object.
Definition: Geometry.h:366
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
static const int kNfaces
Definition: PyrGeom.h:57
static const int kNedges
Definition: PyrGeom.h:54

References Nektar::LibUtilities::ePyramid, kNedges, kNfaces, Nektar::SpatialDomains::Geometry3D::m_eorient, Nektar::SpatialDomains::Geometry3D::m_faces, Nektar::SpatialDomains::Geometry3D::m_forient, Nektar::SpatialDomains::Geometry::m_globalID, Nektar::SpatialDomains::Geometry::m_shapeType, SetUpEdgeOrientation(), SetUpFaceOrientation(), SetUpLocalEdges(), and SetUpLocalVertices().

◆ ~PyrGeom()

Nektar::SpatialDomains::PyrGeom::~PyrGeom ( )

Definition at line 76 of file PyrGeom.cpp.

77{
78}

Member Function Documentation

◆ SetUpEdgeOrientation()

void Nektar::SpatialDomains::PyrGeom::SetUpEdgeOrientation ( )
private

Definition at line 355 of file PyrGeom.cpp.

356{
357 // This 2D array holds the local id's of all the vertices for every
358 // edge. For every edge, they are ordered to what we define as being
359 // Forwards.
360 const unsigned int edgeVerts[kNedges][2] = {{0, 1}, {1, 2}, {3, 2}, {0, 3},
361 {0, 4}, {1, 4}, {2, 4}, {3, 4}};
362
363 int i;
364 for (i = 0; i < kNedges; i++)
365 {
366 if (m_edges[i]->GetVid(0) == m_verts[edgeVerts[i][0]]->GetGlobalID())
367 {
369 }
370 else if (m_edges[i]->GetVid(0) ==
371 m_verts[edgeVerts[i][1]]->GetGlobalID())
372 {
374 }
375 else
376 {
377 ASSERTL0(false, "Could not find matching vertex for the edge");
378 }
379 }
380}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
int GetVid(int i) const
Get the ID of vertex i of this object.
Definition: Geometry.cpp:139
int GetGlobalID(void) const
Get the ID of this object.
Definition: Geometry.h:327

References ASSERTL0, Nektar::StdRegions::eBackwards, Nektar::StdRegions::eForwards, Nektar::SpatialDomains::Geometry::GetGlobalID(), Nektar::SpatialDomains::Geometry::GetVid(), kNedges, Nektar::SpatialDomains::Geometry3D::m_edges, Nektar::SpatialDomains::Geometry3D::m_eorient, and Nektar::SpatialDomains::Geometry3D::m_verts.

Referenced by PyrGeom().

◆ SetUpFaceOrientation()

void Nektar::SpatialDomains::PyrGeom::SetUpFaceOrientation ( )
private

Definition at line 382 of file PyrGeom.cpp.

383{
384 int f, i;
385
386 // These arrays represent the vector of the A and B coordinate of
387 // the local elemental coordinate system where A corresponds with
388 // the coordinate direction xi_i with the lowest index i (for that
389 // particular face) Coordinate 'B' then corresponds to the other
390 // local coordinate (i.e. with the highest index)
391 Array<OneD, NekDouble> elementAaxis(m_coordim);
392 Array<OneD, NekDouble> elementBaxis(m_coordim);
393
394 // These arrays correspond to the local coordinate
395 // system of the face itself (i.e. the Geometry2D)
396 // faceAaxis correspond to the xi_0 axis
397 // faceBaxis correspond to the xi_1 axis
398 Array<OneD, NekDouble> faceAaxis(m_coordim);
399 Array<OneD, NekDouble> faceBaxis(m_coordim);
400
401 // This is the base vertex of the face (i.e. the Geometry2D) This
402 // corresponds to thevertex with local ID 0 of the Geometry2D
403 unsigned int baseVertex;
404
405 // The lenght of the vectors above
406 NekDouble elementAaxis_length;
407 NekDouble elementBaxis_length;
408 NekDouble faceAaxis_length;
409 NekDouble faceBaxis_length;
410
411 // This 2D array holds the local id's of all the vertices for every
412 // face. For every face, they are ordered in such a way that the
413 // implementation below allows a unified approach for all faces.
414 const unsigned int faceVerts[kNfaces][4] = {
415 {0, 1, 2, 3},
416 {0, 1, 4, 0}, // Last four elements are triangles which only
417 {1, 2, 4, 0}, // require three vertices.
418 {3, 2, 4, 0},
419 {0, 3, 4, 0}};
420
421 NekDouble dotproduct1 = 0.0;
422 NekDouble dotproduct2 = 0.0;
423
424 unsigned int orientation;
425
426 // Loop over all the faces to set up the orientation
427 for (f = 0; f < kNqfaces + kNtfaces; f++)
428 {
429 // initialisation
430 elementAaxis_length = 0.0;
431 elementBaxis_length = 0.0;
432 faceAaxis_length = 0.0;
433 faceBaxis_length = 0.0;
434
435 dotproduct1 = 0.0;
436 dotproduct2 = 0.0;
437
438 baseVertex = m_faces[f]->GetVid(0);
439
440 // We are going to construct the vectors representing the A and
441 // B axis of every face. These vectors will be constructed as a
442 // vector-representation of the edges of the face. However, for
443 // both coordinate directions, we can represent the vectors by
444 // two different edges. That's why we need to make sure that we
445 // pick the edge to which the baseVertex of the
446 // Geometry2D-representation of the face belongs...
447
448 // Compute the length of edges on a base-face
449 if (f > 0)
450 {
451 if (baseVertex == m_verts[faceVerts[f][0]]->GetVid())
452 {
453 for (i = 0; i < m_coordim; i++)
454 {
455 elementAaxis[i] = (*m_verts[faceVerts[f][1]])[i] -
456 (*m_verts[faceVerts[f][0]])[i];
457 elementBaxis[i] = (*m_verts[faceVerts[f][2]])[i] -
458 (*m_verts[faceVerts[f][0]])[i];
459 }
460 }
461 else if (baseVertex == m_verts[faceVerts[f][1]]->GetVid())
462 {
463 for (i = 0; i < m_coordim; i++)
464 {
465 elementAaxis[i] = (*m_verts[faceVerts[f][1]])[i] -
466 (*m_verts[faceVerts[f][0]])[i];
467 elementBaxis[i] = (*m_verts[faceVerts[f][2]])[i] -
468 (*m_verts[faceVerts[f][1]])[i];
469 }
470 }
471 else if (baseVertex == m_verts[faceVerts[f][2]]->GetVid())
472 {
473 for (i = 0; i < m_coordim; i++)
474 {
475 elementAaxis[i] = (*m_verts[faceVerts[f][1]])[i] -
476 (*m_verts[faceVerts[f][2]])[i];
477 elementBaxis[i] = (*m_verts[faceVerts[f][2]])[i] -
478 (*m_verts[faceVerts[f][0]])[i];
479 }
480 }
481 else
482 {
483 ASSERTL0(false, "Could not find matching vertex for the face");
484 }
485 }
486 else
487 {
488 if (baseVertex == m_verts[faceVerts[f][0]]->GetGlobalID())
489 {
490 for (i = 0; i < m_coordim; i++)
491 {
492 elementAaxis[i] = (*m_verts[faceVerts[f][1]])[i] -
493 (*m_verts[faceVerts[f][0]])[i];
494 elementBaxis[i] = (*m_verts[faceVerts[f][3]])[i] -
495 (*m_verts[faceVerts[f][0]])[i];
496 }
497 }
498 else if (baseVertex == m_verts[faceVerts[f][1]]->GetGlobalID())
499 {
500 for (i = 0; i < m_coordim; i++)
501 {
502 elementAaxis[i] = (*m_verts[faceVerts[f][1]])[i] -
503 (*m_verts[faceVerts[f][0]])[i];
504 elementBaxis[i] = (*m_verts[faceVerts[f][2]])[i] -
505 (*m_verts[faceVerts[f][1]])[i];
506 }
507 }
508 else if (baseVertex == m_verts[faceVerts[f][2]]->GetGlobalID())
509 {
510 for (i = 0; i < m_coordim; i++)
511 {
512 elementAaxis[i] = (*m_verts[faceVerts[f][2]])[i] -
513 (*m_verts[faceVerts[f][3]])[i];
514 elementBaxis[i] = (*m_verts[faceVerts[f][2]])[i] -
515 (*m_verts[faceVerts[f][1]])[i];
516 }
517 }
518 else if (baseVertex == m_verts[faceVerts[f][3]]->GetGlobalID())
519 {
520 for (i = 0; i < m_coordim; i++)
521 {
522 elementAaxis[i] = (*m_verts[faceVerts[f][2]])[i] -
523 (*m_verts[faceVerts[f][3]])[i];
524 elementBaxis[i] = (*m_verts[faceVerts[f][3]])[i] -
525 (*m_verts[faceVerts[f][0]])[i];
526 }
527 }
528 else
529 {
530 ASSERTL0(false, "Could not find matching vertex for the face");
531 }
532 }
533
534 // Now, construct the edge-vectors of the local coordinates of
535 // the Geometry2D-representation of the face
536 for (i = 0; i < m_coordim; i++)
537 {
538 int v = m_faces[f]->GetNumVerts() - 1;
539 faceAaxis[i] =
540 (*m_faces[f]->GetVertex(1))[i] - (*m_faces[f]->GetVertex(0))[i];
541 faceBaxis[i] =
542 (*m_faces[f]->GetVertex(v))[i] - (*m_faces[f]->GetVertex(0))[i];
543
544 elementAaxis_length += pow(elementAaxis[i], 2);
545 elementBaxis_length += pow(elementBaxis[i], 2);
546 faceAaxis_length += pow(faceAaxis[i], 2);
547 faceBaxis_length += pow(faceBaxis[i], 2);
548 }
549
550 elementAaxis_length = sqrt(elementAaxis_length);
551 elementBaxis_length = sqrt(elementBaxis_length);
552 faceAaxis_length = sqrt(faceAaxis_length);
553 faceBaxis_length = sqrt(faceBaxis_length);
554
555 // Calculate the inner product of both the A-axis
556 // (i.e. Elemental A axis and face A axis)
557 for (i = 0; i < m_coordim; i++)
558 {
559 dotproduct1 += elementAaxis[i] * faceAaxis[i];
560 }
561
562 NekDouble norm =
563 fabs(dotproduct1) / elementAaxis_length / faceAaxis_length;
564 orientation = 0;
565
566 // if the innerproduct is equal to the (absolute value of the ) products
567 // of the lengths of both vectors, then, the coordinate systems will NOT
568 // be transposed
569 if (fabs(norm - 1.0) < NekConstants::kNekZeroTol)
570 {
571 // if the inner product is negative, both A-axis point
572 // in reverse direction
573 if (dotproduct1 < 0.0)
574 {
575 orientation += 2;
576 }
577
578 // calculate the inner product of both B-axis
579 for (i = 0; i < m_coordim; i++)
580 {
581 dotproduct2 += elementBaxis[i] * faceBaxis[i];
582 }
583
584 norm = fabs(dotproduct2) / elementBaxis_length / faceBaxis_length;
585
586 // check that both these axis are indeed parallel
587 ASSERTL1(fabs(norm - 1.0) < NekConstants::kNekZeroTol,
588 "These vectors should be parallel");
589
590 // if the inner product is negative, both B-axis point
591 // in reverse direction
592 if (dotproduct2 < 0.0)
593 {
594 orientation++;
595 }
596 }
597 // The coordinate systems are transposed
598 else
599 {
600 orientation = 4;
601
602 // Calculate the inner product between the elemental A-axis
603 // and the B-axis of the face (which are now the corresponding axis)
604 dotproduct1 = 0.0;
605 for (i = 0; i < m_coordim; i++)
606 {
607 dotproduct1 += elementAaxis[i] * faceBaxis[i];
608 }
609
610 norm = fabs(dotproduct1) / elementAaxis_length / faceBaxis_length;
611 ASSERTL1(fabs(norm - 1.0) < NekConstants::kNekZeroTol,
612 "These vectors should be parallel");
613
614 // if the result is negative, both axis point in reverse
615 // directions
616 if (dotproduct1 < 0.0)
617 {
618 orientation += 2;
619 }
620
621 // Do the same for the other two corresponding axis
622 dotproduct2 = 0.0;
623 for (i = 0; i < m_coordim; i++)
624 {
625 dotproduct2 += elementBaxis[i] * faceAaxis[i];
626 }
627
628 norm = fabs(dotproduct2) / elementBaxis_length / faceAaxis_length;
629
630 // check that both these axis are indeed parallel
631 ASSERTL1(fabs(norm - 1.0) < NekConstants::kNekZeroTol,
632 "These vectors should be parallel");
633
634 if (dotproduct2 < 0.0)
635 {
636 orientation++;
637 }
638 }
639
640 orientation = orientation + 5;
641
642 if (f != 0) // check triangle orientation
643 {
644 ASSERTL0(
646 "Orientation of triangular face (id = " +
647 boost::lexical_cast<string>(m_faces[f]->GetGlobalID()) +
648 ") is inconsistent with face " +
649 boost::lexical_cast<string>(f) +
650 " of pyramid element (id = " +
651 boost::lexical_cast<string>(m_globalID) +
652 ") since Dir2 is aligned with Dir1. Mesh setup "
653 "needs investigation");
654 }
655
656 // Fill the m_forient array
657 m_forient[f] = (StdRegions::Orientation)orientation;
658 }
659}
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:249
int m_coordim
Coordinate dimension of this geometry object.
Definition: Geometry.h:186
static const int kNtfaces
Definition: PyrGeom.h:56
static const int kNqfaces
Definition: PyrGeom.h:55
static const NekDouble kNekZeroTol
double NekDouble
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294

References ASSERTL0, ASSERTL1, Nektar::StdRegions::eDir1FwdDir2_Dir2FwdDir1, Nektar::SpatialDomains::Geometry::GetGlobalID(), Nektar::SpatialDomains::Geometry::GetVertex(), Nektar::SpatialDomains::Geometry::GetVid(), Nektar::NekConstants::kNekZeroTol, kNfaces, kNqfaces, kNtfaces, Nektar::SpatialDomains::Geometry::m_coordim, Nektar::SpatialDomains::Geometry3D::m_faces, Nektar::SpatialDomains::Geometry3D::m_forient, Nektar::SpatialDomains::Geometry::m_globalID, Nektar::SpatialDomains::Geometry3D::m_verts, and tinysimd::sqrt().

Referenced by PyrGeom().

◆ SetUpLocalEdges()

void Nektar::SpatialDomains::PyrGeom::SetUpLocalEdges ( )
private

Definition at line 165 of file PyrGeom.cpp.

166{
167 // find edge 0
168 int i, j;
169 unsigned int check;
170
171 SegGeomSharedPtr edge;
172
173 // First set up the 4 bottom edges
174 int f;
175 for (f = 1; f < 5; f++)
176 {
177 int nEdges = m_faces[f]->GetNumEdges();
178 check = 0;
179 for (i = 0; i < 4; i++)
180 {
181 for (j = 0; j < nEdges; j++)
182 {
183 if (m_faces[0]->GetEid(i) == m_faces[f]->GetEid(j))
184 {
185 edge =
186 dynamic_pointer_cast<SegGeom>((m_faces[0])->GetEdge(i));
187 m_edges.push_back(edge);
188 check++;
189 }
190 }
191 }
192
193 if (check < 1)
194 {
195 std::ostringstream errstrm;
196 errstrm << "Connected faces do not share an edge. Faces ";
197 errstrm << (m_faces[0])->GetGlobalID() << ", "
198 << (m_faces[f])->GetGlobalID();
199 ASSERTL0(false, errstrm.str());
200 }
201 else if (check > 1)
202 {
203 std::ostringstream errstrm;
204 errstrm << "Connected faces share more than one edge. Faces ";
205 errstrm << (m_faces[0])->GetGlobalID() << ", "
206 << (m_faces[f])->GetGlobalID();
207 ASSERTL0(false, errstrm.str());
208 }
209 }
210
211 // Then, set up the 4 vertical edges
212 check = 0;
213 for (i = 0; i < 3; i++) // Set up the vertical edge :face(1) and face(4)
214 {
215 for (j = 0; j < 3; j++)
216 {
217 if ((m_faces[1])->GetEid(i) == (m_faces[4])->GetEid(j))
218 {
219 edge = dynamic_pointer_cast<SegGeom>((m_faces[1])->GetEdge(i));
220 m_edges.push_back(edge);
221 check++;
222 }
223 }
224 }
225 if (check < 1)
226 {
227 std::ostringstream errstrm;
228 errstrm << "Connected faces do not share an edge. Faces ";
229 errstrm << (m_faces[1])->GetGlobalID() << ", "
230 << (m_faces[4])->GetGlobalID();
231 ASSERTL0(false, errstrm.str());
232 }
233 else if (check > 1)
234 {
235 std::ostringstream errstrm;
236 errstrm << "Connected faces share more than one edge. Faces ";
237 errstrm << (m_faces[1])->GetGlobalID() << ", "
238 << (m_faces[4])->GetGlobalID();
239 ASSERTL0(false, errstrm.str());
240 }
241
242 // Set up vertical edges: face(1) through face(4)
243 for (f = 1; f < 4; f++)
244 {
245 check = 0;
246 for (i = 0; i < m_faces[f]->GetNumEdges(); i++)
247 {
248 for (j = 0; j < m_faces[f + 1]->GetNumEdges(); j++)
249 {
250 if ((m_faces[f])->GetEid(i) == (m_faces[f + 1])->GetEid(j))
251 {
252 edge =
253 dynamic_pointer_cast<SegGeom>((m_faces[f])->GetEdge(i));
254 m_edges.push_back(edge);
255 check++;
256 }
257 }
258 }
259
260 if (check < 1)
261 {
262 std::ostringstream errstrm;
263 errstrm << "Connected faces do not share an edge. Faces ";
264 errstrm << (m_faces[f])->GetGlobalID() << ", "
265 << (m_faces[f + 1])->GetGlobalID();
266 ASSERTL0(false, errstrm.str());
267 }
268 else if (check > 1)
269 {
270 std::ostringstream errstrm;
271 errstrm << "Connected faces share more than one edge. Faces ";
272 errstrm << (m_faces[f])->GetGlobalID() << ", "
273 << (m_faces[f + 1])->GetGlobalID();
274 ASSERTL0(false, errstrm.str());
275 }
276 }
277}
int GetEid(int i) const
Get the ID of edge i of this object.
Definition: Geometry.cpp:147
std::shared_ptr< SegGeom > SegGeomSharedPtr
Definition: Geometry2D.h:62

References ASSERTL0, Nektar::SpatialDomains::Geometry::GetEdge(), Nektar::SpatialDomains::Geometry::GetEid(), Nektar::SpatialDomains::Geometry::GetGlobalID(), Nektar::SpatialDomains::Geometry3D::m_edges, and Nektar::SpatialDomains::Geometry3D::m_faces.

Referenced by PyrGeom().

◆ SetUpLocalVertices()

void Nektar::SpatialDomains::PyrGeom::SetUpLocalVertices ( )
private

Definition at line 279 of file PyrGeom.cpp.

280{
281 // Set up the first 2 vertices (i.e. vertex 0,1)
282 if (m_edges[0]->GetVid(0) == m_edges[1]->GetVid(0) ||
283 m_edges[0]->GetVid(0) == m_edges[1]->GetVid(1))
284 {
285 m_verts.push_back(m_edges[0]->GetVertex(1));
286 m_verts.push_back(m_edges[0]->GetVertex(0));
287 }
288 else if (m_edges[0]->GetVid(1) == m_edges[1]->GetVid(0) ||
289 m_edges[0]->GetVid(1) == m_edges[1]->GetVid(1))
290 {
291 m_verts.push_back(m_edges[0]->GetVertex(0));
292 m_verts.push_back(m_edges[0]->GetVertex(1));
293 }
294 else
295 {
296 std::ostringstream errstrm;
297 errstrm << "Connected edges do not share a vertex. Edges ";
298 errstrm << m_edges[0]->GetGlobalID() << ", "
299 << m_edges[1]->GetGlobalID();
300 ASSERTL0(false, errstrm.str());
301 }
302
303 // set up the other bottom vertices (i.e. vertex 2,3)
304 for (int i = 1; i < 3; i++)
305 {
306 if (m_edges[i]->GetVid(0) == m_verts[i]->GetGlobalID())
307 {
308 m_verts.push_back(m_edges[i]->GetVertex(1));
309 }
310 else if (m_edges[i]->GetVid(1) == m_verts[i]->GetGlobalID())
311 {
312 m_verts.push_back(m_edges[i]->GetVertex(0));
313 }
314 else
315 {
316 std::ostringstream errstrm;
317 errstrm << "Connected edges do not share a vertex. Edges ";
318 errstrm << m_edges[i]->GetGlobalID() << ", "
319 << m_edges[i - 1]->GetGlobalID();
320 ASSERTL0(false, errstrm.str());
321 }
322 }
323
324 // set up top vertex
325 if (m_edges[4]->GetVid(0) == m_verts[0]->GetGlobalID())
326 {
327 m_verts.push_back(m_edges[4]->GetVertex(1));
328 }
329 else
330 {
331 m_verts.push_back(m_edges[4]->GetVertex(0));
332 }
333
334 int check = 0;
335 for (int i = 5; i < 8; ++i)
336 {
337 if ((m_edges[i]->GetVid(0) == m_verts[i - 4]->GetGlobalID() &&
338 m_edges[i]->GetVid(1) == m_verts[4]->GetGlobalID()) ||
339 (m_edges[i]->GetVid(1) == m_verts[i - 4]->GetGlobalID() &&
340 m_edges[i]->GetVid(0) == m_verts[4]->GetGlobalID()))
341 {
342 check++;
343 }
344 }
345 if (check != 3)
346 {
347 std::ostringstream errstrm;
348 errstrm << "Connected edges do not share a vertex. Edges ";
349 errstrm << m_edges[3]->GetGlobalID() << ", "
350 << m_edges[2]->GetGlobalID();
351 ASSERTL0(false, errstrm.str());
352 }
353}

References ASSERTL0, Nektar::SpatialDomains::Geometry::GetGlobalID(), Nektar::SpatialDomains::Geometry::GetVertex(), Nektar::SpatialDomains::Geometry::GetVid(), Nektar::SpatialDomains::Geometry3D::m_edges, and Nektar::SpatialDomains::Geometry3D::m_verts.

Referenced by PyrGeom().

◆ SetUpXmap()

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

Set up the m_xmap object by determining the order of each direction from derived faces.

Definition at line 692 of file PyrGeom.cpp.

693{
694 vector<int> tmp;
695 int order0, order1;
696
697 if (m_forient[0] < 9)
698 {
699 tmp.push_back(m_faces[0]->GetXmap()->GetTraceNcoeffs(0));
700 tmp.push_back(m_faces[0]->GetXmap()->GetTraceNcoeffs(2));
701 order0 = *max_element(tmp.begin(), tmp.end());
702 }
703 else
704 {
705 tmp.push_back(m_faces[0]->GetXmap()->GetTraceNcoeffs(1));
706 tmp.push_back(m_faces[0]->GetXmap()->GetTraceNcoeffs(3));
707 order0 = *max_element(tmp.begin(), tmp.end());
708 }
709
710 if (m_forient[0] < 9)
711 {
712 tmp.clear();
713 tmp.push_back(m_faces[0]->GetXmap()->GetTraceNcoeffs(1));
714 tmp.push_back(m_faces[0]->GetXmap()->GetTraceNcoeffs(3));
715 tmp.push_back(m_faces[2]->GetXmap()->GetTraceNcoeffs(2));
716 order1 = *max_element(tmp.begin(), tmp.end());
717 }
718 else
719 {
720 tmp.clear();
721 tmp.push_back(m_faces[0]->GetXmap()->GetTraceNcoeffs(0));
722 tmp.push_back(m_faces[0]->GetXmap()->GetTraceNcoeffs(2));
723 tmp.push_back(m_faces[2]->GetXmap()->GetTraceNcoeffs(2));
724 order1 = *max_element(tmp.begin(), tmp.end());
725 }
726
727 tmp.clear();
728 tmp.push_back(order0);
729 tmp.push_back(order1);
730 tmp.push_back(m_faces[1]->GetXmap()->GetTraceNcoeffs(1));
731 tmp.push_back(m_faces[1]->GetXmap()->GetTraceNcoeffs(2));
732 tmp.push_back(m_faces[3]->GetXmap()->GetTraceNcoeffs(1));
733 tmp.push_back(m_faces[3]->GetXmap()->GetTraceNcoeffs(2));
734 int order2 = *max_element(tmp.begin(), tmp.end());
735
736 const LibUtilities::BasisKey A1(
738 LibUtilities::PointsKey(order0 + 1,
740 const LibUtilities::BasisKey A2(
742 LibUtilities::PointsKey(order1 + 1,
744 const LibUtilities::BasisKey C(
746 LibUtilities::PointsKey(order2, LibUtilities::eGaussRadauMAlpha2Beta0));
747
749}
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
StdRegions::StdExpansionSharedPtr m_xmap
mapping containing isoparametric transformation.
Definition: Geometry.h:192
StdRegions::StdExpansionSharedPtr GetXmap() const
Return the mapping object Geometry::m_xmap that represents the coordinate transformation from standar...
Definition: Geometry.h:436
@ eGaussLobattoLegendre
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:53
@ eModifiedPyr_C
Principle Modified Functions.
Definition: BasisType.h:55
@ eModified_A
Principle Modified Functions .
Definition: BasisType.h:50

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::LibUtilities::eGaussLobattoLegendre, Nektar::LibUtilities::eModified_A, Nektar::LibUtilities::eModifiedPyr_C, Nektar::SpatialDomains::Geometry::GetXmap(), Nektar::SpatialDomains::Geometry3D::m_faces, Nektar::SpatialDomains::Geometry3D::m_forient, and Nektar::SpatialDomains::Geometry::m_xmap.

Referenced by v_Reset(), and v_Setup().

◆ v_GenGeomFactors()

void Nektar::SpatialDomains::PyrGeom::v_GenGeomFactors ( )
protectedvirtual

Implements Nektar::SpatialDomains::Geometry.

Definition at line 80 of file PyrGeom.cpp.

81{
82 if (!m_setupState)
83 {
85 }
86
88 {
89 GeomType Gtype = eRegular;
90
91 v_FillGeom();
92
93 // check to see if expansions are linear
94 m_straightEdge = true;
95 if (m_xmap->GetBasisNumModes(0) != 2 ||
96 m_xmap->GetBasisNumModes(1) != 2 ||
97 m_xmap->GetBasisNumModes(2) != 2)
98 {
99 Gtype = eDeformed;
100 m_straightEdge = false;
101 }
102
103 // check to see if all quadrilateral faces are parallelograms
104 if (Gtype == eRegular)
105 {
106 m_isoParameter = Array<OneD, Array<OneD, NekDouble>>(3);
107 for (int i = 0; i < 3; ++i)
108 {
109 m_isoParameter[i] = Array<OneD, NekDouble>(5, 0.);
110 NekDouble A = (*m_verts[0])(i);
111 NekDouble B = (*m_verts[1])(i);
112 NekDouble C = (*m_verts[2])(i);
113 NekDouble D = (*m_verts[3])(i);
114 NekDouble E = (*m_verts[4])(i);
115 m_isoParameter[i][0] = 0.25 * (-A + B + C + D + E + E);
116
117 m_isoParameter[i][1] = 0.25 * (-A + B + C - D); // xi1
118 m_isoParameter[i][2] = 0.25 * (-A - B + C + D); // xi2
119 m_isoParameter[i][3] = 0.5 * (-A + E); // xi3
120
121 m_isoParameter[i][4] = 0.25 * (A - B + C - D); // xi1*xi2
122 NekDouble tmp = fabs(m_isoParameter[i][1]) +
123 fabs(m_isoParameter[i][2]) +
124 fabs(m_isoParameter[i][3]);
125 if (fabs(m_isoParameter[i][4]) >
127 {
128 Gtype = eDeformed;
129 }
130 }
131 }
132
133 if (Gtype == eRegular)
134 {
136 }
137
139 Gtype, m_coordim, m_xmap, m_coeffs);
141 }
142}
virtual void v_CalculateInverseIsoParam() override
Definition: Geometry3D.cpp:572
virtual void v_FillGeom() override
Put all quadrature information into face/edge structure and backward transform.
Definition: Geometry3D.cpp:449
bool m_setupState
Wether or not the setup routines have been run.
Definition: Geometry.h:196
Array< OneD, Array< OneD, NekDouble > > m_isoParameter
Definition: Geometry.h:207
Array< OneD, Array< OneD, NekDouble > > m_coeffs
Array containing expansion coefficients of m_xmap.
Definition: Geometry.h:204
GeomState m_geomFactorsState
State of the geometric factors.
Definition: Geometry.h:190
GeomFactorsSharedPtr m_geomFactors
Geometric factors.
Definition: Geometry.h:188
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.
@ ePtsFilled
Geometric information has been generated.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), 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::Geometry::m_setupState, Nektar::SpatialDomains::Geometry::m_straightEdge, Nektar::SpatialDomains::Geometry3D::m_verts, Nektar::SpatialDomains::Geometry::m_xmap, Nektar::SpatialDomains::Geometry3D::v_CalculateInverseIsoParam(), Nektar::SpatialDomains::Geometry3D::v_FillGeom(), and v_Setup().

◆ v_GetDir()

int Nektar::SpatialDomains::PyrGeom::v_GetDir ( const int  i,
const int  j 
) const
protectedvirtual

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

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 144 of file PyrGeom.cpp.

145{
146 if (faceidx == 0)
147 {
148 return facedir;
149 }
150 else if (faceidx == 1 || faceidx == 3)
151 {
152 return 2 * facedir;
153 }
154 else
155 {
156 return 1 + facedir;
157 }
158}

◆ v_GetEdgeNormalToFaceVert()

int Nektar::SpatialDomains::PyrGeom::v_GetEdgeNormalToFaceVert ( const int  i,
const int  j 
) const
protectedvirtual

Returns the standard lement edge IDs that are normal to a given face vertex.

For example, on a hexahedron, on face 0 at vertices 0,1,2,3 the edges normal to that face are 4,5,6,7, ; so GetEdgeNormalToFaceVert(0,j) would therefore return the values 4, 5, 6 and 7 respectively. We assume that j runs between 0 and 3 inclusive on a quadrilateral face and between 0 and 2 inclusive on a triangular face.

This is used to help set up a length scale normal to an face

Parameters
iThe face to query for the normal edge
jThe local vertex index between 0 and nverts on this face

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 160 of file PyrGeom.cpp.

161{
162 return EdgeNormalToFaceVert[i][j];
163}
static const unsigned int EdgeNormalToFaceVert[5][4]
Definition: PyrGeom.h:74

References EdgeNormalToFaceVert.

◆ v_Reset()

void Nektar::SpatialDomains::PyrGeom::v_Reset ( CurveMap curvedEdges,
CurveMap curvedFaces 
)
protectedvirtual

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

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 661 of file PyrGeom.cpp.

662{
663 Geometry::v_Reset(curvedEdges, curvedFaces);
664
665 for (int i = 0; i < 5; ++i)
666 {
667 m_faces[i]->Reset(curvedEdges, curvedFaces);
668 }
669
670 SetUpXmap();
671 SetUpCoeffs(m_xmap->GetNcoeffs());
672}
void SetUpCoeffs(const int nCoeffs)
Initialise the Geometry::m_coeffs array.
Definition: Geometry.h:690
virtual void v_Reset(CurveMap &curvedEdges, CurveMap &curvedFaces)
Reset this geometry object: unset the current state, zero Geometry::m_coeffs and remove allocated Geo...
Definition: Geometry.cpp:390
void SetUpXmap()
Set up the m_xmap object by determining the order of each direction from derived faces.
Definition: PyrGeom.cpp:692

References Nektar::SpatialDomains::Geometry3D::m_faces, Nektar::SpatialDomains::Geometry::m_xmap, Nektar::SpatialDomains::Geometry::SetUpCoeffs(), SetUpXmap(), and Nektar::SpatialDomains::Geometry::v_Reset().

◆ v_Setup()

void Nektar::SpatialDomains::PyrGeom::v_Setup ( )
protectedvirtual

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 674 of file PyrGeom.cpp.

675{
676 if (!m_setupState)
677 {
678 for (int i = 0; i < 5; ++i)
679 {
680 m_faces[i]->Setup();
681 }
682 SetUpXmap();
683 SetUpCoeffs(m_xmap->GetNcoeffs());
684 m_setupState = true;
685 }
686}

References Nektar::SpatialDomains::Geometry3D::m_faces, Nektar::SpatialDomains::Geometry::m_setupState, Nektar::SpatialDomains::Geometry::m_xmap, Nektar::SpatialDomains::Geometry::SetUpCoeffs(), and SetUpXmap().

Referenced by v_GenGeomFactors().

Member Data Documentation

◆ EdgeNormalToFaceVert

const unsigned int Nektar::SpatialDomains::PyrGeom::EdgeNormalToFaceVert
staticprivate
Initial value:
= {
{4, 5, 6, 7}, {1, 3, 6, 7}, {0, 2, 4, 7}, {1, 3, 4, 5}, {0, 2, 5, 6}}

Definition at line 74 of file PyrGeom.h.

Referenced by v_GetEdgeNormalToFaceVert().

◆ kNedges

const int Nektar::SpatialDomains::PyrGeom::kNedges = 8
static

Definition at line 54 of file PyrGeom.h.

Referenced by PyrGeom(), and SetUpEdgeOrientation().

◆ kNfaces

const int Nektar::SpatialDomains::PyrGeom::kNfaces = kNqfaces + kNtfaces
static

◆ kNqfaces

const int Nektar::SpatialDomains::PyrGeom::kNqfaces = 1
static

◆ kNtfaces

const int Nektar::SpatialDomains::PyrGeom::kNtfaces = 4
static

◆ kNverts

const int Nektar::SpatialDomains::PyrGeom::kNverts = 5
static

Definition at line 53 of file PyrGeom.h.

◆ XMLElementType

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

Definition at line 58 of file PyrGeom.h.