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...
 
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 = 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 () 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 [8][3]
 
static const unsigned int VertexFaceConnectivity [8][3]
 
static const unsigned int EdgeFaceConnectivity [12][2]
 
static const unsigned int EdgeNormalToFaceVert [6][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 49 of file HexGeom.h.

Constructor & Destructor Documentation

◆ HexGeom() [1/2]

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

Definition at line 62 of file HexGeom.cpp.

63 {
65 }
LibUtilities::ShapeType m_shapeType
Type of shape.
Definition: Geometry.h:198

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 67 of file HexGeom.cpp.

68  : Geometry3D(faces[0]->GetEdge(0)->GetVertex(0)->GetCoordim())
69 {
71  m_globalID = id;
72 
73  /// Copy the face shared pointers
74  m_faces.insert(m_faces.begin(), faces, faces + HexGeom::kNfaces);
75 
76  /// Set up orientation vectors with correct amount of elements.
77  m_eorient.resize(kNedges);
78  m_forient.resize(kNfaces);
79 
84 }
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
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 86 of file HexGeom.cpp.

87 {
88 }

Member Function Documentation

◆ SetUpEdgeOrientation()

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

Definition at line 652 of file HexGeom.cpp.

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

◆ SetUpFaceOrientation()

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

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

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

◆ SetUpLocalVertices()

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

Definition at line 335 of file HexGeom.cpp.

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

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 712 of file HexGeom.cpp.

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

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

Implements Nektar::SpatialDomains::Geometry.

Definition at line 90 of file HexGeom.cpp.

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

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

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 171 of file HexGeom.cpp.

172 {
173  if (faceidx == 0 || faceidx == 5)
174  {
175  return facedir;
176  }
177  else if (faceidx == 1 || faceidx == 3)
178  {
179  return 2 * facedir;
180  }
181  else
182  {
183  return 1 + facedir;
184  }
185 }

◆ v_GetEdgeFaceMap()

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

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

References EdgeFaceConnectivity.

◆ v_GetEdgeNormalToFaceVert()

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

167 {
168  return EdgeNormalToFaceVert[i][j];
169 }
static const unsigned int EdgeNormalToFaceVert[6][4]
Definition: HexGeom.h:84

References EdgeNormalToFaceVert.

◆ v_GetVertexEdgeMap()

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

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

References VertexEdgeConnectivity.

◆ v_GetVertexFaceMap()

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

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

References VertexFaceConnectivity.

◆ v_Reset()

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

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

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

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 694 of file HexGeom.cpp.

695 {
696  if (!m_setupState)
697  {
698  for (int i = 0; i < 6; ++i)
699  {
700  m_faces[i]->Setup();
701  }
702  SetUpXmap();
703  SetUpCoeffs(m_xmap->GetNcoeffs());
704  m_setupState = true;
705  }
706 }

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 83 of file HexGeom.h.

Referenced by v_GetEdgeFaceMap().

◆ EdgeNormalToFaceVert

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

Definition at line 84 of file HexGeom.h.

Referenced by v_GetEdgeNormalToFaceVert().

◆ 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 81 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 82 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.