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 () override=default
 
- Public Member Functions inherited from Nektar::SpatialDomains::Geometry3D
 Geometry3D ()
 
 Geometry3D (const int coordim)
 
 ~Geometry3D () override=default
 
- 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
 
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 ResetNonRecursive (CurveMap &curvedEdges, CurveMap &curvedFaces)
 Reset this geometry object non-recursively: unset the current state, zero Geometry::m_coeffs and remove allocated GeomFactors. More...
 
void Setup ()
 
void GenGeomFactors ()
 Handles generation of geometry factors. More...
 

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

void v_GenGeomFactors () override
 
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...
 
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...
 
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...
 
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...
 
int v_GetDir (const int faceidx, const int facedir) const override
 Returns the element coordinate direction corresponding to a given face coordinate direction. More...
 
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...
 
void v_Setup () override
 
- Protected Member Functions inherited from Nektar::SpatialDomains::Geometry3D
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)
 
void v_FillGeom () override
 Put all quadrature information into face/edge structure and backward transform. More...
 
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...
 
void v_CalculateInverseIsoParam () override
 
int v_AllLeftCheck (const Array< OneD, const NekDouble > &gloCoord) override
 
int v_GetShapeDim () const override
 Get the object's shape dimension. More...
 
int v_GetNumVerts () const override
 Get the number of vertices of this object. More...
 
int v_GetNumEdges () const override
 Get the number of edges of this object. More...
 
int v_GetNumFaces () const override
 Get the number of faces of this object. More...
 
PointGeomSharedPtr v_GetVertex (int i) const override
 
Geometry1DSharedPtr v_GetEdge (int i) const override
 Returns edge i of this object. More...
 
Geometry2DSharedPtr v_GetFace (int i) const override
 Returns face i of this object. More...
 
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...
 
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
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
 
int m_straightEdge
 
- Static Protected Attributes inherited from Nektar::SpatialDomains::Geometry
static GeomFactorsVector m_regGeomFactorsManager
 

Detailed Description

Definition at line 47 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 63 of file HexGeom.cpp.

