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 ()
 
- Public Member Functions inherited from Nektar::SpatialDomains::Geometry3D
 Geometry3D ()
 
 Geometry3D (const int coordim)
 
virtual ~Geometry3D ()
 
- Public Member Functions inherited from Nektar::SpatialDomains::Geometry
 Geometry ()
 Default constructor. More...
 
 Geometry (int coordim)
 Constructor when supplied a coordinate dimension. More...
 
virtual ~Geometry ()
 Default destructor. More...
 
int GetCoordim () const
 Return the coordinate dimension of this object (i.e. the dimension of the space in which this object is embedded). More...
 
void SetCoordim (int coordim)
 Sets the coordinate dimension of this object (i.e. the dimension of the space in which this object is embedded). More...
 
GeomFactorsSharedPtr GetGeomFactors ()
 Get the geometric factors for this object, generating them if required. More...
 
GeomFactorsSharedPtr GetRefGeomFactors (const Array< OneD, const LibUtilities::BasisSharedPtr > &tbasis)
 
GeomFactorsSharedPtr GetMetricInfo ()
 Get the geometric factors for this object. More...
 
LibUtilities::ShapeType GetShapeType (void)
 Get the geometric shape type of this object. More...
 
int GetGlobalID (void) const
 Get the ID of this object. More...
 
void SetGlobalID (int globalid)
 Set the ID of this object. More...
 
int GetVid (int i) const
 Get the ID of vertex i of this object. More...
 
int GetEid (int i) const
 Get the ID of edge i of this object. More...
 
int GetFid (int i) const
 Get the ID of face i of this object. More...
 
int GetTid (int i) const
 Get the ID of trace i of this object. More...
 
PointGeomSharedPtr GetVertex (int i) const
 Returns vertex i of this object. More...
 
Geometry1DSharedPtr GetEdge (int i) const
 Returns edge i of this object. More...
 
Geometry2DSharedPtr GetFace (int i) const
 Returns face i of this object. More...
 
StdRegions::Orientation GetEorient (const int i) const
 Returns the orientation of edge i with respect to the ordering of edges in the standard element. More...
 
StdRegions::Orientation GetForient (const int i) const
 Returns the orientation of face i with respect to the ordering of faces in the standard element. More...
 
int GetNumVerts () const
 Get the number of vertices of this object. More...
 
int GetNumEdges () const
 Get the number of edges of this object. More...
 
int GetNumFaces () const
 Get the number of faces of this object. More...
 
int GetShapeDim () const
 Get the object's shape dimension. More...
 
StdRegions::StdExpansionSharedPtr GetXmap () const
 Return the mapping object Geometry::m_xmap that represents the coordinate transformation from standard element to physical element. More...
 
const Array< OneD, const NekDouble > & GetCoeffs (const int i) const
 Return the coefficients of the transformation Geometry::m_xmap in coordinate direction i. More...
 
void FillGeom ()
 Populate the coordinate mapping Geometry::m_coeffs information from any children geometry elements. More...
 
std::array< NekDouble, 6 > GetBoundingBox ()
 Generates the bounding box for the element. More...
 
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...
 
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)
 Clamp local coords to be within standard regions [-1, 1]^dim. More...
 
int GetVertexEdgeMap (int i, int j) const
 Returns the standard element edge IDs that are connected to a given vertex. More...
 
int GetVertexFaceMap (int i, int j) const
 Returns the standard element face IDs that are connected to a given vertex. More...
 
int GetEdgeFaceMap (int i, int j) const
 Returns the standard element edge IDs that are connected to a given face. More...
 
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

virtual void v_GenGeomFactors ()
 
virtual int v_GetVertexEdgeMap (const int i, const int j) const
 Returns the standard element edge IDs that are connected to a given vertex. More...
 
virtual int v_GetVertexFaceMap (const int i, const int j) const
 Returns the standard element face IDs that are connected to a given vertex. More...
 
virtual int v_GetEdgeFaceMap (const int i, const int j) const
 Returns the standard element edge IDs that are connected to a given face. More...
 
virtual int v_GetDir (const int faceidx, const int facedir) const
 Returns the element coordinate direction corresponding to a given face coordinate direction. More...
 
