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 () override
 
- Public Member Functions inherited from Nektar::SpatialDomains::Geometry3D
 Geometry3D ()
 
 Geometry3D (const int coordim)
 
 ~Geometry3D () override
 
- Public Member Functions inherited from Nektar::SpatialDomains::Geometry
 Geometry ()
 Default constructor. More...
 
 Geometry (int coordim)
 Constructor when supplied a coordinate dimension. More...
 
virtual ~Geometry ()
 Default destructor. More...
 
int GetCoordim () const
 Return the coordinate dimension of this object (i.e. the dimension of the space in which this object is embedded). More...
 
void SetCoordim (int coordim)
 Sets the coordinate dimension of this object (i.e. the dimension of the space in which this object is embedded). More...
 
GeomFactorsSharedPtr GetGeomFactors ()
 Get the geometric factors for this object, generating them if required. More...
 
GeomFactorsSharedPtr GetRefGeomFactors (const Array< OneD, const LibUtilities::BasisSharedPtr > &tbasis)
 
GeomFactorsSharedPtr GetMetricInfo ()
 Get the geometric factors for this object. More...
 
LibUtilities::ShapeType GetShapeType (void)
 Get the geometric shape type of this object. More...
 
int GetGlobalID (void) const
 Get the ID of this object. More...
 
void SetGlobalID (int globalid)
 Set the ID of this object. More...
 
int GetVid (int i) const
 Get the ID of vertex i of this object. More...
 
int GetEid (int i) const
 Get the ID of edge i of this object. More...
 
int GetFid (int i) const
 Get the ID of face i of this object. More...
 
int GetTid (int i) const
 Get the ID of trace i of this object. More...
 
PointGeomSharedPtr GetVertex (int i) const
 Returns vertex i of this object. More...
 
Geometry1DSharedPtr GetEdge (int i) const
 Returns edge i of this object. More...
 
Geometry2DSharedPtr GetFace (int i) const
 Returns face i of this object. More...
 
StdRegions::Orientation GetEorient (const int i) const
 Returns the orientation of edge i with respect to the ordering of edges in the standard element. More...
 
StdRegions::Orientation GetForient (const int i) const
 Returns the orientation of face i with respect to the ordering of faces in the standard element. More...
 
int GetNumVerts () const
 Get the number of vertices of this object. More...
 
int GetNumEdges () const
 Get the number of edges of this object. More...
 
int GetNumFaces () const
 Get the number of faces of this object. More...
 
int GetShapeDim () const
 Get the object's shape dimension. More...
 
StdRegions::StdExpansionSharedPtr GetXmap () const
 Return the mapping object Geometry::m_xmap that represents the coordinate transformation from standard element to physical element. More...
 
const Array< OneD, const NekDouble > & GetCoeffs (const int i) const
 Return the coefficients of the transformation Geometry::m_xmap in coordinate direction i. More...
 
void FillGeom ()
 Populate the coordinate mapping Geometry::m_coeffs information from any children geometry elements. More...
 
std::array< NekDouble, 6 > GetBoundingBox ()
 Generates the bounding box for the element. More...
 
void ClearBoundingBox ()
 
bool ContainsPoint (const Array< OneD, const NekDouble > &gloCoord, NekDouble tol=0.0)
 Determine whether an element contains a particular Cartesian coordinate \((x,y,z)\). More...
 
bool ContainsPoint (const Array< OneD, const NekDouble > &gloCoord, Array< OneD, NekDouble > &locCoord, NekDouble tol)
 Determine whether an element contains a particular Cartesian coordinate \((x,y,z)\). More...
 
bool ContainsPoint (const Array< OneD, const NekDouble > &gloCoord, Array< OneD, NekDouble > &locCoord, NekDouble tol, NekDouble &dist)
 Determine whether an element contains a particular Cartesian coordinate \(\vec{x} = (x,y,z)\). More...
 
NekDouble GetLocCoords (const Array< OneD, const NekDouble > &coords, Array< OneD, NekDouble > &Lcoords)
 Determine the local collapsed coordinates that correspond to a given Cartesian coordinate for this geometry object. More...
 
NekDouble GetCoord (const int i, const Array< OneD, const NekDouble > &Lcoord)
 Given local collapsed coordinate Lcoord, return the value of physical coordinate in direction i. More...
 
int PreliminaryCheck (const Array< OneD, const NekDouble > &gloCoord)
 A fast and robust check if a given global coord is outside of a deformed element. For regular elements, this check is unnecessary. More...
 
bool MinMaxCheck (const Array< OneD, const NekDouble > &gloCoord)
 Check if given global coord is within the BoundingBox of the element. More...
 
bool ClampLocCoords (Array< OneD, NekDouble > &locCoord, NekDouble tol=std::numeric_limits< NekDouble >::epsilon())
 Clamp local coords to be within standard regions [-1, 1]^dim. More...
 
NekDouble FindDistance (const Array< OneD, const NekDouble > &xs, Array< OneD, NekDouble > &xi)
 
int GetVertexEdgeMap (int i, int j) const
 Returns the standard element edge IDs that are connected to a given vertex. More...
 
int GetVertexFaceMap (int i, int j) const
 Returns the standard element face IDs that are connected to a given vertex. More...
 