64 : Geometry3D(faces[0]->GetEdge(0)->GetVertex(0)->GetCoordim())
65{
67 m_globalID = id;
68
69 /// Copy the face shared pointers
70 m_faces.insert(m_faces.begin(), faces, faces + HexGeom::kNfaces);
71
72 /// Set up orientation vectors with correct amount of elements.
73 m_eorient.resize(kNedges);
74 m_forient.resize(kNfaces);
75
80}
std::vector< StdRegions::Orientation > m_forient
Definition: Geometry3D.h:81
std::vector< StdRegions::Orientation > m_eorient
Definition: Geometry3D.h:80
PointGeomSharedPtr GetVertex(int i) const
Returns vertex i of this object.
Definition: Geometry.h:357
Geometry1DSharedPtr GetEdge(int i) const
Returns edge i of this object.
Definition: Geometry.h:365
int GetCoordim() const
Return the coordinate dimension of this object (i.e. the dimension of the space in which this object ...
Definition: Geometry.h:283
static const int kNfaces
Definition: HexGeom.h:58
static const int kNedges
Definition: HexGeom.h:55

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

Member Function Documentation

◆ SetUpEdgeOrientation()

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

Definition at line 665 of file HexGeom.cpp.

666{
667
668 // This 2D array holds the local id's of all the vertices
669 // for every edge. For every edge, they are ordered to what we
670 // define as being Forwards
671 const unsigned int edgeVerts[kNedges][2] = {{0, 1}, {1, 2}, {2, 3}, {3, 0},
672 {0, 4}, {1, 5}, {2, 6}, {3, 7},
673 {4, 5}, {5, 6}, {6, 7}, {7, 4}};
674
675 int i;
676 for (i = 0; i < kNedges; i++)
677 {
678 if (m_edges[i]->GetVid(0) == m_verts[edgeVerts[i][0]]->GetGlobalID())
679 {
681 }
682 else if (m_edges[i]->GetVid(0) ==
683 m_verts[edgeVerts[i][1]]->GetGlobalID())
684 {
686 }
687 else
688 {
689 ASSERTL0(false, "Could not find matching vertex for the edge");
690 }
691 }
692}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
int GetVid(int i) const
Get the ID of vertex i of this object.
Definition: Geometry.cpp:104
int GetGlobalID(void) const
Get the ID of this object.
Definition: Geometry.h:326

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

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

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

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

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

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

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

726{
727 // Determine necessary order for standard region. This can almost certainly
728 // be simplified but works for now!
729 std::vector<int> tmp1;
730
731 if (m_forient[0] < 9)
732 {
733 tmp1.push_back(m_faces[0]->GetXmap()->GetTraceNcoeffs(0));
734 tmp1.push_back(m_faces[0]->GetXmap()->GetTraceNcoeffs(2));
735 }
736 else
737 {
738 tmp1.push_back(m_faces[0]->GetXmap()->GetTraceNcoeffs(1));
739 tmp1.push_back(m_faces[0]->GetXmap()->GetTraceNcoeffs(3));
740 }
741
742 if (m_forient[5] < 9)
743 {
744 tmp1.push_back(m_faces[5]->GetXmap()->GetTraceNcoeffs(0));
745 tmp1.push_back(m_faces[5]->GetXmap()->GetTraceNcoeffs(2));
746 }
747 else
748 {
749 tmp1.push_back(m_faces[5]->GetXmap()->GetTraceNcoeffs(1));
750 tmp1.push_back(m_faces[5]->GetXmap()->GetTraceNcoeffs(3));
751 }
752
753 int order0 = *std::max_element(tmp1.begin(), tmp1.end());
754
755 tmp1.clear();
756
757 if (m_forient[0] < 9)
758 {
759 tmp1.push_back(m_faces[0]->GetXmap()->GetTraceNcoeffs(1));
760 tmp1.push_back(m_faces[0]->GetXmap()->GetTraceNcoeffs(3));
761 }
762 else
763 {
764 tmp1.push_back(m_faces[0]->GetXmap()->GetTraceNcoeffs(0));
765 tmp1.push_back(m_faces[0]->GetXmap()->GetTraceNcoeffs(2));
766 }
767
768 if (m_forient[5] < 9)
769 {
770 tmp1.push_back(m_faces[5]->GetXmap()->GetTraceNcoeffs(1));
771 tmp1.push_back(m_faces[5]->GetXmap()->GetTraceNcoeffs(3));
772 }
773 else
774 {
775 tmp1.push_back(m_faces[5]->GetXmap()->GetTraceNcoeffs(0));
776 tmp1.push_back(m_faces[5]->GetXmap()->GetTraceNcoeffs(2));
777 }
778
779 int order1 = *std::max_element(tmp1.begin(), tmp1.end());
780
781 tmp1.clear();
782
783 if (m_forient[1] < 9)
784 {
785 tmp1.push_back(m_faces[1]->GetXmap()->GetTraceNcoeffs(1));
786 tmp1.push_back(m_faces[1]->GetXmap()->GetTraceNcoeffs(3));
787 }
788 else
789 {
790 tmp1.push_back(m_faces[1]->GetXmap()->GetTraceNcoeffs(0));
791 tmp1.push_back(m_faces[1]->GetXmap()->GetTraceNcoeffs(2));
792 }
793
794 if (m_forient[3] < 9)
795 {
796 tmp1.push_back(m_faces[3]->GetXmap()->GetTraceNcoeffs(1));
797 tmp1.push_back(m_faces[3]->GetXmap()->GetTraceNcoeffs(3));
798 }
799 else
800 {
801 tmp1.push_back(m_faces[3]->GetXmap()->GetTraceNcoeffs(0));
802 tmp1.push_back(m_faces[3]->GetXmap()->GetTraceNcoeffs(2));
803 }
804
805 int order2 = *std::max_element(tmp1.begin(), tmp1.end());
806
807 const LibUtilities::BasisKey A(
809 LibUtilities::PointsKey(order0 + 1,
811 const LibUtilities::BasisKey B(
813 LibUtilities::PointsKey(order1 + 1,
815 const LibUtilities::BasisKey C(
817 LibUtilities::PointsKey(order2 + 1,
819
821}
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:194
StdRegions::StdExpansionSharedPtr GetXmap() const
Return the mapping object Geometry::m_xmap that represents the coordinate transformation from standar...
Definition: Geometry.h:435
@ eGaussLobattoLegendre
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:51
@ eModified_A
Principle Modified Functions .
Definition: BasisType.h:48

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

Referenced by v_Reset(), and v_Setup().

◆ v_GenGeomFactors()

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

Implements Nektar::SpatialDomains::Geometry.

Definition at line 82 of file HexGeom.cpp.

83{
84 if (!m_setupState)
85 {
87 }
88
90 {
91 GeomType Gtype = eRegular;
92
93 v_FillGeom();
94
95 // check to see if expansions are linear
97 if (m_xmap->GetBasisNumModes(0) != 2 ||
98 m_xmap->GetBasisNumModes(1) != 2 ||
99 m_xmap->GetBasisNumModes(2) != 2)
100 {
101 Gtype = eDeformed;
102 m_straightEdge = 0;
103 }
104
105 // check to see if all faces are parallelograms
106 if (Gtype == eRegular)
107 {
108 m_isoParameter = Array<OneD, Array<OneD, NekDouble>>(3);
109 for (int i = 0; i < 3; ++i)
110 {
111 m_isoParameter[i] = Array<OneD, NekDouble>(8, 0.);
112 NekDouble A = (*m_verts[0])(i);
113 NekDouble B = (*m_verts[1])(i);
114 NekDouble C = (*m_verts[2])(i);
115 NekDouble D = (*m_verts[3])(i);
116 NekDouble E = (*m_verts[4])(i);
117 NekDouble F = (*m_verts[5])(i);
118 NekDouble G = (*m_verts[6])(i);
119 NekDouble H = (*m_verts[7])(i);
120 m_isoParameter[i][0] =
121 0.125 * (A + B + C + D + E + F + G + H); // 1
122
123 m_isoParameter[i][1] =
124 0.125 * (-A + B + C - D - E + F + G - H); // xi1
125 m_isoParameter[i][2] =
126 0.125 * (-A - B + C + D - E - F + G + H); // xi2
127 m_isoParameter[i][3] =
128 0.125 * (-A - B - C - D + E + F + G + H); // xi3
129
130 m_isoParameter[i][4] =
131 0.125 * (A - B + C - D + E - F + G - H); // xi1*xi2
132 m_isoParameter[i][5] =
133 0.125 * (A + B - C - D - E - F + G + H); // xi2*xi3
134 m_isoParameter[i][6] =
135 0.125 * (A - B - C + D - E + F + G - H); // xi1*xi3
136
137 m_isoParameter[i][7] =
138 0.125 * (-A + B - C + D + E - F + G - H); // xi1*xi2*xi3
139 NekDouble tmp = fabs(m_isoParameter[i][1]) +
140 fabs(m_isoParameter[i][2]) +
141 fabs(m_isoParameter[i][3]);
143 for (int d = 4; d < 8; ++d)
144 {
145 if (fabs(m_isoParameter[i][d]) > tmp)
146 {
147 Gtype = eDeformed;
148 }
149 }
150 }
151 }
152
153 if (Gtype == eRegular)
154 {
156 }
157
159 Gtype, m_coordim, m_xmap, m_coeffs);
161 }
162}
void v_CalculateInverseIsoParam() override
Definition: Geometry3D.cpp:565
void v_FillGeom() override
Put all quadrature information into face/edge structure and backward transform.
Definition: Geometry3D.cpp:442
bool m_setupState
Wether or not the setup routines have been run.
Definition: Geometry.h:198
Array< OneD, Array< OneD, NekDouble > > m_isoParameter
Definition: Geometry.h:209
Array< OneD, Array< OneD, NekDouble > > m_coeffs
Array containing expansion coefficients of m_xmap.
Definition: Geometry.h:206
GeomState m_geomFactorsState
State of the geometric factors.
Definition: Geometry.h:192
GeomFactorsSharedPtr m_geomFactors
Geometric factors.
Definition: Geometry.h:190
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 184 of file HexGeom.cpp.

185{
186 if (faceidx == 0 || faceidx == 5)
187 {
188 return facedir;
189 }
190 else if (faceidx == 1 || faceidx == 3)
191 {
192 return 2 * facedir;
193 }
194 else
195 {
196 return 1 + facedir;
197 }
198}

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

175{
176 return EdgeFaceConnectivity[i][j];
177}
static const unsigned int EdgeFaceConnectivity[12][2]
Definition: HexGeom.h:80

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

180{
181 return EdgeNormalToFaceVert[i][j];
182}
static const unsigned int EdgeNormalToFaceVert[6][4]
Definition: HexGeom.h:81

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

165{
166 return VertexEdgeConnectivity[i][j];
167}
static const unsigned int VertexEdgeConnectivity[8][3]
Definition: HexGeom.h:78

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

170{
171 return VertexFaceConnectivity[i][j];
172}
static const unsigned int VertexFaceConnectivity[8][3]
Definition: HexGeom.h:79

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

695{
696 Geometry::v_Reset(curvedEdges, curvedFaces);
697
698 for (int i = 0; i < 6; ++i)
699 {
700 m_faces[i]->Reset(curvedEdges, curvedFaces);
701 }
702
703 SetUpXmap();
704 SetUpCoeffs(m_xmap->GetNcoeffs());
705}
void SetUpCoeffs(const int nCoeffs)
Initialise the Geometry::m_coeffs array.
Definition: Geometry.h:700
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:364
void SetUpXmap()
Set up the m_xmap object by determining the order of each direction from derived faces.
Definition: HexGeom.cpp:725

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

708{
709 if (!m_setupState)
710 {
711 for (int i = 0; i < 6; ++i)
712 {
713 m_faces[i]->Setup();
714 }
715 SetUpXmap();
716 SetUpCoeffs(m_xmap->GetNcoeffs());
717 m_setupState = true;
718 }
719}

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

Referenced by v_GetEdgeNormalToFaceVert().

◆ kNedges

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

Definition at line 55 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 54 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 78 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 79 of file HexGeom.h.

Referenced by v_GetVertexFaceMap().

◆ XMLElementType

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

Definition at line 59 of file HexGeom.h.