virtual void v_Reset (CurveMap &curvedEdges, CurveMap &curvedFaces)
 Reset this geometry object: unset the current state, zero Geometry::m_coeffs and remove allocated GeomFactors. More...
 
virtual void v_Setup ()
 
- Protected Member Functions inherited from Nektar::SpatialDomains::Geometry3D
virtual NekDouble v_GetLocCoords (const Array< OneD, const NekDouble > &coords, Array< OneD, NekDouble > &Lcoords)
 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)
 
virtual void v_FillGeom ()
 Put all quadrature information into face/edge structure and backward transform. More...
 
virtual NekDouble v_GetCoord (const int i, const Array< OneD, const NekDouble > &Lcoord)
 Given local collapsed coordinate Lcoord return the value of physical coordinate in direction i. More...
 
virtual int v_GetShapeDim () const
 Get the object's shape dimension. More...
 
virtual int v_GetNumVerts () const
 Get the number of vertices of this object. More...
 
virtual int v_GetNumEdges () const
 Get the number of edges of this object. More...
 
virtual int v_GetNumFaces () const
 Get the number of faces of this object. More...
 
virtual PointGeomSharedPtr v_GetVertex (int i) const
 
virtual Geometry1DSharedPtr v_GetEdge (int i) const
 Returns edge i of this object. More...
 
virtual Geometry2DSharedPtr v_GetFace (int i) const
 Returns face i of this object. More...
 
virtual StdRegions::Orientation v_GetEorient (const int i) const
 Returns the orientation of edge i with respect to the ordering of edges in the standard element. More...
 
virtual StdRegions::Orientation v_GetForient (const int i) const
 Returns the orientation of face i with respect to the ordering of faces in the standard element. More...
 
- Protected Member Functions inherited from Nektar::SpatialDomains::Geometry
void GenGeomFactors ()
 Handles generation of geometry factors. More...
 
virtual StdRegions::StdExpansionSharedPtr v_GetXmap () const
 Return the mapping object Geometry::m_xmap that represents the coordinate transformation from standard element to physical element. More...
 
virtual bool v_ContainsPoint (const Array< OneD, const NekDouble > &gloCoord, Array< OneD, NekDouble > &locCoord, NekDouble tol, NekDouble &dist)
 
void SetUpCoeffs (const int nCoeffs)
 Initialise the Geometry::m_coeffs array. More...
 

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]
 

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

Detailed Description

Definition at line 47 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:199

References Nektar::LibUtilities::ePrism.

