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...
 
int PreliminaryCheck (const Array< OneD, const NekDouble > &gloCoord)
 A fast and robust check if a given global coord is outside of a deformed element. For regular elements, this check is unnecessary. More...
 
bool MinMaxCheck (const Array< OneD, const NekDouble > &gloCoord)
 Check if given global coord is within the BoundingBox of the element. More...
 
bool ClampLocCoords (Array< OneD, NekDouble > &locCoord, NekDouble tol=std::numeric_limits< NekDouble >::epsilon())
 Clamp local coords to be within standard regions [-1, 1]^dim. More...
 
NekDouble FindDistance (const Array< OneD, const NekDouble > &xs, Array< OneD, NekDouble > &xi)
 
int GetVertexEdgeMap (int i, int j) const
 Returns the standard element edge IDs that are connected to a given vertex. More...
 
int GetVertexFaceMap (int i, int j) const
 Returns the standard element face IDs that are connected to a given vertex. More...
 
int GetEdgeFaceMap (int i, int j) const
 Returns the standard element edge IDs that are connected to a given face. More...
 
int GetEdgeNormalToFaceVert (int i, int j) const
 Returns the standard lement edge IDs that are normal to a given face vertex. More...
 
int GetDir (const int i, const int j=0) const
 Returns the element coordinate direction corresponding to a given face coordinate direction. More...
 
void Reset (CurveMap &curvedEdges, CurveMap &curvedFaces)
 Reset this geometry object: unset the current state, zero Geometry::m_coeffs and remove allocated GeomFactors. More...
 
void Setup ()
 

Static Public Attributes

static const int kNverts = 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)
 
void NewtonIterationForLocCoord (const Array< OneD, const NekDouble > &coords, Array< OneD, NekDouble > &Lcoords)
 
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 void v_CalculateInverseIsoParam () override
 
virtual int v_AllLeftCheck (const Array< OneD, const NekDouble > &gloCoord) override
 
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 PointGeomSharedPtr v_GetVertex (int i) const =0
 
virtual Geometry1DSharedPtr v_GetEdge (int i) const
 Returns edge i of this object. More...
 
virtual Geometry2DSharedPtr v_GetFace (int i) const
 Returns face i of this object. More...
 
virtual StdRegions::Orientation v_GetEorient (const int i) const
 Returns the orientation of edge i with respect to the ordering of edges in the standard element. More...
 
virtual StdRegions::Orientation v_GetForient (const int i) const
 Returns the orientation of face i with respect to the ordering of faces in the standard element. More...
 
virtual int v_GetNumVerts () const
 Get the number of vertices of this object. More...
 
virtual int v_GetNumEdges () const
 Get the number of edges of this object. More...
 
virtual int v_GetNumFaces () const
 Get the number of faces of this object. More...
 
virtual int v_GetShapeDim () const
 Get the object's shape dimension. More...
 
virtual StdRegions::StdExpansionSharedPtr v_GetXmap () const
 Return the mapping object Geometry::m_xmap that represents the coordinate transformation from standard element to physical element. More...
 
virtual void v_FillGeom ()
 Populate the coordinate mapping Geometry::m_coeffs information from any children geometry elements. More...
 
virtual bool v_ContainsPoint (const Array< OneD, const NekDouble > &gloCoord, Array< OneD, NekDouble > &locCoord, NekDouble tol, NekDouble &dist)
 Determine whether an element contains a particular Cartesian coordinate \(\vec{x} = (x,y,z)\). More...
 
virtual int v_AllLeftCheck (const Array< OneD, const NekDouble > &gloCoord)
 
virtual NekDouble v_GetCoord (const int i, const Array< OneD, const NekDouble > &Lcoord)
 Given local collapsed coordinate Lcoord, return the value of physical coordinate in direction i. More...
 
virtual NekDouble v_GetLocCoords (const Array< OneD, const NekDouble > &coords, Array< OneD, NekDouble > &Lcoords)
 Determine the local collapsed coordinates that correspond to a given Cartesian coordinate for this geometry object. More...
 
virtual NekDouble v_FindDistance (const Array< OneD, const NekDouble > &xs, Array< OneD, NekDouble > &xi)
 
virtual int v_GetVertexEdgeMap (int i, int j) const
 Returns the standard element edge IDs that are connected to a given vertex. More...
 
