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

virtual void v_GenGeomFactors () override
 
virtual 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...
 
virtual 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...
 
virtual 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...
 
virtual 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...
 
virtual int v_GetDir (const int faceidx, const int facedir) const override
 Returns the element coordinate direction corresponding to a given face coordinate direction. More...
 
virtual 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...
 
virtual void v_Setup () override
 
- Protected Member Functions inherited from Nektar::SpatialDomains::Geometry3D
virtual NekDouble v_GetLocCoords (const Array< OneD, const NekDouble > &coords, Array< OneD, NekDouble > &Lcoords) override
 Determine the local collapsed coordinates that correspond to a given Cartesian coordinate for this geometry object. More...
 
void NewtonIterationForLocCoord (const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &ptsx, const Array< OneD, const NekDouble > &ptsy, const Array< OneD, const NekDouble > &ptsz, Array< OneD, NekDouble > &Lcoords, NekDouble &dist)
 
virtual void v_FillGeom () override
 Put all quadrature information into face/edge structure and backward transform. More...
 
virtual NekDouble v_GetCoord (const int i, const Array< OneD, const NekDouble > &Lcoord) override
 Given local collapsed coordinate Lcoord return the value of physical coordinate in direction i. More...
 
virtual int v_GetShapeDim () const override
 Get the object's shape dimension. More...
 
virtual int v_GetNumVerts () const override
 Get the number of vertices of this object. More...
 
virtual int v_GetNumEdges () const override
 Get the number of edges of this object. More...
 
virtual int v_GetNumFaces () const override
 Get the number of faces of this object. More...
 
virtual PointGeomSharedPtr v_GetVertex (int i) const override
 
virtual Geometry1DSharedPtr v_GetEdge (int i) const override
 Returns edge i of this object. More...
 
virtual Geometry2DSharedPtr v_GetFace (int i) const override
 Returns face i of this object. More...
 
virtual StdRegions::Orientation v_GetEorient (const int i) const override
 Returns the orientation of edge i with respect to the ordering of edges in the standard element. More...
 
virtual StdRegions::Orientation v_GetForient (const int i) const override
 Returns the orientation of face i with respect to the ordering of faces in the standard element. More...
 
- Protected Member Functions inherited from Nektar::SpatialDomains::Geometry
void GenGeomFactors ()
 Handles generation of geometry factors. More...
 
virtual 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)
 
virtual NekDouble v_FindDistance (const Array< OneD, const NekDouble > &xs, Array< OneD, NekDouble > &xi)
 
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]
 
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...
 
- 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 58 of file PrismGeom.cpp.

59 {
61 }
LibUtilities::ShapeType m_shapeType
Type of shape.
Definition: Geometry.h:198

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

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

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

84 {
85 }

Member Function Documentation

◆ SetUpEdgeOrientation()

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

Definition at line 398 of file PrismGeom.cpp.

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

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

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

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

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

Referenced by PrismGeom().

◆ SetUpLocalVertices()

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

Definition at line 329 of file PrismGeom.cpp.

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

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

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

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

104 {
105  if (!m_setupState)
106  {
108  }
109 
111  {
112  int i, f;
113  GeomType Gtype = eRegular;
114 
115  v_FillGeom();
116 
117  // check to see if expansions are linear
118  for (i = 0; i < m_coordim; ++i)
119  {
120  if (m_xmap->GetBasisNumModes(0) != 2 ||
121  m_xmap->GetBasisNumModes(1) != 2 ||
122  m_xmap->GetBasisNumModes(2) != 2)
123  {
124  Gtype = eDeformed;
125  }
126  }
127 
128  // check to see if all quadrilateral faces are parallelograms
129  if (Gtype == eRegular)
130  {
131  // Vertex ids of quad faces
132  const unsigned int faceVerts[3][4] = {
133  {0, 1, 2, 3}, {1, 2, 5, 4}, {0, 3, 5, 4}};
134 
135  for (f = 0; f < 3; f++)
136  {
137  // Ensure each face is a parallelogram? Check this.
138  for (i = 0; i < m_coordim; i++)
139  {
140  if (fabs((*m_verts[faceVerts[f][0]])(i) -
141  (*m_verts[faceVerts[f][1]])(i) +
142  (*m_verts[faceVerts[f][2]])(i) -
143  (*m_verts[faceVerts[f][3]])(i)) >
145  {
146  Gtype = eDeformed;
147  break;
148  }
149  }
150 
151  if (Gtype == eDeformed)
152  {
153  break;
154  }
155  }
156  }
157 
159  Gtype, m_coordim, m_xmap, m_coeffs);
161  }
162 }
virtual void v_FillGeom() override
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:194
Array< OneD, Array< OneD, NekDouble > > m_coeffs
Array containing expansion coefficients of m_xmap.
Definition: Geometry.h:202
GeomState m_geomFactorsState
State of the geometric factors.
Definition: Geometry.h:188
GeomFactorsSharedPtr m_geomFactors
Geometric factors.
Definition: Geometry.h:186
virtual void v_Setup() override
Definition: PrismGeom.cpp:720
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
overrideprotectedvirtual

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

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 87 of file PrismGeom.cpp.

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

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

175 {
176  return EdgeFaceConnectivity[i][j];
177 }
static const unsigned int EdgeFaceConnectivity[9][2]
Definition: PrismGeom.h:81

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

180 {
181  return EdgeNormalToFaceVert[i][j];
182 }
static const unsigned int EdgeNormalToFaceVert[5][4]
Definition: PrismGeom.h:82

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

165 {
166  return VertexEdgeConnectivity[i][j];
167 }
static const unsigned int VertexEdgeConnectivity[6][3]
Definition: PrismGeom.h:79

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

170 {
171  return VertexFaceConnectivity[i][j];
172 }
static const unsigned int VertexFaceConnectivity[6][3]
Definition: PrismGeom.h:80

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

711 {
712  Geometry::v_Reset(curvedEdges, curvedFaces);
713 
714  for (int i = 0; i < 5; ++i)
715  {
716  m_faces[i]->Reset(curvedEdges, curvedFaces);
717  }
718 }
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:376

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

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

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

Referenced by v_GetEdgeNormalToFaceVert().

◆ 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 79 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 80 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.