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

#include <PrismGeom.h>

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

Public Member Functions

 PrismGeom ()
 
 PrismGeom (int id, const Geometry2DSharedPtr faces[])
 
 ~PrismGeom () 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 Setup ()
 

Static Public Attributes

static const int kNverts = 6
 
static const int kNedges = 9
 
static const int kNqfaces = 3
 
static const int kNtfaces = 2
 
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_GetVertexEdgeMap (const int i, const int j) const override
 Returns the standard element edge IDs that are connected to a given vertex. More...
 
int v_GetVertexFaceMap (const int i, const int j) const override
 Returns the standard element face IDs that are connected to a given vertex. More...
 
int v_GetEdgeFaceMap (const int i, const int j) const override
 Returns the standard element edge IDs that are connected to a given face. More...
 
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
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 VertexEdgeConnectivity [6][3]
 
static const unsigned int VertexFaceConnectivity [6][3]
 
static const unsigned int EdgeFaceConnectivity [9][2]
 
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 45 of file PrismGeom.h.

Constructor & Destructor Documentation

◆ PrismGeom() [1/2]

Nektar::SpatialDomains::PrismGeom::PrismGeom ( )

Definition at line 56 of file PrismGeom.cpp.

57{
59}
LibUtilities::ShapeType m_shapeType
Type of shape.
Definition: Geometry.h:197

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

◆ PrismGeom() [2/2]

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

Copy the face shared pointers.

Set up orientation vectors with correct amount of elements.

Set up local objects.

Definition at line 61 of file PrismGeom.cpp.