int GetEdgeFaceMap (int i, int j) const
 Returns the standard element edge IDs that are connected to a given face. More...
 
int GetEdgeNormalToFaceVert (int i, int j) const
 Returns the standard lement edge IDs that are normal to a given face vertex. More...
 
int GetDir (const int i, const int j=0) const
 Returns the element coordinate direction corresponding to a given face coordinate direction. More...
 
void Reset (CurveMap &curvedEdges, CurveMap &curvedFaces)
 Reset this geometry object: unset the current state, zero Geometry::m_coeffs and remove allocated GeomFactors. More...
 
void ResetNonRecursive (CurveMap &curvedEdges, CurveMap &curvedFaces)
 Reset this geometry object non-recursively: unset the current state, zero Geometry::m_coeffs and remove allocated GeomFactors. More...
 
void Setup ()
 
void GenGeomFactors ()
 Handles generation of geometry factors. More...
 

Static Public 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

void v_GenGeomFactors () override
 
int v_GetEdgeNormalToFaceVert (const int i, const int j) const override
 Returns the standard lement edge IDs that are normal to a given face vertex. More...
 
int v_GetDir (const int faceidx, const int facedir) const override
 Returns the element coordinate direction corresponding to a given face coordinate direction. More...
 
void v_Reset (CurveMap &curvedEdges, CurveMap &curvedFaces) override
 Reset this geometry object: unset the current state, zero Geometry::m_coeffs and remove allocated GeomFactors. More...
 
void v_Setup () override
 
- Protected Member Functions inherited from Nektar::SpatialDomains::Geometry3D
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)
 
void v_FillGeom () override
 Put all quadrature information into face/edge structure and backward transform. More...
 
NekDouble v_GetCoord (const int i, const Array< OneD, const NekDouble > &Lcoord) override
 Given local collapsed coordinate Lcoord return the value of physical coordinate in direction i. More...
 
void v_CalculateInverseIsoParam () override
 
int v_AllLeftCheck (const Array< OneD, const NekDouble > &gloCoord) override
 
int v_GetShapeDim () const override
 Get the object's shape dimension. More...
 
int v_GetNumVerts () const override
 Get the number of vertices of this object. More...
 
int v_GetNumEdges () const override
 Get the number of edges of this object. More...
 
int v_GetNumFaces () const override
 Get the number of faces of this object. More...
 
PointGeomSharedPtr v_GetVertex (int i) const override
 
Geometry1DSharedPtr v_GetEdge (int i) const override
 Returns edge i of this object. More...
 
Geometry2DSharedPtr v_GetFace (int i) const override
 Returns face i of this object. More...
 
StdRegions::Orientation v_GetEorient (const int i) const override
 Returns the orientation of edge i with respect to the ordering of edges in the standard element. More...
 
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
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
 
int m_straightEdge
 
- Static Protected Attributes inherited from Nektar::SpatialDomains::Geometry
static GeomFactorsVector m_regGeomFactorsManager
 

Detailed Description

Definition at line 44 of file PyrGeom.h.

Constructor & Destructor Documentation

◆ PyrGeom() [1/2]

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

Definition at line 50 of file PyrGeom.cpp.

51{
53}
LibUtilities::ShapeType m_shapeType
Type of shape.
Definition: Geometry.h:203

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 55 of file PyrGeom.cpp.

56 : Geometry3D(faces[0]->GetEdge(0)->GetVertex(0)->GetCoordim())
57{
59 m_globalID = id;
60
61 /// Copy the face shared pointers
62 m_faces.insert(m_faces.begin(), faces, faces + PyrGeom::kNfaces);
63
64 /// Set up orientation vectors with correct amount of elements.
65 m_eorient.resize(kNedges);
66 m_forient.resize(kNfaces);
67
72}
std::vector< StdRegions::Orientation > m_forient
Definition: Geometry3D.h:81
std::vector< StdRegions::Orientation > m_eorient
Definition: Geometry3D.h:80
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:55
static const int kNedges
Definition: PyrGeom.h:52

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

Definition at line 74 of file PyrGeom.cpp.

75{
76}

Member Function Documentation

◆ SetUpEdgeOrientation()

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

Definition at line 353 of file PyrGeom.cpp.

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

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

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

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 277 of file PyrGeom.cpp.

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

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 690 of file PyrGeom.cpp.

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

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

Implements Nektar::SpatialDomains::Geometry.

Definition at line 78 of file PyrGeom.cpp.

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

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

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 142 of file PyrGeom.cpp.

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

◆ v_GetEdgeNormalToFaceVert()

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

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 158 of file PyrGeom.cpp.

159{
160 return EdgeNormalToFaceVert[i][j];
161}
static const unsigned int EdgeNormalToFaceVert[5][4]
Definition: PyrGeom.h:72

References EdgeNormalToFaceVert.

◆ v_Reset()

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

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

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 659 of file PyrGeom.cpp.

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

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

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 672 of file PyrGeom.cpp.

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

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 72 of file PyrGeom.h.

Referenced by v_GetEdgeNormalToFaceVert().

◆ kNedges

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

Definition at line 52 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 51 of file PyrGeom.h.

◆ XMLElementType

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

Definition at line 56 of file PyrGeom.h.