◆ 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:84
std::vector< StdRegions::Orientation > m_eorient
Definition: Geometry3D.h:83
PointGeomSharedPtr GetVertex(int i) const
Returns vertex i of this object.
Definition: Geometry.h:348
Geometry1DSharedPtr GetEdge(int i) const
Returns edge i of this object.
Definition: Geometry.h:356
int GetCoordim() const
Return the coordinate dimension of this object (i.e. the dimension of the space in which this object ...
Definition: Geometry.h:276

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

Definition at line 81 of file PrismGeom.cpp.

82 {
83 }

Member Function Documentation

◆ SetUpEdgeOrientation()

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

Definition at line 393 of file PrismGeom.cpp.

394 {
395 
396  // This 2D array holds the local id's of all the vertices
397  // for every edge. For every edge, they are ordered to what we
398  // define as being Forwards
399  const unsigned int edgeVerts[kNedges][2] = {
400  {0, 1}, {1, 2}, {3, 2}, {0, 3}, {0, 4}, {1, 4}, {2, 5}, {3, 5}, {4, 5}};
401 
402  int i;
403  for (i = 0; i < kNedges; i++)
404  {
405  if (m_edges[i]->GetVid(0) == m_verts[edgeVerts[i][0]]->GetGlobalID())
406  {
408  }
409  else if (m_edges[i]->GetVid(0) == m_verts[edgeVerts[i][1]]->GetGlobalID())
410  {
412  }
413  else
414  {
415  ASSERTL0(false, "Could not find matching vertex for the edge");
416  }
417  }
418 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
int GetVid(int i) const
Get the ID of vertex i of this object.
Definition: Geometry.cpp:137
int GetGlobalID(void) const
Get the ID of this object.
Definition: Geometry.h:319

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 420 of file PrismGeom.cpp.

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

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 177 of file PrismGeom.cpp.

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

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

Referenced by PrismGeom().

◆ SetUpLocalVertices()

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

Definition at line 324 of file PrismGeom.cpp.

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

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 731 of file PrismGeom.cpp.

732 {
733  vector<int> tmp;
734  int order0, order1;
735 
736  if (m_forient[0] < 9)
737  {
738  tmp.push_back(m_faces[0]->GetXmap()->GetTraceNcoeffs(0));
739  tmp.push_back(m_faces[0]->GetXmap()->GetTraceNcoeffs(2));
740  order0 = *max_element(tmp.begin(), tmp.end());
741  }
742  else
743  {
744  tmp.push_back(m_faces[0]->GetXmap()->GetTraceNcoeffs(1));
745  tmp.push_back(m_faces[0]->GetXmap()->GetTraceNcoeffs(3));
746  order0 = *max_element(tmp.begin(), tmp.end());
747  }
748 
749  if (m_forient[0] < 9)
750  {
751  tmp.clear();
752  tmp.push_back(m_faces[0]->GetXmap()->GetTraceNcoeffs(1));
753  tmp.push_back(m_faces[0]->GetXmap()->GetTraceNcoeffs(3));
754  tmp.push_back(m_faces[2]->GetXmap()->GetTraceNcoeffs(2));
755  order1 = *max_element(tmp.begin(), tmp.end());
756  }
757  else
758  {
759  tmp.clear();
760  tmp.push_back(m_faces[0]->GetXmap()->GetTraceNcoeffs(0));
761  tmp.push_back(m_faces[0]->GetXmap()->GetTraceNcoeffs(2));
762  tmp.push_back(m_faces[2]->GetXmap()->GetTraceNcoeffs(2));
763  order1 = *max_element(tmp.begin(), tmp.end());
764  }
765 
766  tmp.clear();
767  tmp.push_back(order0);
768  tmp.push_back(order1);
769  tmp.push_back(m_faces[1]->GetXmap()->GetTraceNcoeffs(1));
770  tmp.push_back(m_faces[1]->GetXmap()->GetTraceNcoeffs(2));
771  tmp.push_back(m_faces[3]->GetXmap()->GetTraceNcoeffs(1));
772  tmp.push_back(m_faces[3]->GetXmap()->GetTraceNcoeffs(2));
773  int order2 = *max_element(tmp.begin(), tmp.end());
774 
775  const LibUtilities::BasisKey A(
777  order0,
778  LibUtilities::PointsKey(order0+1, LibUtilities::eGaussLobattoLegendre));
779  const LibUtilities::BasisKey B(
781  order1,
782  LibUtilities::PointsKey(order1+1, LibUtilities::eGaussLobattoLegendre));
783  const LibUtilities::BasisKey C(
785  order2,
786  LibUtilities::PointsKey(order2,
788 
790 }
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:191
StdRegions::StdExpansionSharedPtr GetXmap() const
Return the mapping object Geometry::m_xmap that represents the coordinate transformation from standar...
Definition: Geometry.h:426
@ eGaussLobattoLegendre
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:51
@ eGaussRadauMAlpha1Beta0
Gauss Radau pinned at x=-1, .
Definition: PointsType.h:58
@ 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::eGaussRadauMAlpha1Beta0, 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 ( )
protectedvirtual

Implements Nektar::SpatialDomains::Geometry.

Definition at line 101 of file PrismGeom.cpp.

102 {
103  if(!m_setupState)
104  {
106  }
107 
109  {
110  int i, f;
111  GeomType Gtype = eRegular;
112 
113  v_FillGeom();
114 
115  // check to see if expansions are linear
116  for (i = 0; i < m_coordim; ++i)
117  {
118  if (m_xmap->GetBasisNumModes(0) != 2 ||
119  m_xmap->GetBasisNumModes(1) != 2 ||
120  m_xmap->GetBasisNumModes(2) != 2)
121  {
122  Gtype = eDeformed;
123  }
124  }
125 
126  // check to see if all quadrilateral faces are parallelograms
127  if (Gtype == eRegular)
128  {
129  // Vertex ids of quad faces
130  const unsigned int faceVerts[3][4] = {
131  {0, 1, 2, 3}, {1, 2, 5, 4}, {0, 3, 5, 4}};
132 
133  for (f = 0; f < 3; f++)
134  {
135  // Ensure each face is a parallelogram? Check this.
136  for (i = 0; i < m_coordim; i++)
137  {
138  if (fabs((*m_verts[faceVerts[f][0]])(i) -
139  (*m_verts[faceVerts[f][1]])(i) +
140  (*m_verts[faceVerts[f][2]])(i) -
141  (*m_verts[faceVerts[f][3]])(i)) >
143  {
144  Gtype = eDeformed;
145  break;
146  }
147  }
148 
149  if (Gtype == eDeformed)
150  {
151  break;
152  }
153  }
154  }
155 
157  Gtype, m_coordim, m_xmap, m_coeffs);
159  }
160 }
virtual void v_FillGeom()
Put all quadrature information into face/edge structure and backward transform.
Definition: Geometry3D.cpp:365
bool m_setupState
Wether or not the setup routines have been run.
Definition: Geometry.h:195
Array< OneD, Array< OneD, NekDouble > > m_coeffs
Array containing expansion coefficients of m_xmap.
Definition: Geometry.h:203
GeomState m_geomFactorsState
State of the geometric factors.
Definition: Geometry.h:189
GeomFactorsSharedPtr m_geomFactors
Geometric factors.
Definition: Geometry.h:187
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_setupState, Nektar::SpatialDomains::Geometry3D::m_verts, Nektar::SpatialDomains::Geometry::m_xmap, Nektar::SpatialDomains::Geometry3D::v_FillGeom(), and v_Setup().

◆ v_GetDir()

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

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

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 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
protectedvirtual

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 172 of file PrismGeom.cpp.

173 {
174  return EdgeFaceConnectivity[i][j];
175 }
static const unsigned int EdgeFaceConnectivity[9][2]
Definition: PrismGeom.h:79

References EdgeFaceConnectivity.

◆ v_GetVertexEdgeMap()

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

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 162 of file PrismGeom.cpp.

163 {
164  return VertexEdgeConnectivity[i][j];
165 }
static const unsigned int VertexEdgeConnectivity[6][3]
Definition: PrismGeom.h:77

References VertexEdgeConnectivity.

◆ v_GetVertexFaceMap()

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

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 167 of file PrismGeom.cpp.

168 {
169  return VertexFaceConnectivity[i][j];
170 }
static const unsigned int VertexFaceConnectivity[6][3]
Definition: PrismGeom.h:78

References VertexFaceConnectivity.

◆ v_Reset()

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

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

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 703 of file PrismGeom.cpp.

704 {
705  Geometry::v_Reset(curvedEdges, curvedFaces);
706 
707  for (int i = 0; i < 5; ++i)
708  {
709  m_faces[i]->Reset(curvedEdges, curvedFaces);
710  }
711 }
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:355

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

◆ v_Setup()

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

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 713 of file PrismGeom.cpp.

714 {
715  if(!m_setupState)
716  {
717  for (int i = 0; i < 5; ++i)
718  {
719  m_faces[i]->Setup();
720  }
721  SetUpXmap();
722  SetUpCoeffs(m_xmap->GetNcoeffs());
723  m_setupState = true;
724  }
725 }
void SetUpCoeffs(const int nCoeffs)
Initialise the Geometry::m_coeffs array.
Definition: Geometry.h:658
void SetUpXmap()
Set up the m_xmap object by determining the order of each direction from derived faces.
Definition: PrismGeom.cpp:731

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 79 of file PrismGeom.h.

Referenced by v_GetEdgeFaceMap().

◆ kNedges

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

Definition at line 55 of file PrismGeom.h.

Referenced by PrismGeom(), and SetUpEdgeOrientation().

◆ kNfaces

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

Definition at line 58 of file PrismGeom.h.

Referenced by PrismGeom(), and SetUpFaceOrientation().

◆ kNqfaces

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

Definition at line 56 of file PrismGeom.h.

Referenced by SetUpFaceOrientation().

◆ kNtfaces

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

Definition at line 57 of file PrismGeom.h.

Referenced by SetUpFaceOrientation().

◆ kNverts

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

Definition at line 54 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 77 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 78 of file PrismGeom.h.

Referenced by v_GetVertexFaceMap().

◆ XMLElementType

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

Definition at line 59 of file PrismGeom.h.