62 : Geometry3D(faces[0]->GetEdge(0)->GetVertex(0)->GetCoordim())
63{
65 m_globalID = id;
66
67 /// Copy the face shared pointers.
68 m_faces.insert(m_faces.begin(), faces, faces + PrismGeom::kNfaces);
69
70 /// Set up orientation vectors with correct amount of elements.
71 m_eorient.resize(kNedges);
72 m_forient.resize(kNfaces);
73
74 /// Set up local objects.
79}
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:355
Geometry1DSharedPtr GetEdge(int i) const
Returns edge i of this object.
Definition: Geometry.h:363
int GetCoordim() const
Return the coordinate dimension of this object (i.e. the dimension of the space in which this object ...
Definition: Geometry.h:281

References Nektar::LibUtilities::ePrism, 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().

◆ ~PrismGeom()

Nektar::SpatialDomains::PrismGeom::~PrismGeom ( )
override

Definition at line 81 of file PrismGeom.cpp.

82{
83}

Member Function Documentation

◆ SetUpEdgeOrientation()

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

Definition at line 404 of file PrismGeom.cpp.

405{
406
407 // This 2D array holds the local id's of all the vertices
408 // for every edge. For every edge, they are ordered to what we
409 // define as being Forwards
410 const unsigned int edgeVerts[kNedges][2] = {
411 {0, 1}, {1, 2}, {3, 2}, {0, 3}, {0, 4}, {1, 4}, {2, 5}, {3, 5}, {4, 5}};
412
413 int i;
414 for (i = 0; i < kNedges; i++)
415 {
416 if (m_edges[i]->GetVid(0) == m_verts[edgeVerts[i][0]]->GetGlobalID())
417 {
419 }
420 else if (m_edges[i]->GetVid(0) ==
421 m_verts[edgeVerts[i][1]]->GetGlobalID())
422 {
424 }
425 else
426 {
427 ASSERTL0(false, "Could not find matching vertex for the edge");
428 }
429 }
430}
#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:324

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 PrismGeom().

◆ SetUpFaceOrientation()

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

Definition at line 432 of file PrismGeom.cpp.

433{
434 int f, i;
435
436 // These arrays represent the vector of the A and B
437 // coordinate of the local elemental coordinate system
438 // where A corresponds with the coordinate direction xi_i
439 // with the lowest index i (for that particular face)
440 // Coordinate 'B' then corresponds to the other local
441 // coordinate (i.e. with the highest index)
442 Array<OneD, NekDouble> elementAaxis(m_coordim);
443 Array<OneD, NekDouble> elementBaxis(m_coordim);
444
445 // These arrays correspond to the local coordinate
446 // system of the face itself (i.e. the Geometry2D)
447 // faceAaxis correspond to the xi_0 axis
448 // faceBaxis correspond to the xi_1 axis
449 Array<OneD, NekDouble> faceAaxis(m_coordim);
450 Array<OneD, NekDouble> faceBaxis(m_coordim);
451
452 // This is the base vertex of the face (i.e. the Geometry2D)
453 // This corresponds to thevertex with local ID 0 of the
454 // Geometry2D
455 unsigned int baseVertex;
456
457 // The lenght of the vectors above
458 NekDouble elementAaxis_length;
459 NekDouble elementBaxis_length;
460 NekDouble faceAaxis_length;
461 NekDouble faceBaxis_length;
462
463 // This 2D array holds the local id's of all the vertices
464 // for every face. For every face, they are ordered in such
465 // a way that the implementation below allows a unified approach
466 // for all faces.
467 const unsigned int faceVerts[kNfaces][QuadGeom::kNverts] = {
468 {0, 1, 2, 3},
469 {0, 1, 4, 0}, // This is triangle requires only three vertices
470 {1, 2, 5, 4},
471 {3, 2, 5, 0}, // This is triangle requires only three vertices
472 {0, 3, 5, 4},
473 };
474
475 NekDouble dotproduct1 = 0.0;
476 NekDouble dotproduct2 = 0.0;
477
478 unsigned int orientation;
479
480 // Loop over all the faces to set up the orientation
481 for (f = 0; f < kNqfaces + kNtfaces; f++)
482 {
483 // initialisation
484 elementAaxis_length = 0.0;
485 elementBaxis_length = 0.0;
486 faceAaxis_length = 0.0;
487 faceBaxis_length = 0.0;
488
489 dotproduct1 = 0.0;
490 dotproduct2 = 0.0;
491
492 baseVertex = m_faces[f]->GetVid(0);
493
494 // We are going to construct the vectors representing the A and B axis
495 // of every face. These vectors will be constructed as a
496 // vector-representation
497 // of the edges of the face. However, for both coordinate directions, we
498 // can
499 // represent the vectors by two different edges. That's why we need to
500 // make sure that
501 // we pick the edge to which the baseVertex of the
502 // Geometry2D-representation of the face
503 // belongs...
504
505 // Compute the length of edges on a base-face
506 if (f == 1 || f == 3)
507 { // Face is a Triangle
508 if (baseVertex == m_verts[faceVerts[f][0]]->GetVid())
509 {
510 for (i = 0; i < m_coordim; i++)
511 {
512 elementAaxis[i] = (*m_verts[faceVerts[f][1]])[i] -
513 (*m_verts[faceVerts[f][0]])[i];
514 elementBaxis[i] = (*m_verts[faceVerts[f][2]])[i] -
515 (*m_verts[faceVerts[f][0]])[i];
516 }
517 }
518 else if (baseVertex == m_verts[faceVerts[f][1]]->GetVid())
519 {
520 for (i = 0; i < m_coordim; i++)
521 {
522 elementAaxis[i] = (*m_verts[faceVerts[f][1]])[i] -
523 (*m_verts[faceVerts[f][0]])[i];
524 elementBaxis[i] = (*m_verts[faceVerts[f][2]])[i] -
525 (*m_verts[faceVerts[f][1]])[i];
526 }
527 }
528 else if (baseVertex == m_verts[faceVerts[f][2]]->GetVid())
529 {
530 for (i = 0; i < m_coordim; i++)
531 {
532 elementAaxis[i] = (*m_verts[faceVerts[f][1]])[i] -
533 (*m_verts[faceVerts[f][2]])[i];
534 elementBaxis[i] = (*m_verts[faceVerts[f][2]])[i] -
535 (*m_verts[faceVerts[f][0]])[i];
536 }
537 }
538 else
539 {
540 ASSERTL0(false, "Could not find matching vertex for the face");
541 }
542 }
543 else
544 { // Face is a Quad
545 if (baseVertex == m_verts[faceVerts[f][0]]->GetGlobalID())
546 {
547 for (i = 0; i < m_coordim; i++)
548 {
549 elementAaxis[i] = (*m_verts[faceVerts[f][1]])[i] -
550 (*m_verts[faceVerts[f][0]])[i];
551 elementBaxis[i] = (*m_verts[faceVerts[f][3]])[i] -
552 (*m_verts[faceVerts[f][0]])[i];
553 }
554 }
555 else if (baseVertex == m_verts[faceVerts[f][1]]->GetGlobalID())
556 {
557 for (i = 0; i < m_coordim; i++)
558 {
559 elementAaxis[i] = (*m_verts[faceVerts[f][1]])[i] -
560 (*m_verts[faceVerts[f][0]])[i];
561 elementBaxis[i] = (*m_verts[faceVerts[f][2]])[i] -
562 (*m_verts[faceVerts[f][1]])[i];
563 }
564 }
565 else if (baseVertex == m_verts[faceVerts[f][2]]->GetGlobalID())
566 {
567 for (i = 0; i < m_coordim; i++)
568 {
569 elementAaxis[i] = (*m_verts[faceVerts[f][2]])[i] -
570 (*m_verts[faceVerts[f][3]])[i];
571 elementBaxis[i] = (*m_verts[faceVerts[f][2]])[i] -
572 (*m_verts[faceVerts[f][1]])[i];
573 }
574 }
575 else if (baseVertex == m_verts[faceVerts[f][3]]->GetGlobalID())
576 {
577 for (i = 0; i < m_coordim; i++)
578 {
579 elementAaxis[i] = (*m_verts[faceVerts[f][2]])[i] -
580 (*m_verts[faceVerts[f][3]])[i];
581 elementBaxis[i] = (*m_verts[faceVerts[f][3]])[i] -
582 (*m_verts[faceVerts[f][0]])[i];
583 }
584 }
585 else
586 {
587 ASSERTL0(false, "Could not find matching vertex for the face");
588 }
589 }
590 // Now, construct the edge-vectors of the local coordinates of
591 // the Geometry2D-representation of the face
592 for (i = 0; i < m_coordim; i++)
593 {
594 int v = m_faces[f]->GetNumVerts() - 1;
595 faceAaxis[i] =
596 (*m_faces[f]->GetVertex(1))[i] - (*m_faces[f]->GetVertex(0))[i];
597 faceBaxis[i] =
598 (*m_faces[f]->GetVertex(v))[i] - (*m_faces[f]->GetVertex(0))[i];
599
600 elementAaxis_length += pow(elementAaxis[i], 2);
601 elementBaxis_length += pow(elementBaxis[i], 2);
602 faceAaxis_length += pow(faceAaxis[i], 2);
603 faceBaxis_length += pow(faceBaxis[i], 2);
604 }
605
606 elementAaxis_length = sqrt(elementAaxis_length);
607 elementBaxis_length = sqrt(elementBaxis_length);
608 faceAaxis_length = sqrt(faceAaxis_length);
609 faceBaxis_length = sqrt(faceBaxis_length);
610
611 // Calculate the inner product of both the A-axis
612 // (i.e. Elemental A axis and face A axis)
613 for (i = 0; i < m_coordim; i++)
614 {
615 dotproduct1 += elementAaxis[i] * faceAaxis[i];
616 }
617
618 orientation = 0;
619
620 // if the innerproduct is equal to the (absolute value of the ) products
621 // of the lengths of both vectors, then, the coordinate systems will NOT
622 // be transposed
623 if (fabs(elementAaxis_length * faceAaxis_length - fabs(dotproduct1)) <
625 {
626 // if the inner product is negative, both A-axis point
627 // in reverse direction
628 if (dotproduct1 < 0.0)
629 {
630 orientation += 2;
631 }
632
633 // calculate the inner product of both B-axis
634 for (i = 0; i < m_coordim; i++)
635 {
636 dotproduct2 += elementBaxis[i] * faceBaxis[i];
637 }
638
639 ASSERTL1(fabs(fabs(dotproduct2 / elementBaxis_length /
640 faceBaxis_length) -
642 "These vectors should be parallel");
643
644 // if the inner product is negative, both B-axis point
645 // in reverse direction
646 if (dotproduct2 < 0.0)
647 {
648 orientation++;
649 }
650 }
651 // The coordinate systems are transposed
652 else
653 {
654 orientation = 4;
655
656 // Calculate the inner product between the elemental A-axis
657 // and the B-axis of the face (which are now the corresponding axis)
658 dotproduct1 = 0.0;
659 for (i = 0; i < m_coordim; i++)
660 {
661 dotproduct1 += elementAaxis[i] * faceBaxis[i];
662 }
663
664 // check that both these axis are indeed parallel
665 ASSERTL1(fabs(fabs(dotproduct1) / elementAaxis_length /
666 faceBaxis_length -
668 "These vectors should be parallel");
669
670 // if the result is negative, both axis point in reverse
671 // directions
672 if (dotproduct1 < 0.0)
673 {
674 orientation += 2;
675 }
676
677 // Do the same for the other two corresponding axis
678 dotproduct2 = 0.0;
679 for (i = 0; i < m_coordim; i++)
680 {
681 dotproduct2 += elementBaxis[i] * faceAaxis[i];
682 }
683
684 ASSERTL1(fabs(fabs(dotproduct2) / elementBaxis_length /
685 faceAaxis_length -
687 "These vectors should be parallel");
688
689 if (dotproduct2 < 0.0)
690 {
691 orientation++;
692 }
693 }
694
695 orientation = orientation + 5;
696
697 if ((f == 1) || (f == 3)) // check triange orientation
698 {
699 ASSERTL0(
701 "Orientation of triangular face (id = " +
702 boost::lexical_cast<string>(m_faces[f]->GetGlobalID()) +
703 ") is inconsistent with face " +
704 boost::lexical_cast<string>(f) +
705 " of prism element (id = " +
706 boost::lexical_cast<string>(m_globalID) +
707 ") since Dir2 is aligned with Dir1. Mesh setup "
708 "needs investigation");
709 }
710
711 // Fill the m_forient array
712 m_forient[f] = (StdRegions::Orientation)orientation;
713 }
714}
#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:183
static const int kNverts
Definition: QuadGeom.h:73
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::QuadGeom::kNverts, 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 PrismGeom().

◆ SetUpLocalEdges()

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

Definition at line 190 of file PrismGeom.cpp.

191{
192 // find edge 0
193 int i, j;
194 unsigned int check;
195
196 SegGeomSharedPtr edge;
197
198 // First set up the 4 bottom edges
199 int f; // Connected face index
200 for (f = 1; f < 5; f++)
201 {
202 int nEdges = m_faces[f]->GetNumEdges();
203 check = 0;
204 for (i = 0; i < 4; i++)
205 {
206 for (j = 0; j < nEdges; j++)
207 {
208 if (m_faces[0]->GetEid(i) == m_faces[f]->GetEid(j))
209 {
210 edge =
211 dynamic_pointer_cast<SegGeom>((m_faces[0])->GetEdge(i));
212 m_edges.push_back(edge);
213 check++;
214 }
215 }
216 }
217
218 if (check < 1)
219 {
220 std::ostringstream errstrm;
221 errstrm << "Connected faces do not share an edge. Faces ";
222 errstrm << (m_faces[0])->GetGlobalID() << ", "
223 << (m_faces[f])->GetGlobalID();
224 ASSERTL0(false, errstrm.str());
225 }
226 else if (check > 1)
227 {
228 std::ostringstream errstrm;
229 errstrm << "Connected faces share more than one edge. Faces ";
230 errstrm << (m_faces[0])->GetGlobalID() << ", "
231 << (m_faces[f])->GetGlobalID();
232 ASSERTL0(false, errstrm.str());
233 }
234 }
235
236 // Then, set up the 4 vertical edges
237 check = 0;
238 for (i = 0; i < 3; i++) // Set up the vertical edge :face(1) and face(4)
239 {
240 for (j = 0; j < 4; j++)
241 {
242 if ((m_faces[1])->GetEid(i) == (m_faces[4])->GetEid(j))
243 {
244 edge = dynamic_pointer_cast<SegGeom>((m_faces[1])->GetEdge(i));
245 m_edges.push_back(edge);
246 check++;
247 }
248 }
249 }
250 if (check < 1)
251 {
252 std::ostringstream errstrm;
253 errstrm << "Connected faces do not share an edge. Faces ";
254 errstrm << (m_faces[1])->GetGlobalID() << ", "
255 << (m_faces[4])->GetGlobalID();
256 ASSERTL0(false, errstrm.str());
257 }
258 else if (check > 1)
259 {
260 std::ostringstream errstrm;
261 errstrm << "Connected faces share more than one edge. Faces ";
262 errstrm << (m_faces[1])->GetGlobalID() << ", "
263 << (m_faces[4])->GetGlobalID();
264 ASSERTL0(false, errstrm.str());
265 }
266 // Set up vertical edges: face(1) through face(4)
267 for (f = 1; f < 4; f++)
268 {
269 check = 0;
270 for (i = 0; i < m_faces[f]->GetNumEdges(); i++)
271 {
272 for (j = 0; j < m_faces[f + 1]->GetNumEdges(); j++)
273 {
274 if ((m_faces[f])->GetEid(i) == (m_faces[f + 1])->GetEid(j))
275 {
276 edge =
277 dynamic_pointer_cast<SegGeom>((m_faces[f])->GetEdge(i));
278 m_edges.push_back(edge);
279 check++;
280 }
281 }
282 }
283
284 if (check < 1)
285 {
286 std::ostringstream errstrm;
287 errstrm << "Connected faces do not share an edge. Faces ";
288 errstrm << (m_faces[f])->GetGlobalID() << ", "
289 << (m_faces[f + 1])->GetGlobalID();
290 ASSERTL0(false, errstrm.str());
291 }
292 else if (check > 1)
293 {
294 std::ostringstream errstrm;
295 errstrm << "Connected faces share more than one edge. Faces ";
296 errstrm << (m_faces[f])->GetGlobalID() << ", "
297 << (m_faces[f + 1])->GetGlobalID();
298 ASSERTL0(false, errstrm.str());
299 }
300 }
301
302 // Finally, set up the 1 top edge
303 check = 0;
304 for (i = 0; i < 4; i++)
305 {
306 for (j = 0; j < 4; j++)
307 {
308 if ((m_faces[2])->GetEid(i) == (m_faces[4])->GetEid(j))
309 {
310 edge = dynamic_pointer_cast<SegGeom>((m_faces[2])->GetEdge(i));
311 m_edges.push_back(edge);
312 check++;
313 }
314 }
315 }
316
317 if (check < 1)
318 {
319 std::ostringstream errstrm;
320 errstrm << "Connected faces do not share an edge. Faces ";
321 errstrm << (m_faces[1])->GetGlobalID() << ", "
322 << (m_faces[3])->GetGlobalID();
323 ASSERTL0(false, errstrm.str());
324 }
325 else if (check > 1)
326 {
327 std::ostringstream errstrm;
328 errstrm << "Connected faces share more than one edge. Faces ";
329 errstrm << (m_faces[1])->GetGlobalID() << ", "
330 << (m_faces[3])->GetGlobalID();
331 ASSERTL0(false, errstrm.str());
332 }
333}
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 PrismGeom().

◆ SetUpLocalVertices()

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

Definition at line 335 of file PrismGeom.cpp.

336{
337
338 // Set up the first 2 vertices (i.e. vertex 0,1)
339 if ((m_edges[0]->GetVid(0) == m_edges[1]->GetVid(0)) ||
340 (m_edges[0]->GetVid(0) == m_edges[1]->GetVid(1)))
341 {
342 m_verts.push_back(m_edges[0]->GetVertex(1));
343 m_verts.push_back(m_edges[0]->GetVertex(0));
344 }
345 else if ((m_edges[0]->GetVid(1) == m_edges[1]->GetVid(0)) ||
346 (m_edges[0]->GetVid(1) == m_edges[1]->GetVid(1)))
347 {
348 m_verts.push_back(m_edges[0]->GetVertex(0));
349 m_verts.push_back(m_edges[0]->GetVertex(1));
350 }
351 else
352 {
353 std::ostringstream errstrm;
354 errstrm << "Connected edges do not share a vertex. Edges ";
355 errstrm << m_edges[0]->GetGlobalID() << ", "
356 << m_edges[1]->GetGlobalID();
357 ASSERTL0(false, errstrm.str());
358 }
359
360 // set up the other bottom vertices (i.e. vertex 2,3)
361 for (int i = 1; i < 3; i++)
362 {
363 if (m_edges[i]->GetVid(0) == m_verts[i]->GetGlobalID())
364 {
365 m_verts.push_back(m_edges[i]->GetVertex(1));
366 }
367 else if (m_edges[i]->GetVid(1) == m_verts[i]->GetGlobalID())
368 {
369 m_verts.push_back(m_edges[i]->GetVertex(0));
370 }
371 else
372 {
373 std::ostringstream errstrm;
374 errstrm << "Connected edges do not share a vertex. Edges ";
375 errstrm << m_edges[i]->GetGlobalID() << ", "
376 << m_edges[i - 1]->GetGlobalID();
377 ASSERTL0(false, errstrm.str());
378 }
379 }
380
381 // set up top vertices
382 // First, set up vertices 4,5
383 if ((m_edges[8]->GetVid(0) == m_edges[4]->GetVid(0)) ||
384 (m_edges[8]->GetVid(0) == m_edges[4]->GetVid(1)))
385 {
386 m_verts.push_back(m_edges[8]->GetVertex(0));
387 m_verts.push_back(m_edges[8]->GetVertex(1));
388 }
389 else if ((m_edges[8]->GetVid(1) == m_edges[4]->GetVid(0)) ||
390 (m_edges[8]->GetVid(1) == m_edges[4]->GetVid(1)))
391 {
392 m_verts.push_back(m_edges[8]->GetVertex(1));
393 m_verts.push_back(m_edges[8]->GetVertex(0));
394 }
395 else
396 {
397 std::ostringstream errstrm;
398 errstrm << "Connected edges do not share a vertex. Edges ";
399 errstrm << m_edges[8]->GetGlobalID();
400 ASSERTL0(false, errstrm.str());
401 }
402}

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 PrismGeom().

◆ SetUpXmap()

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

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

Definition at line 744 of file PrismGeom.cpp.

745{
746 vector<int> tmp;
747 int order0, order1;
748
749 if (m_forient[0] < 9)
750 {
751 tmp.push_back(m_faces[0]->GetXmap()->GetTraceNcoeffs(0));
752 tmp.push_back(m_faces[0]->GetXmap()->GetTraceNcoeffs(2));
753 order0 = *max_element(tmp.begin(), tmp.end());
754 }
755 else
756 {
757 tmp.push_back(m_faces[0]->GetXmap()->GetTraceNcoeffs(1));
758 tmp.push_back(m_faces[0]->GetXmap()->GetTraceNcoeffs(3));
759 order0 = *max_element(tmp.begin(), tmp.end());
760 }
761
762 if (m_forient[0] < 9)
763 {
764 tmp.clear();
765 tmp.push_back(m_faces[0]->GetXmap()->GetTraceNcoeffs(1));
766 tmp.push_back(m_faces[0]->GetXmap()->GetTraceNcoeffs(3));
767 tmp.push_back(m_faces[2]->GetXmap()->GetTraceNcoeffs(2));
768 order1 = *max_element(tmp.begin(), tmp.end());
769 }
770 else
771 {
772 tmp.clear();
773 tmp.push_back(m_faces[0]->GetXmap()->GetTraceNcoeffs(0));
774 tmp.push_back(m_faces[0]->GetXmap()->GetTraceNcoeffs(2));
775 tmp.push_back(m_faces[2]->GetXmap()->GetTraceNcoeffs(2));
776 order1 = *max_element(tmp.begin(), tmp.end());
777 }
778
779 tmp.clear();
780 tmp.push_back(order0);
781 tmp.push_back(order1);
782 tmp.push_back(m_faces[1]->GetXmap()->GetTraceNcoeffs(1));
783 tmp.push_back(m_faces[1]->GetXmap()->GetTraceNcoeffs(2));
784 tmp.push_back(m_faces[3]->GetXmap()->GetTraceNcoeffs(1));
785 tmp.push_back(m_faces[3]->GetXmap()->GetTraceNcoeffs(2));
786 int order2 = *max_element(tmp.begin(), tmp.end());
787
788 const LibUtilities::BasisKey A(
790 LibUtilities::PointsKey(order0 + 1,
792 const LibUtilities::BasisKey B(
794 LibUtilities::PointsKey(order1 + 1,
796 const LibUtilities::BasisKey C(
798 LibUtilities::PointsKey(order2, LibUtilities::eGaussRadauMAlpha1Beta0));
799
801}
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:189
StdRegions::StdExpansionSharedPtr GetXmap() const
Return the mapping object Geometry::m_xmap that represents the coordinate transformation from standar...
Definition: Geometry.h:433
@ eGaussLobattoLegendre
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:51
@ eModified_B
Principle Modified Functions .
Definition: BasisType.h:49
@ eModified_A
Principle Modified Functions .
Definition: BasisType.h:48

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

Referenced by v_Setup().

◆ v_GenGeomFactors()

void Nektar::SpatialDomains::PrismGeom::v_GenGeomFactors ( )
overrideprotectedvirtual

Implements Nektar::SpatialDomains::Geometry.

Definition at line 101 of file PrismGeom.cpp.

102{
103 if (!m_setupState)
104 {
106 }
107
109 {
110 GeomType Gtype = eRegular;
111
112 v_FillGeom();
113
114 // check to see if expansions are linear
115 m_straightEdge = 1;
116 if (m_xmap->GetBasisNumModes(0) != 2 ||
117 m_xmap->GetBasisNumModes(1) != 2 ||
118 m_xmap->GetBasisNumModes(2) != 2)
119 {
120 Gtype = eDeformed;
121 m_straightEdge = 0;
122 }
123
124 // check to see if all quadrilateral faces are parallelograms
125 if (Gtype == eRegular)
126 {
127 m_isoParameter = Array<OneD, Array<OneD, NekDouble>>(3);
128 for (int i = 0; i < 3; ++i)
129 {
130 m_isoParameter[i] = Array<OneD, NekDouble>(6, 0.);
131 NekDouble A = (*m_verts[0])(i);
132 NekDouble B = (*m_verts[1])(i);
133 NekDouble C = (*m_verts[2])(i);
134 NekDouble D = (*m_verts[3])(i);
135 NekDouble E = (*m_verts[4])(i);
136 NekDouble F = (*m_verts[5])(i);
137 m_isoParameter[i][0] = 0.25 * (B + C + E + F);
138
139 m_isoParameter[i][1] = 0.25 * (-A + B + C - D); // xi1
140 m_isoParameter[i][2] = 0.25 * (-B + C - E + F); // xi2
141 m_isoParameter[i][3] = 0.25 * (-A - D + E + F); // xi3
142
143 m_isoParameter[i][4] = 0.25 * (A - B + C - D); // xi1*xi2
144 m_isoParameter[i][5] = 0.25 * (A - D - E + F); // xi2*xi3
145 NekDouble tmp = fabs(m_isoParameter[i][1]) +
146 fabs(m_isoParameter[i][2]) +
147 fabs(m_isoParameter[i][3]);
149 for (int d = 4; d < 6; ++d)
150 {
151 if (fabs(m_isoParameter[i][d]) > tmp)
152 {
153 Gtype = eDeformed;
154 }
155 }
156 }
157 }
158
159 if (Gtype == eRegular)
160 {
162 }
163
165 Gtype, m_coordim, m_xmap, m_coeffs);
167 }
168}
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:193
Array< OneD, Array< OneD, NekDouble > > m_isoParameter
Definition: Geometry.h:204
Array< OneD, Array< OneD, NekDouble > > m_coeffs
Array containing expansion coefficients of m_xmap.
Definition: Geometry.h:201
GeomState m_geomFactorsState
State of the geometric factors.
Definition: Geometry.h:187
GeomFactorsSharedPtr m_geomFactors
Geometric factors.
Definition: Geometry.h:185
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.
std::vector< double > d(NPUPPER *NPUPPER)

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::UnitTests::d(), Nektar::SpatialDomains::eDeformed, Nektar::SpatialDomains::ePtsFilled, Nektar::SpatialDomains::eRegular, Nektar::NekConstants::kNekZeroTol, Nektar::SpatialDomains::Geometry::m_coeffs, Nektar::SpatialDomains::Geometry::m_coordim, Nektar::SpatialDomains::Geometry::m_geomFactors, Nektar::SpatialDomains::Geometry::m_geomFactorsState, Nektar::SpatialDomains::Geometry::m_isoParameter, Nektar::SpatialDomains::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::PrismGeom::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 85 of file PrismGeom.cpp.

86{
87 if (faceidx == 0)
88 {
89 return facedir;
90 }
91 else if (faceidx == 1 || faceidx == 3)
92 {
93 return 2 * facedir;
94 }
95 else
96 {
97 return 1 + facedir;
98 }
99}

◆ v_GetEdgeFaceMap()

int Nektar::SpatialDomains::PrismGeom::v_GetEdgeFaceMap ( const int  i,
const int  j 
) const
overrideprotectedvirtual

Returns the standard element edge IDs that are connected to a given face.

For example, on a prism, edge 0 is connnected to faces 0 and 1; GetEdgeFaceMap(0,j) would therefore return the values 0 and 1 respectively. We assume that j runs between 0 and 1 inclusive, since every face is connected to precisely two faces for all 3D elements.

This function is used in the construction of the low-energy preconditioner.

Parameters
iThe edge to query connectivity for.
jThe local face index between 0 and 1 connected to this element.
See also
MultiRegions::PreconditionerLowEnergy

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 180 of file PrismGeom.cpp.

181{
182 return EdgeFaceConnectivity[i][j];
183}
static const unsigned int EdgeFaceConnectivity[9][2]
Definition: PrismGeom.h:78

References EdgeFaceConnectivity.

◆ v_GetEdgeNormalToFaceVert()

int Nektar::SpatialDomains::PrismGeom::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 185 of file PrismGeom.cpp.

186{
187 return EdgeNormalToFaceVert[i][j];
188}
static const unsigned int EdgeNormalToFaceVert[5][4]
Definition: PrismGeom.h:79

References EdgeNormalToFaceVert.

◆ v_GetVertexEdgeMap()

int Nektar::SpatialDomains::PrismGeom::v_GetVertexEdgeMap ( const int  i,
const int  j 
) const
overrideprotectedvirtual

Returns the standard element edge IDs that are connected to a given vertex.

For example, on a prism, vertex 0 is connnected to edges 0, 3, and 4; GetVertexEdgeMap(0,j) would therefore return the values 0, 1 and 4 respectively. We assume that j runs between 0 and 2 inclusive, which is true for every 3D element asides from the pyramid.

This function is used in the construction of the low-energy preconditioner.

Parameters
iThe vertex to query connectivity for.
jThe local edge index between 0 and 2 connected to this element.
Todo:
Expand to work with pyramid elements.
See also
MultiRegions::PreconditionerLowEnergy

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 170 of file PrismGeom.cpp.

171{
172 return VertexEdgeConnectivity[i][j];
173}
static const unsigned int VertexEdgeConnectivity[6][3]
Definition: PrismGeom.h:76

References VertexEdgeConnectivity.

◆ v_GetVertexFaceMap()

int Nektar::SpatialDomains::PrismGeom::v_GetVertexFaceMap ( const int  i,
const int  j 
) const
overrideprotectedvirtual

Returns the standard element face IDs that are connected to a given vertex.

For example, on a hexahedron, vertex 0 is connnected to faces 0, 1, and 4; GetVertexFaceMap(0,j) would therefore return the values 0, 1 and 4 respectively. We assume that j runs between 0 and 2 inclusive, which is true for every 3D element asides from the pyramid.

This is used in the construction of the low-energy preconditioner.

Parameters
iThe vertex to query connectivity for.
jThe local face index between 0 and 2 connected to this element.
Todo:
Expand to work with pyramid elements.
See also
MultiRegions::PreconditionerLowEnergy

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 175 of file PrismGeom.cpp.

176{
177 return VertexFaceConnectivity[i][j];
178}
static const unsigned int VertexFaceConnectivity[6][3]
Definition: PrismGeom.h:77

References VertexFaceConnectivity.

◆ v_Reset()

void Nektar::SpatialDomains::PrismGeom::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 716 of file PrismGeom.cpp.

717{
718 Geometry::v_Reset(curvedEdges, curvedFaces);
719
720 for (int i = 0; i < 5; ++i)
721 {
722 m_faces[i]->Reset(curvedEdges, curvedFaces);
723 }
724}
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:384

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

◆ v_Setup()

void Nektar::SpatialDomains::PrismGeom::v_Setup ( )
overrideprotectedvirtual

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 726 of file PrismGeom.cpp.

727{
728 if (!m_setupState)
729 {
730 for (int i = 0; i < 5; ++i)
731 {
732 m_faces[i]->Setup();
733 }
734 SetUpXmap();
735 SetUpCoeffs(m_xmap->GetNcoeffs());
736 m_setupState = true;
737 }
738}
void SetUpCoeffs(const int nCoeffs)
Initialise the Geometry::m_coeffs array.
Definition: Geometry.h:687
void SetUpXmap()
Set up the m_xmap object by determining the order of each direction from derived faces.
Definition: PrismGeom.cpp:744

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

◆ EdgeFaceConnectivity

const unsigned int Nektar::SpatialDomains::PrismGeom::EdgeFaceConnectivity
staticprivate
Initial value:
= {
{0, 1}, {0, 2}, {0, 3}, {0, 4}, {1, 4}, {1, 2}, {2, 3}, {3, 4}, {2, 4}}

Definition at line 78 of file PrismGeom.h.

Referenced by v_GetEdgeFaceMap().

◆ EdgeNormalToFaceVert

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

Definition at line 79 of file PrismGeom.h.

Referenced by v_GetEdgeNormalToFaceVert().

◆ kNedges

const int Nektar::SpatialDomains::PrismGeom::kNedges = 9
static

Definition at line 53 of file PrismGeom.h.

Referenced by PrismGeom(), and SetUpEdgeOrientation().

◆ kNfaces

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

◆ kNqfaces

const int Nektar::SpatialDomains::PrismGeom::kNqfaces = 3
static

◆ kNtfaces

const int Nektar::SpatialDomains::PrismGeom::kNtfaces = 2
static

◆ kNverts

const int Nektar::SpatialDomains::PrismGeom::kNverts = 6
static

Definition at line 52 of file PrismGeom.h.

◆ VertexEdgeConnectivity

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

Definition at line 76 of file PrismGeom.h.

Referenced by v_GetVertexEdgeMap().

◆ VertexFaceConnectivity

const unsigned int Nektar::SpatialDomains::PrismGeom::VertexFaceConnectivity
staticprivate
Initial value:
= {
{0, 1, 4}, {0, 1, 2}, {0, 2, 3}, {0, 3, 4}, {1, 2, 4}, {2, 3, 4}}

Definition at line 77 of file PrismGeom.h.

Referenced by v_GetVertexFaceMap().

◆ XMLElementType

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

Definition at line 57 of file PrismGeom.h.