virtual int v_GetVertexFaceMap (int i, int j) const
 Returns the standard element face IDs that are connected to a given vertex. More...
 
virtual int v_GetEdgeFaceMap (int i, int j) const
 Returns the standard element edge IDs that are connected to a given face. More...
 
virtual int v_GetEdgeNormalToFaceVert (const int i, const int j) const
 Returns the standard lement edge IDs that are normal to a given face vertex. More...
 
virtual int v_GetDir (const int faceidx, const int facedir) const
 Returns the element coordinate direction corresponding to a given face coordinate direction. More...
 
virtual void v_Reset (CurveMap &curvedEdges, CurveMap &curvedFaces)
 Reset this geometry object: unset the current state, zero Geometry::m_coeffs and remove allocated GeomFactors. More...
 
virtual void v_Setup ()
 
virtual void v_GenGeomFactors ()=0
 
void SetUpCoeffs (const int nCoeffs)
 Initialise the Geometry::m_coeffs array. More...
 
virtual void v_CalculateInverseIsoParam ()
 

Private Member Functions

void SetUpLocalEdges ()
 
void SetUpLocalVertices ()
 
void SetUpEdgeOrientation ()
 
void SetUpFaceOrientation ()
 
void SetUpXmap ()
 Set up the m_xmap object by determining the order of each direction from derived faces. More...
 

Static Private Attributes

