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

#include <HexGeom.h>

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

Public Member Functions

 HexGeom ()
 
 HexGeom (int id, const QuadGeomSharedPtr faces[])
 
 ~HexGeom ()
 
- 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 = 8
 
static const int kNedges = 12
 
static const int kNqfaces = 6
 
static const int kNtfaces = 0
 
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 [8][3]
 
static const unsigned int VertexFaceConnectivity [8][3]
 
static const unsigned int EdgeFaceConnectivity [12][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 49 of file HexGeom.h.

Constructor & Destructor Documentation

◆ HexGeom() [1/2]

Nektar::SpatialDomains::HexGeom::HexGeom ( )

Definition at line 59 of file HexGeom.cpp.

60 {
62 }
LibUtilities::ShapeType m_shapeType
Type of shape.
Definition: Geometry.h:199

References Nektar::LibUtilities::eHexahedron.

◆ HexGeom() [2/2]

Nektar::SpatialDomains::HexGeom::HexGeom ( int  id,
const QuadGeomSharedPtr  faces[] 
)

Copy the face shared pointers

Set up orientation vectors with correct amount of elements.

Definition at line 64 of file HexGeom.cpp.

65  : Geometry3D(faces[0]->GetEdge(0)->GetVertex(0)->GetCoordim())
66 {
68  m_globalID = id;
69 
70  /// Copy the face shared pointers
71  m_faces.insert(m_faces.begin(), faces, faces + HexGeom::kNfaces);
72 
73  /// Set up orientation vectors with correct amount of elements.
74  m_eorient.resize(kNedges);
75  m_forient.resize(kNfaces);
76 
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: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
static const int kNfaces
Definition: HexGeom.h:60
static const int kNedges
Definition: HexGeom.h:57

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

◆ ~HexGeom()

Nektar::SpatialDomains::HexGeom::~HexGeom ( )

Definition at line 83 of file HexGeom.cpp.

84 {
85 }

Member Function Documentation

◆ SetUpEdgeOrientation()

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

Definition at line 650 of file HexGeom.cpp.

651 {
652 
653  // This 2D array holds the local id's of all the vertices
654  // for every edge. For every edge, they are ordered to what we
655  // define as being Forwards
656  const unsigned int edgeVerts[kNedges][2] = {{0, 1},
657  {1, 2},
658  {2, 3},
659  {3, 0},
660  {0, 4},
661  {1, 5},
662  {2, 6},
663  {3, 7},
664  {4, 5},
665  {5, 6},
666  {6, 7},
667  {7, 4}};
668 
669  int i;
670  for (i = 0; i < kNedges; i++)
671  {
672  if (m_edges[i]->GetVid(0) == m_verts[edgeVerts[i][0]]->GetGlobalID())
673  {
675  }
676  else if (m_edges[i]->GetVid(0) == m_verts[edgeVerts[i][1]]->GetGlobalID())
677  {
679  }
680  else
681  {
682  ASSERTL0(false, "Could not find matching vertex for the edge");
683  }
684  }
685 }
#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 HexGeom().

◆ SetUpFaceOrientation()

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

Definition at line 422 of file HexGeom.cpp.

423 {
424  int f, i;
425 
426  // These arrays represent the vector of the A and B
427  // coordinate of the local elemental coordinate system
428  // where A corresponds with the coordinate direction xi_i
429  // with the lowest index i (for that particular face)
430  // Coordinate 'B' then corresponds to the other local
431  // coordinate (i.e. with the highest index)
432  Array<OneD, NekDouble> elementAaxis(m_coordim);
433  Array<OneD, NekDouble> elementBaxis(m_coordim);
434 
435  // These arrays correspond to the local coordinate
436  // system of the face itself (i.e. the Geometry2D)
437  // faceAaxis correspond to the xi_0 axis
438  // faceBaxis correspond to the xi_1 axis
439  Array<OneD, NekDouble> faceAaxis(m_coordim);
440  Array<OneD, NekDouble> faceBaxis(m_coordim);
441 
442  // This is the base vertex of the face (i.e. the Geometry2D)
443  // This corresponds to thevertex with local ID 0 of the
444  // Geometry2D
445  unsigned int baseVertex;
446 
447  // The lenght of the vectors above
448  NekDouble elementAaxis_length;
449  NekDouble elementBaxis_length;
450  NekDouble faceAaxis_length;
451  NekDouble faceBaxis_length;
452 
453  // This 2D array holds the local id's of all the vertices
454  // for every face. For every face, they are ordered in such
455  // a way that the implementation below allows a unified approach
456  // for all faces.
457  const unsigned int faceVerts[kNfaces][QuadGeom::kNverts] = {{0, 1, 2, 3},
458  {0, 1, 5, 4},
459  {1, 2, 6, 5},
460  {3, 2, 6, 7},
461  {0, 3, 7, 4},
462  {4, 5, 6, 7}};
463 
464  NekDouble dotproduct1 = 0.0;
465  NekDouble dotproduct2 = 0.0;
466 
467  unsigned int orientation;
468 
469  // Loop over all the faces to set up the orientation
470  for (f = 0; f < kNqfaces + kNtfaces; f++)
471  {
472  // initialisation
473  elementAaxis_length = 0.0;
474  elementBaxis_length = 0.0;
475  faceAaxis_length = 0.0;
476  faceBaxis_length = 0.0;
477 
478  dotproduct1 = 0.0;
479  dotproduct2 = 0.0;
480 
481  baseVertex = m_faces[f]->GetVid(0);
482 
483  // We are going to construct the vectors representing the A and B axis
484  // of every face. These vectors will be constructed as a
485  // vector-representation
486  // of the edges of the face. However, for both coordinate directions, we
487  // can
488  // represent the vectors by two different edges. That's why we need to
489  // make sure that
490  // we pick the edge to which the baseVertex of the
491  // Geometry2D-representation of the face
492  // belongs...
493  if (baseVertex == m_verts[faceVerts[f][0]]->GetGlobalID())
494  {
495  for (i = 0; i < m_coordim; i++)
496  {
497  elementAaxis[i] = (*m_verts[faceVerts[f][1]])[i] -
498  (*m_verts[faceVerts[f][0]])[i];
499  elementBaxis[i] = (*m_verts[faceVerts[f][3]])[i] -
500  (*m_verts[faceVerts[f][0]])[i];
501  }
502  }
503  else if (baseVertex == m_verts[faceVerts[f][1]]->GetGlobalID())
504  {
505  for (i = 0; i < m_coordim; i++)
506  {
507  elementAaxis[i] = (*m_verts[faceVerts[f][1]])[i] -
508  (*m_verts[faceVerts[f][0]])[i];
509  elementBaxis[i] = (*m_verts[faceVerts[f][2]])[i] -
510  (*m_verts[faceVerts[f][1]])[i];
511  }
512  }
513  else if (baseVertex == m_verts[faceVerts[f][2]]->GetGlobalID())
514  {
515  for (i = 0; i < m_coordim; i++)
516  {
517  elementAaxis[i] = (*m_verts[faceVerts[f][2]])[i] -
518  (*m_verts[faceVerts[f][3]])[i];
519  elementBaxis[i] = (*m_verts[faceVerts[f][2]])[i] -
520  (*m_verts[faceVerts[f][1]])[i];
521  }
522  }
523  else if (baseVertex == m_verts[faceVerts[f][3]]->GetGlobalID())
524  {
525  for (i = 0; i < m_coordim; i++)
526  {
527  elementAaxis[i] = (*m_verts[faceVerts[f][2]])[i] -
528  (*m_verts[faceVerts[f][3]])[i];
529  elementBaxis[i] = (*m_verts[faceVerts[f][3]])[i] -
530  (*m_verts[faceVerts[f][0]])[i];
531  }
532  }
533  else
534  {
535  ASSERTL0(false, "Could not find matching vertex for the face");
536  }
537 
538  // Now, construct the edge-vectors of the local coordinates of
539  // the Geometry2D-representation of the face
540  for (i = 0; i < m_coordim; i++)
541  {
542  faceAaxis[i] =
543  (*m_faces[f]->GetVertex(1))[i] - (*m_faces[f]->GetVertex(0))[i];
544  faceBaxis[i] =
545  (*m_faces[f]->GetVertex(3))[i] - (*m_faces[f]->GetVertex(0))[i];
546 
547  elementAaxis_length += pow(elementAaxis[i], 2);
548  elementBaxis_length += pow(elementBaxis[i], 2);
549  faceAaxis_length += pow(faceAaxis[i], 2);
550  faceBaxis_length += pow(faceBaxis[i], 2);
551  }
552 
553  elementAaxis_length = sqrt(elementAaxis_length);
554  elementBaxis_length = sqrt(elementBaxis_length);
555  faceAaxis_length = sqrt(faceAaxis_length);
556  faceBaxis_length = sqrt(faceBaxis_length);
557 
558  // Calculate the inner product of both the A-axis
559  // (i.e. Elemental A axis and face A axis)
560  for (i = 0; i < m_coordim; i++)
561  {
562  dotproduct1 += elementAaxis[i] * faceAaxis[i];
563  }
564 
565  NekDouble norm = fabs(dotproduct1) / elementAaxis_length / faceAaxis_length;
566  orientation = 0;
567 
568  // if the innerproduct is equal to the (absolute value of the ) products
569  // of the lengths of both vectors, then, the coordinate systems will NOT
570  // be transposed
571  if (fabs(norm - 1.0) < NekConstants::kNekZeroTol)
572  {
573  // if the inner product is negative, both A-axis point
574  // in reverse direction
575  if (dotproduct1 < 0.0)
576  {
577  orientation += 2;
578  }
579 
580  // calculate the inner product of both B-axis
581  for (i = 0; i < m_coordim; i++)
582  {
583  dotproduct2 += elementBaxis[i] * faceBaxis[i];
584  }
585 
586  norm = fabs(dotproduct2) / elementBaxis_length / faceBaxis_length;
587 
588  // check that both these axis are indeed parallel
589  ASSERTL1(fabs(norm - 1.0) < NekConstants::kNekZeroTol,
590  "These vectors should be parallel");
591 
592  // if the inner product is negative, both B-axis point
593  // in reverse direction
594  if (dotproduct2 < 0.0)
595  {
596  orientation++;
597  }
598  }
599  // The coordinate systems are transposed
600  else
601  {
602  orientation = 4;
603 
604  // Calculate the inner product between the elemental A-axis
605  // and the B-axis of the face (which are now the corresponding axis)
606  dotproduct1 = 0.0;
607  for (i = 0; i < m_coordim; i++)
608  {
609  dotproduct1 += elementAaxis[i] * faceBaxis[i];
610  }
611 
612  norm = fabs(dotproduct1) / elementAaxis_length / faceBaxis_length;
613 
614  // check that both these axis are indeed parallel
615  ASSERTL1(fabs(norm - 1.0) < NekConstants::kNekZeroTol,
616  "These vectors should be parallel");
617 
618  // if the result is negative, both axis point in reverse
619  // directions
620  if (dotproduct1 < 0.0)
621  {
622  orientation += 2;
623  }
624 
625  // Do the same for the other two corresponding axis
626  dotproduct2 = 0.0;
627  for (i = 0; i < m_coordim; i++)
628  {
629  dotproduct2 += elementBaxis[i] * faceAaxis[i];
630  }
631 
632  norm = fabs(dotproduct2) / elementBaxis_length / faceAaxis_length;
633 
634  // check that both these axis are indeed parallel
635  ASSERTL1(fabs(norm - 1.0) < NekConstants::kNekZeroTol,
636  "These vectors should be parallel");
637 
638  if (dotproduct2 < 0.0)
639  {
640  orientation++;
641  }
642  }
643 
644  orientation = orientation + 5;
645  // Fill the m_forient array
646  m_forient[f] = (StdRegions::Orientation)orientation;
647  }
648 }
#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 kNqfaces
Definition: HexGeom.h:58
static const int kNtfaces
Definition: HexGeom.h:59
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::SpatialDomains::Geometry::GetGlobalID(), Nektar::SpatialDomains::Geometry::GetVertex(), 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::Geometry3D::m_verts, and tinysimd::sqrt().

Referenced by HexGeom().

◆ SetUpLocalEdges()

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

Definition at line 183 of file HexGeom.cpp.

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

◆ SetUpLocalVertices()

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

Definition at line 331 of file HexGeom.cpp.

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

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

◆ SetUpXmap()

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

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

Definition at line 718 of file HexGeom.cpp.

719 {
720  // Determine necessary order for standard region. This can almost certainly
721  // be simplified but works for now!
722  vector<int> tmp1;
723 
724  if (m_forient[0] < 9)
725  {
726  tmp1.push_back(m_faces[0]->GetXmap()->GetTraceNcoeffs(0));
727  tmp1.push_back(m_faces[0]->GetXmap()->GetTraceNcoeffs(2));
728  }
729  else
730  {
731  tmp1.push_back(m_faces[0]->GetXmap()->GetTraceNcoeffs(1));
732  tmp1.push_back(m_faces[0]->GetXmap()->GetTraceNcoeffs(3));
733  }
734 
735  if (m_forient[5] < 9)
736  {
737  tmp1.push_back(m_faces[5]->GetXmap()->GetTraceNcoeffs(0));
738  tmp1.push_back(m_faces[5]->GetXmap()->GetTraceNcoeffs(2));
739  }
740  else
741  {
742  tmp1.push_back(m_faces[5]->GetXmap()->GetTraceNcoeffs(1));
743  tmp1.push_back(m_faces[5]->GetXmap()->GetTraceNcoeffs(3));
744  }
745 
746  int order0 = *max_element(tmp1.begin(), tmp1.end());
747 
748  tmp1.clear();
749 
750  if (m_forient[0] < 9)
751  {
752  tmp1.push_back(m_faces[0]->GetXmap()->GetTraceNcoeffs(1));
753  tmp1.push_back(m_faces[0]->GetXmap()->GetTraceNcoeffs(3));
754  }
755  else
756  {
757  tmp1.push_back(m_faces[0]->GetXmap()->GetTraceNcoeffs(0));
758  tmp1.push_back(m_faces[0]->GetXmap()->GetTraceNcoeffs(2));
759  }
760 
761  if (m_forient[5] < 9)
762  {
763  tmp1.push_back(m_faces[5]->GetXmap()->GetTraceNcoeffs(1));
764  tmp1.push_back(m_faces[5]->GetXmap()->GetTraceNcoeffs(3));
765  }
766  else
767  {
768  tmp1.push_back(m_faces[5]->GetXmap()->GetTraceNcoeffs(0));
769  tmp1.push_back(m_faces[5]->GetXmap()->GetTraceNcoeffs(2));
770  }
771 
772  int order1 = *max_element(tmp1.begin(), tmp1.end());
773 
774  tmp1.clear();
775 
776  if (m_forient[1] < 9)
777  {
778  tmp1.push_back(m_faces[1]->GetXmap()->GetTraceNcoeffs(1));
779  tmp1.push_back(m_faces[1]->GetXmap()->GetTraceNcoeffs(3));
780  }
781  else
782  {
783  tmp1.push_back(m_faces[1]->GetXmap()->GetTraceNcoeffs(0));
784  tmp1.push_back(m_faces[1]->GetXmap()->GetTraceNcoeffs(2));
785  }
786 
787  if (m_forient[3] < 9)
788  {
789  tmp1.push_back(m_faces[3]->GetXmap()->GetTraceNcoeffs(1));
790  tmp1.push_back(m_faces[3]->GetXmap()->GetTraceNcoeffs(3));
791  }
792  else
793  {
794  tmp1.push_back(m_faces[3]->GetXmap()->GetTraceNcoeffs(0));
795  tmp1.push_back(m_faces[3]->GetXmap()->GetTraceNcoeffs(2));
796  }
797 
798  int order2 = *max_element(tmp1.begin(), tmp1.end());
799 
800  const LibUtilities::BasisKey A(
802  order0,
803  LibUtilities::PointsKey(order0+1, LibUtilities::eGaussLobattoLegendre));
804  const LibUtilities::BasisKey B(
806  order1,
807  LibUtilities::PointsKey(order1+1, LibUtilities::eGaussLobattoLegendre));
808  const LibUtilities::BasisKey C(
810  order2,
811  LibUtilities::PointsKey(order2+1, LibUtilities::eGaussLobattoLegendre));
812 
814 }
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
@ eModified_A
Principle Modified Functions .
Definition: BasisType.h:48

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

Referenced by v_Reset(), and v_Setup().

◆ v_GenGeomFactors()

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

Implements Nektar::SpatialDomains::Geometry.

Definition at line 87 of file HexGeom.cpp.

88 {
89  if(!m_setupState)
90  {
92  }
93 
95  {
96  int i, f;
97  GeomType Gtype = eRegular;
98 
99  v_FillGeom();
100 
101  // check to see if expansions are linear
102  for (i = 0; i < m_coordim; ++i)
103  {
104  if (m_xmap->GetBasisNumModes(0) != 2 ||
105  m_xmap->GetBasisNumModes(1) != 2 ||
106  m_xmap->GetBasisNumModes(2) != 2)
107  {
108  Gtype = eDeformed;
109  }
110  }
111 
112  // check to see if all faces are parallelograms
113  if (Gtype == eRegular)
114  {
115  const unsigned int faceVerts[kNfaces][QuadGeom::kNverts] = {
116  {0, 1, 2, 3},
117  {0, 1, 5, 4},
118  {1, 2, 6, 5},
119  {3, 2, 6, 7},
120  {0, 3, 7, 4},
121  {4, 5, 6, 7}};
122 
123  for (f = 0; f < kNfaces; f++)
124  {
125  // Ensure each face is a parallelogram? Check this.
126  for (i = 0; i < m_coordim; i++)
127  {
128  if (fabs((*m_verts[faceVerts[f][0]])(i) -
129  (*m_verts[faceVerts[f][1]])(i) +
130  (*m_verts[faceVerts[f][2]])(i) -
131  (*m_verts[faceVerts[f][3]])(i)) >
133  {
134  Gtype = eDeformed;
135  break;
136  }
137  }
138 
139  if (Gtype == eDeformed)
140  {
141  break;
142  }
143  }
144  }
145 
147  Gtype, m_coordim, m_xmap, m_coeffs);
149  }
150 }
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, kNfaces, Nektar::SpatialDomains::QuadGeom::kNverts, 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::HexGeom::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 167 of file HexGeom.cpp.

168 {
169  if (faceidx == 0 || faceidx == 5)
170  {
171  return facedir;
172  }
173  else if (faceidx == 1 || faceidx == 3)
174  {
175  return 2 * facedir;
176  }
177  else
178  {
179  return 1 + facedir;
180  }
181 }

◆ v_GetEdgeFaceMap()

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

163 {
164  return EdgeFaceConnectivity[i][j];
165 }
static const unsigned int EdgeFaceConnectivity[12][2]
Definition: HexGeom.h:81

References EdgeFaceConnectivity.

◆ v_GetVertexEdgeMap()

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

153 {
154  return VertexEdgeConnectivity[i][j];
155 }
static const unsigned int VertexEdgeConnectivity[8][3]
Definition: HexGeom.h:79

References VertexEdgeConnectivity.

◆ v_GetVertexFaceMap()

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

158 {
159  return VertexFaceConnectivity[i][j];
160 }
static const unsigned int VertexFaceConnectivity[8][3]
Definition: HexGeom.h:80

References VertexFaceConnectivity.

◆ v_Reset()

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

688 {
689  Geometry::v_Reset(curvedEdges, curvedFaces);
690 
691  for (int i = 0; i < 6; ++i)
692  {
693  m_faces[i]->Reset(curvedEdges, curvedFaces);
694  }
695 
696  SetUpXmap();
697  SetUpCoeffs(m_xmap->GetNcoeffs());
698 }
void SetUpCoeffs(const int nCoeffs)
Initialise the Geometry::m_coeffs array.
Definition: Geometry.h:658
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
void SetUpXmap()
Set up the m_xmap object by determining the order of each direction from derived faces.
Definition: HexGeom.cpp:718

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

◆ v_Setup()

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

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 700 of file HexGeom.cpp.

701 {
702  if(!m_setupState)
703  {
704  for (int i = 0; i < 6; ++i)
705  {
706  m_faces[i]->Setup();
707  }
708  SetUpXmap();
709  SetUpCoeffs(m_xmap->GetNcoeffs());
710  m_setupState = true;
711  }
712 }

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::HexGeom::EdgeFaceConnectivity
staticprivate
Initial value:
= {
{0, 1}, {0, 2}, {0, 3}, {0, 4}, {1, 4}, {1, 2},
{2, 3}, {3, 4}, {1, 5}, {2, 5}, {3, 5}, {4, 5}}

Definition at line 81 of file HexGeom.h.

Referenced by v_GetEdgeFaceMap().

◆ kNedges

const int Nektar::SpatialDomains::HexGeom::kNedges = 12
static

Definition at line 57 of file HexGeom.h.

Referenced by HexGeom(), and SetUpEdgeOrientation().

◆ kNfaces

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

Definition at line 60 of file HexGeom.h.

Referenced by HexGeom(), SetUpFaceOrientation(), and v_GenGeomFactors().

◆ kNqfaces

const int Nektar::SpatialDomains::HexGeom::kNqfaces = 6
static

Definition at line 58 of file HexGeom.h.

Referenced by SetUpFaceOrientation().

◆ kNtfaces

const int Nektar::SpatialDomains::HexGeom::kNtfaces = 0
static

Definition at line 59 of file HexGeom.h.

Referenced by SetUpFaceOrientation().

◆ kNverts

const int Nektar::SpatialDomains::HexGeom::kNverts = 8
static

Definition at line 56 of file HexGeom.h.

◆ VertexEdgeConnectivity

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

Definition at line 79 of file HexGeom.h.

Referenced by v_GetVertexEdgeMap().

◆ VertexFaceConnectivity

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

Definition at line 80 of file HexGeom.h.

Referenced by v_GetVertexFaceMap().

◆ XMLElementType

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

Definition at line 61 of file HexGeom.h.