static const unsigned int VertexEdgeConnectivity [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...
 
Array< OneD, Array< OneD, NekDouble > > m_isoParameter
 
Array< OneD, Array< OneD, NekDouble > > m_invIsoParam
 
bool m_straightEdge
 
- 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 ( )

◆ 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:358
Geometry1DSharedPtr GetEdge(int i) const
Returns edge i of this object.
Definition: Geometry.h:366
int GetCoordim() const
Return the coordinate dimension of this object (i.e. the dimension of the space in which this object ...
Definition: Geometry.h:284
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 673 of file HexGeom.cpp.

674{
675
676 // This 2D array holds the local id's of all the vertices
677 // for every edge. For every edge, they are ordered to what we
678 // define as being Forwards
679 const unsigned int edgeVerts[kNedges][2] = {{0, 1}, {1, 2}, {2, 3}, {3, 0},
680 {0, 4}, {1, 5}, {2, 6}, {3, 7},
681 {4, 5}, {5, 6}, {6, 7}, {7, 4}};
682
683 int i;
684 for (i = 0; i < kNedges; i++)
685 {
686 if (m_edges[i]->GetVid(0) == m_verts[edgeVerts[i][0]]->GetGlobalID())
687 {
689 }
690 else if (m_edges[i]->GetVid(0) ==
691 m_verts[edgeVerts[i][1]]->GetGlobalID())
692 {
694 }
695 else
696 {
697 ASSERTL0(false, "Could not find matching vertex for the edge");
698 }
699 }
700}
#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:327

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

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

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

357{
358 // Set up the first 2 vertices (i.e. vertex 0,1)
359 if ((m_edges[0]->GetVid(0) == m_edges[1]->GetVid(0)) ||
360 (m_edges[0]->GetVid(0) == m_edges[1]->GetVid(1)))
361 {
362 m_verts.push_back(m_edges[0]->GetVertex(1));
363 m_verts.push_back(m_edges[0]->GetVertex(0));
364 }
365 else if ((m_edges[0]->GetVid(1) == m_edges[1]->GetVid(0)) ||
366 (m_edges[0]->GetVid(1) == m_edges[1]->GetVid(1)))
367 {
368 m_verts.push_back(m_edges[0]->GetVertex(0));
369 m_verts.push_back(m_edges[0]->GetVertex(1));
370 }
371 else
372 {
373 std::ostringstream errstrm;
374 errstrm << "Connected edges do not share a vertex. Edges ";
375 errstrm << m_edges[0]->GetGlobalID() << ", "
376 << m_edges[1]->GetGlobalID();
377 ASSERTL0(false, errstrm.str());
378 }
379
380 // set up the other bottom vertices (i.e. vertex 2,3)
381 int i;
382 for (i = 1; i < 3; i++)
383 {
384 if (m_edges[i]->GetVid(0) == m_verts[i]->GetGlobalID())
385 {
386 m_verts.push_back(m_edges[i]->GetVertex(1));
387 }
388 else if (m_edges[i]->GetVid(1) == m_verts[i]->GetGlobalID())
389 {
390 m_verts.push_back(m_edges[i]->GetVertex(0));
391 }
392 else
393 {
394 std::ostringstream errstrm;
395 errstrm << "Connected edges do not share a vertex. Edges ";
396 errstrm << m_edges[i]->GetGlobalID() << ", "
397 << m_edges[i - 1]->GetGlobalID();
398 ASSERTL0(false, errstrm.str());
399 }
400 }
401
402 // set up top vertices
403 // First, set up vertices 4,5
404 if ((m_edges[8]->GetVid(0) == m_edges[9]->GetVid(0)) ||
405 (m_edges[8]->GetVid(0) == m_edges[9]->GetVid(1)))
406 {
407 m_verts.push_back(m_edges[8]->GetVertex(1));
408 m_verts.push_back(m_edges[8]->GetVertex(0));
409 }
410 else if ((m_edges[8]->GetVid(1) == m_edges[9]->GetVid(0)) ||
411 (m_edges[8]->GetVid(1) == m_edges[9]->GetVid(1)))
412 {
413 m_verts.push_back(m_edges[8]->GetVertex(0));
414 m_verts.push_back(m_edges[8]->GetVertex(1));
415 }
416 else
417 {
418 std::ostringstream errstrm;
419 errstrm << "Connected edges do not share a vertex. Edges ";
420 errstrm << m_edges[8]->GetGlobalID() << ", "
421 << m_edges[9]->GetGlobalID();
422 ASSERTL0(false, errstrm.str());
423 }
424
425 // set up the other top vertices (i.e. vertex 6,7)
426 for (i = 9; i < 11; i++)
427 {
428 if (m_edges[i]->GetVid(0) == m_verts[i - 4]->GetGlobalID())
429 {
430 m_verts.push_back(m_edges[i]->GetVertex(1));
431 }
432 else if (m_edges[i]->GetVid(1) == m_verts[i - 4]->GetGlobalID())
433 {
434 m_verts.push_back(m_edges[i]->GetVertex(0));
435 }
436 else
437 {
438 std::ostringstream errstrm;
439 errstrm << "Connected edges do not share a vertex. Edges ";
440 errstrm << m_edges[i]->GetGlobalID() << ", "
441 << m_edges[i - 1]->GetGlobalID();
442 ASSERTL0(false, errstrm.str());
443 }
444 }
445}

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

734{
735 // Determine necessary order for standard region. This can almost certainly
736 // be simplified but works for now!
737 vector<int> tmp1;
738
739 if (m_forient[0] < 9)
740 {
741 tmp1.push_back(m_faces[0]->GetXmap()->GetTraceNcoeffs(0));
742 tmp1.push_back(m_faces[0]->GetXmap()->GetTraceNcoeffs(2));
743 }
744 else
745 {
746 tmp1.push_back(m_faces[0]->GetXmap()->GetTraceNcoeffs(1));
747 tmp1.push_back(m_faces[0]->GetXmap()->GetTraceNcoeffs(3));
748 }
749
750 if (m_forient[5] < 9)
751 {
752 tmp1.push_back(m_faces[5]->GetXmap()->GetTraceNcoeffs(0));
753 tmp1.push_back(m_faces[5]->GetXmap()->GetTraceNcoeffs(2));
754 }
755 else
756 {
757 tmp1.push_back(m_faces[5]->GetXmap()->GetTraceNcoeffs(1));
758 tmp1.push_back(m_faces[5]->GetXmap()->GetTraceNcoeffs(3));
759 }
760
761 int order0 = *max_element(tmp1.begin(), tmp1.end());
762
763 tmp1.clear();
764
765 if (m_forient[0] < 9)
766 {
767 tmp1.push_back(m_faces[0]->GetXmap()->GetTraceNcoeffs(1));
768 tmp1.push_back(m_faces[0]->GetXmap()->GetTraceNcoeffs(3));
769 }
770 else
771 {
772 tmp1.push_back(m_faces[0]->GetXmap()->GetTraceNcoeffs(0));
773 tmp1.push_back(m_faces[0]->GetXmap()->GetTraceNcoeffs(2));
774 }
775
776 if (m_forient[5] < 9)
777 {
778 tmp1.push_back(m_faces[5]->GetXmap()->GetTraceNcoeffs(1));
779 tmp1.push_back(m_faces[5]->GetXmap()->GetTraceNcoeffs(3));
780 }
781 else
782 {
783 tmp1.push_back(m_faces[5]->GetXmap()->GetTraceNcoeffs(0));
784 tmp1.push_back(m_faces[5]->GetXmap()->GetTraceNcoeffs(2));
785 }
786
787 int order1 = *max_element(tmp1.begin(), tmp1.end());
788
789 tmp1.clear();
790
791 if (m_forient[1] < 9)
792 {
793 tmp1.push_back(m_faces[1]->GetXmap()->GetTraceNcoeffs(1));
794 tmp1.push_back(m_faces[1]->GetXmap()->GetTraceNcoeffs(3));
795 }
796 else
797 {
798 tmp1.push_back(m_faces[1]->GetXmap()->GetTraceNcoeffs(0));
799 tmp1.push_back(m_faces[1]->GetXmap()->GetTraceNcoeffs(2));
800 }
801
802 if (m_forient[3] < 9)
803 {
804 tmp1.push_back(m_faces[3]->GetXmap()->GetTraceNcoeffs(1));
805 tmp1.push_back(m_faces[3]->GetXmap()->GetTraceNcoeffs(3));
806 }
807 else
808 {
809 tmp1.push_back(m_faces[3]->GetXmap()->GetTraceNcoeffs(0));
810 tmp1.push_back(m_faces[3]->GetXmap()->GetTraceNcoeffs(2));
811 }
812
813 int order2 = *max_element(tmp1.begin(), tmp1.end());
814
815 const LibUtilities::BasisKey A(
817 LibUtilities::PointsKey(order0 + 1,
819 const LibUtilities::BasisKey B(
821 LibUtilities::PointsKey(order1 + 1,
823 const LibUtilities::BasisKey C(
825 LibUtilities::PointsKey(order2 + 1,
827
829}
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:192
StdRegions::StdExpansionSharedPtr GetXmap() const
Return the mapping object Geometry::m_xmap that represents the coordinate transformation from standar...
Definition: Geometry.h:436
@ 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 GeomType Gtype = eRegular;
100
101 v_FillGeom();
102
103 // check to see if expansions are linear
104 m_straightEdge = true;
105 if (m_xmap->GetBasisNumModes(0) != 2 ||
106 m_xmap->GetBasisNumModes(1) != 2 ||
107 m_xmap->GetBasisNumModes(2) != 2)
108 {
109 Gtype = eDeformed;
110 m_straightEdge = false;
111 }
112
113 // check to see if all faces are parallelograms
114 if (Gtype == eRegular)
115 {
116 m_isoParameter = Array<OneD, Array<OneD, NekDouble>>(3);
117 for (int i = 0; i < 3; ++i)
118 {
119 m_isoParameter[i] = Array<OneD, NekDouble>(8, 0.);
120 NekDouble A = (*m_verts[0])(i);
121 NekDouble B = (*m_verts[1])(i);
122 NekDouble C = (*m_verts[2])(i);
123 NekDouble D = (*m_verts[3])(i);
124 NekDouble E = (*m_verts[4])(i);
125 NekDouble F = (*m_verts[5])(i);
126 NekDouble G = (*m_verts[6])(i);
127 NekDouble H = (*m_verts[7])(i);
128 m_isoParameter[i][0] =
129 0.125 * (A + B + C + D + E + F + G + H); // 1
130
131 m_isoParameter[i][1] =
132 0.125 * (-A + B + C - D - E + F + G - H); // xi1
133 m_isoParameter[i][2] =
134 0.125 * (-A - B + C + D - E - F + G + H); // xi2
135 m_isoParameter[i][3] =
136 0.125 * (-A - B - C - D + E + F + G + H); // xi3
137
138 m_isoParameter[i][4] =
139 0.125 * (A - B + C - D + E - F + G - H); // xi1*xi2
140 m_isoParameter[i][5] =
141 0.125 * (A + B - C - D - E - F + G + H); // xi2*xi3
142 m_isoParameter[i][6] =
143 0.125 * (A - B - C + D - E + F + G - H); // xi1*xi3
144
145 m_isoParameter[i][7] =
146 0.125 * (-A + B - C + D + E - F + G - H); // xi1*xi2*xi3
147 NekDouble tmp = fabs(m_isoParameter[i][1]) +
148 fabs(m_isoParameter[i][2]) +
149 fabs(m_isoParameter[i][3]);
151 for (int d = 4; d < 8; ++d)
152 {
153 if (fabs(m_isoParameter[i][d]) > tmp)
154 {
155 Gtype = eDeformed;
156 }
157 }
158 }
159 }
160
161 if (Gtype == eRegular)
162 {
164 }
165
167 Gtype, m_coordim, m_xmap, m_coeffs);
169 }
170}
virtual void v_CalculateInverseIsoParam() override
Definition: Geometry3D.cpp:572
virtual void v_FillGeom() override
Put all quadrature information into face/edge structure and backward transform.
Definition: Geometry3D.cpp:449
bool m_setupState
Wether or not the setup routines have been run.
Definition: Geometry.h:196
Array< OneD, Array< OneD, NekDouble > > m_isoParameter
Definition: Geometry.h:207
Array< OneD, Array< OneD, NekDouble > > m_coeffs
Array containing expansion coefficients of m_xmap.
Definition: Geometry.h:204
GeomState m_geomFactorsState
State of the geometric factors.
Definition: Geometry.h:190
GeomFactorsSharedPtr m_geomFactors
Geometric factors.
Definition: Geometry.h:188
virtual void v_Setup() override
Definition: HexGeom.cpp:715
GeomType
Indicates the type of element geometry.
@ eRegular
Geometry is straight-sided with constant geometric factors.
@ eDeformed
Geometry is curved or has non-constant factors.
@ ePtsFilled
Geometric information has been generated.
std::vector< double > d(NPUPPER *NPUPPER)

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::UnitTests::d(), Nektar::SpatialDomains::eDeformed, Nektar::SpatialDomains::ePtsFilled, Nektar::SpatialDomains::eRegular, Nektar::NekConstants::kNekZeroTol, Nektar::SpatialDomains::Geometry::m_coeffs, Nektar::SpatialDomains::Geometry::m_coordim, Nektar::SpatialDomains::Geometry::m_geomFactors, Nektar::SpatialDomains::Geometry::m_geomFactorsState, Nektar::SpatialDomains::Geometry::m_isoParameter, Nektar::SpatialDomains::Geometry::m_setupState, Nektar::SpatialDomains::Geometry::m_straightEdge, Nektar::SpatialDomains::Geometry3D::m_verts, Nektar::SpatialDomains::Geometry::m_xmap, Nektar::SpatialDomains::Geometry3D::v_CalculateInverseIsoParam(), Nektar::SpatialDomains::Geometry3D::v_FillGeom(), and v_Setup().

◆ v_GetDir()

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

193{
194 if (faceidx == 0 || faceidx == 5)
195 {
196 return facedir;
197 }
198 else if (faceidx == 1 || faceidx == 3)
199 {
200 return 2 * facedir;
201 }
202 else
203 {
204 return 1 + facedir;
205 }
206}

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

183{
184 return EdgeFaceConnectivity[i][j];
185}
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 187 of file HexGeom.cpp.

188{
189 return EdgeNormalToFaceVert[i][j];
190}
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 172 of file HexGeom.cpp.

173{
174 return VertexEdgeConnectivity[i][j];
175}
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 177 of file HexGeom.cpp.

178{
179 return VertexFaceConnectivity[i][j];
180}
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 702 of file HexGeom.cpp.

703{
704 Geometry::v_Reset(curvedEdges, curvedFaces);
705
706 for (int i = 0; i < 6; ++i)
707 {
708 m_faces[i]->Reset(curvedEdges, curvedFaces);
709 }
710
711 SetUpXmap();
712 SetUpCoeffs(m_xmap->GetNcoeffs());
713}
void SetUpCoeffs(const int nCoeffs)
Initialise the Geometry::m_coeffs array.
Definition: Geometry.h:690
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:390
void SetUpXmap()
Set up the m_xmap object by determining the order of each direction from derived faces.
Definition: HexGeom.cpp:733

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

716{
717 if (!m_setupState)
718 {
719 for (int i = 0; i < 6; ++i)
720 {
721 m_faces[i]->Setup();
722 }
723 SetUpXmap();
724 SetUpCoeffs(m_xmap->GetNcoeffs());
725 m_setupState = true;
726 }
727}

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

◆ kNqfaces

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

◆ kNtfaces

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

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