Nektar++
Loading...
Searching...
No Matches
Public Member Functions | Static Public Attributes | Protected Member Functions | Protected Attributes | Private Member Functions | Static Private Attributes | List of all members
Nektar::SpatialDomains::TetGeom Class Reference

#include <TetGeom.h>

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

Public Member Functions

 TetGeom ()
 
 TetGeom (int id, std::array< TriGeom *, kNfaces > faces)
 
 ~TetGeom () 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.
 
 Geometry (int coordim)
 Constructor when supplied a coordinate dimension.
 
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).
 
void SetCoordim (int coordim)
 Sets the coordinate dimension of this object (i.e. the dimension of the space in which this object is embedded).
 
GeomFactorsSharedPtr GetGeomFactors ()
 Get the geometric factors for this object, generating them if required.
 
GeomFactorsSharedPtr GetRefGeomFactors (const Array< OneD, const LibUtilities::BasisSharedPtr > &tbasis)
 
GeomFactorsSharedPtr GetMetricInfo ()
 Get the geometric factors for this object.
 
LibUtilities::ShapeType GetShapeType (void)
 Get the geometric shape type of this object.
 
int GetGlobalID (void) const
 Get the ID of this object.
 
void SetGlobalID (int globalid)
 Set the ID of this object.
 
int GetVid (int i) const
 Returns global id of vertex i of this object.
 
int GetEid (int i) const
 Get the ID of edge i of this object.
 
int GetFid (int i) const
 Get the ID of face i of this object.
 
int GetTid (int i) const
 Get the ID of trace i of this object.
 
PointGeomGetVertex (int i) const
 Returns vertex i of this object.
 
Geometry1DGetEdge (int i) const
 Returns edge i of this object.
 
Geometry2DGetFace (int i) const
 Returns face i of this object.
 
StdRegions::Orientation GetEorient (const int i) const
 Returns the orientation of edge i with respect to the ordering of edges in the standard element.
 
StdRegions::Orientation GetForient (const int i) const
 Returns the orientation of face i with respect to the ordering of faces in the standard element.
 
int GetNumVerts () const
 Get the number of vertices of this object.
 
int GetNumEdges () const
 Get the number of edges of this object.
 
int GetNumFaces () const
 Get the number of faces of this object.
 
int GetShapeDim () const
 Get the object's shape dimension.
 
StdRegions::StdExpansionSharedPtr GetXmap () const
 Return the mapping object Geometry::m_xmap that represents the coordinate transformation from standard element to physical element.
 
const Array< OneD, const NekDouble > & GetCoeffs (const int i) const
 Return the coefficients of the transformation Geometry::m_xmap in coordinate direction i.
 
void FillGeom ()
 Populate the coordinate mapping Geometry::m_coeffs information from any children geometry elements.
 
std::array< NekDouble, 6 > GetBoundingBox ()
 Generates the bounding box for the element.
 
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)\).
 
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)\).
 
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)\).
 
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.
 
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.
 
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.
 
bool MinMaxCheck (const Array< OneD, const NekDouble > &gloCoord)
 Check if given global coord is within the BoundingBox of the element.
 
bool ClampLocCoords (Array< OneD, NekDouble > &locCoord, NekDouble tol=std::numeric_limits< NekDouble >::epsilon())
 Clamp local coords to be within standard regions [-1, 1]^dim.
 
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.
 
int GetVertexFaceMap (int i, int j) const
 Returns the standard element face IDs that are connected to a given vertex.
 
int GetEdgeFaceMap (int i, int j) const
 Returns the standard element edge IDs that are connected to a given face.
 
int GetEdgeNormalToFaceVert (int i, int j) const
 Returns the standard lement edge IDs that are normal to a given face vertex.
 
int GetDir (const int i, const int j=0) const
 Returns the element coordinate direction corresponding to a given face coordinate direction.
 
void Reset (CurveMap &curvedEdges, CurveMap &curvedFaces)
 Reset this geometry object: unset the current state, zero Geometry::m_coeffs and remove allocated GeomFactors.
 
void ResetNonRecursive (CurveMap &curvedEdges, CurveMap &curvedFaces)
 Reset this geometry object non-recursively: unset the current state, zero Geometry::m_coeffs and remove allocated GeomFactors.
 
void Setup ()
 
void GenGeomFactors ()
 Handles generation of geometry factors.
 

Static Public Attributes

static const int kNverts = 4
 
static const int kNedges = 6
 
static const int kNqfaces = 0
 
static const int kNtfaces = 4
 
static const int kNfaces = kNqfaces + kNtfaces
 
static const int kNfacets = kNfaces
 
static const std::string XMLElementType
 
- Static Public Attributes inherited from Nektar::SpatialDomains::Geometry3D
static const int kDim = 3
 

Protected Member Functions

int v_GetVertexEdgeMap (const int i, const int j) const override
 Returns the standard element edge IDs that are connected to a given vertex.
 
int v_GetVertexFaceMap (const int i, const int j) const override
 Returns the standard element face IDs that are connected to a given vertex.
 
int v_GetEdgeFaceMap (const int i, const int j) const override
 Returns the standard element edge IDs that are connected to a given face.
 
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.
 
int v_GetDir (const int faceidx, const int facedir) const override
 Returns the element coordinate direction corresponding to a given face coordinate direction.
 
void v_Reset (CurveMap &curvedEdges, CurveMap &curvedFaces) override
 Reset this geometry object: unset the current state, zero Geometry::m_coeffs and remove allocated GeomFactors.
 
void v_Setup () override
 
void v_GenGeomFactors () override
 
void v_FillGeom () override
 Put all quadrature information into face/edge structure and backward transform.
 
int v_GetNumVerts () const final
 Get the number of vertices of this object.
 
int v_GetNumEdges () const final
 Get the number of edges of this object.
 
int v_GetNumFaces () const final
 Get the number of faces of this object.
 
PointGeomv_GetVertex (const int i) const final
 Returns vertex i of this object.
 
Geometry1Dv_GetEdge (const int i) const final
 Returns edge i of this object.
 
Geometry2Dv_GetFace (const int i) const final
 Returns face i of this object.
 
StdRegions::Orientation v_GetEorient (const int i) const final
 Returns the orientation of edge i with respect to the ordering of edges in the standard element.
 
StdRegions::Orientation v_GetForient (const int i) const final
 Returns the orientation of face i with respect to the ordering of faces in the standard element.
 
- 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.
 
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)
 
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.
 
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.
 
- Protected Member Functions inherited from Nektar::SpatialDomains::Geometry
virtual int v_GetVid (int i) const
 Get the ID of vertex i of this object.
 
virtual StdRegions::StdExpansionSharedPtr v_GetXmap () const
 Return the mapping object Geometry::m_xmap that represents the coordinate transformation from standard element to physical element.
 
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)\).
 
virtual NekDouble v_FindDistance (const Array< OneD, const NekDouble > &xs, Array< OneD, NekDouble > &xi)
 
void SetUpCoeffs (const int nCoeffs)
 Initialise the Geometry::m_coeffs array.
 

Protected Attributes

std::array< PointGeom *, kNvertsm_verts
 
std::array< SegGeom *, kNedgesm_edges
 
std::array< TriGeom *, kNfacesm_faces
 
std::array< StdRegions::Orientation, kNedgesm_eorient
 
std::array< StdRegions::Orientation, kNfacesm_forient
 
- Protected Attributes inherited from Nektar::SpatialDomains::Geometry3D
int m_eid
 
bool m_ownverts
 
- Protected Attributes inherited from Nektar::SpatialDomains::Geometry
int m_coordim
 Coordinate dimension of this geometry object.
 
GeomFactorsSharedPtr m_geomFactors
 Geometric factors.
 
GeomState m_geomFactorsState
 State of the geometric factors.
 
StdRegions::StdExpansionSharedPtr m_xmap
 \(\chi\) mapping containing isoparametric transformation.
 
GeomState m_state
 Enumeration to dictate whether coefficients are filled.
 
bool m_setupState
 Wether or not the setup routines have been run.
 
GeomType m_geomType
 Type of geometry.
 
LibUtilities::ShapeType m_shapeType
 Type of shape.
 
int m_globalID
 Global ID.
 
Array< OneD, Array< OneD, NekDouble > > m_coeffs
 Array containing expansion coefficients of m_xmap.
 
Array< OneD, NekDoublem_boundingBox
 Array containing bounding box.
 
Array< OneD, Array< OneD, NekDouble > > m_isoParameter
 
Array< OneD, Array< OneD, NekDouble > > m_invIsoParam
 
int m_straightEdge
 

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.
 

Static Private Attributes

static const unsigned int VertexEdgeConnectivity [4][3]
 
static const unsigned int VertexFaceConnectivity [4][3]
 
static const unsigned int EdgeFaceConnectivity [6][2]
 
static const unsigned int EdgeNormalToFaceVert [4][3]
 

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.
 
- Static Protected Attributes inherited from Nektar::SpatialDomains::Geometry
static GeomFactorsVector m_regGeomFactorsManager
 

Detailed Description

Definition at line 45 of file TetGeom.h.

Constructor & Destructor Documentation

◆ TetGeom() [1/2]

Nektar::SpatialDomains::TetGeom::TetGeom ( )

◆ TetGeom() [2/2]

Nektar::SpatialDomains::TetGeom::TetGeom ( int  id,
std::array< TriGeom *, kNfaces faces 
)

Definition at line 64 of file TetGeom.cpp.

65 : Geometry3D(faces[0]->GetEdge(0)->GetVertex(0)->GetCoordim())
66{
68 m_globalID = id;
69 m_faces = faces;
70
75}
PointGeom * GetVertex(int i) const
Returns vertex i of this object.
Definition Geometry.h:361
int GetCoordim() const
Return the coordinate dimension of this object (i.e. the dimension of the space in which this object ...
Definition Geometry.h:279
Geometry1D * GetEdge(int i) const
Returns edge i of this object.
Definition Geometry.h:369
std::array< TriGeom *, kNfaces > m_faces
Definition TetGeom.h:114

References Nektar::LibUtilities::eTetrahedron, m_faces, Nektar::SpatialDomains::Geometry::m_globalID, Nektar::SpatialDomains::Geometry::m_shapeType, SetUpEdgeOrientation(), SetUpFaceOrientation(), SetUpLocalEdges(), and SetUpLocalVertices().

◆ ~TetGeom()

Nektar::SpatialDomains::TetGeom::~TetGeom ( )
overridedefault

Member Function Documentation

◆ SetUpEdgeOrientation()

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

Definition at line 310 of file TetGeom.cpp.

311{
312
313 // This 2D array holds the local id's of all the vertices
314 // for every edge. For every edge, they are ordered to what we
315 // define as being Forwards
316 const unsigned int edgeVerts[kNedges][2] = {{0, 1}, {1, 2}, {0, 2},
317 {0, 3}, {1, 3}, {2, 3}};
318
319 int i;
320 for (i = 0; i < kNedges; i++)
321 {
322 if (m_edges[i]->GetVid(0) == m_verts[edgeVerts[i][0]]->GetGlobalID())
323 {
325 }
326 else if (m_edges[i]->GetVid(0) ==
327 m_verts[edgeVerts[i][1]]->GetGlobalID())
328 {
330 }
331 else
332 {
333 ASSERTL0(false, "Could not find matching vertex for the edge");
334 }
335 }
336}
#define ASSERTL0(condition, msg)
int GetVid(int i) const
Returns global id of vertex i of this object.
Definition Geometry.h:353
int GetGlobalID(void) const
Get the ID of this object.
Definition Geometry.h:322
std::array< SegGeom *, kNedges > m_edges
Definition TetGeom.h:113
static const int kNedges
Definition TetGeom.h:49
std::array< PointGeom *, kNverts > m_verts
Definition TetGeom.h:112
std::array< StdRegions::Orientation, kNedges > m_eorient
Definition TetGeom.h:115

References ASSERTL0, Nektar::StdRegions::eBackwards, Nektar::StdRegions::eForwards, Nektar::SpatialDomains::Geometry::GetGlobalID(), Nektar::SpatialDomains::Geometry::GetVid(), kNedges, m_edges, m_eorient, and m_verts.

Referenced by TetGeom().

◆ SetUpFaceOrientation()

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

Definition at line 338 of file TetGeom.cpp.

339{
340
341 int f, i;
342
343 // These arrays represent the vector of the A and B
344 // coordinate of the local elemental coordinate system
345 // where A corresponds with the coordinate direction xi_i
346 // with the lowest index i (for that particular face)
347 // Coordinate 'B' then corresponds to the other local
348 // coordinate (i.e. with the highest index)
349 Array<OneD, NekDouble> elementAaxis(m_coordim);
350 Array<OneD, NekDouble> elementBaxis(m_coordim);
351
352 // These arrays correspond to the local coordinate
353 // system of the face itself (i.e. the Geometry2D)
354 // faceAaxis correspond to the xi_0 axis
355 // faceBaxis correspond to the xi_1 axis
356 Array<OneD, NekDouble> faceAaxis(m_coordim);
357 Array<OneD, NekDouble> faceBaxis(m_coordim);
358
359 // This is the base vertex of the face (i.e. the Geometry2D)
360 // This corresponds to thevertex with local ID 0 of the
361 // Geometry2D
362 unsigned int baseVertex;
363
364 // The lenght of the vectors above
365 NekDouble elementAaxis_length;
366 NekDouble elementBaxis_length;
367 NekDouble faceAaxis_length;
368 NekDouble faceBaxis_length;
369
370 // This 2D array holds the local id's of all the vertices
371 // for every face. For every face, they are ordered in such
372 // a way that the implementation below allows a unified approach
373 // for all faces.
374 const unsigned int faceVerts[kNfaces][TriGeom::kNverts] = {
375 {0, 1, 2}, {0, 1, 3}, {1, 2, 3}, {0, 2, 3}};
376
377 NekDouble dotproduct1 = 0.0;
378 NekDouble dotproduct2 = 0.0;
379
380 unsigned int orientation;
381
382 // Loop over all the faces to set up the orientation
383 for (f = 0; f < kNqfaces + kNtfaces; f++)
384 {
385 // initialisation
386 elementAaxis_length = 0.0;
387 elementBaxis_length = 0.0;
388 faceAaxis_length = 0.0;
389 faceBaxis_length = 0.0;
390
391 dotproduct1 = 0.0;
392 dotproduct2 = 0.0;
393
394 baseVertex = m_faces[f]->GetVid(0);
395
396 // We are going to construct the vectors representing the A and B axis
397 // of every face. These vectors will be constructed as a
398 // vector-representation
399 // of the edges of the face. However, for both coordinate directions, we
400 // can
401 // represent the vectors by two different edges. That's why we need to
402 // make sure that
403 // we pick the edge to which the baseVertex of the
404 // Geometry2D-representation of the face
405 // belongs...
406
407 // Compute the length of edges on a base-face
408 if (baseVertex == m_verts[faceVerts[f][0]]->GetGlobalID())
409 {
410 for (i = 0; i < m_coordim; i++)
411 {
412 elementAaxis[i] = (*m_verts[faceVerts[f][1]])[i] -
413 (*m_verts[faceVerts[f][0]])[i];
414 elementBaxis[i] = (*m_verts[faceVerts[f][2]])[i] -
415 (*m_verts[faceVerts[f][0]])[i];
416 }
417 }
418 else if (baseVertex == m_verts[faceVerts[f][1]]->GetGlobalID())
419 {
420 for (i = 0; i < m_coordim; i++)
421 {
422 elementAaxis[i] = (*m_verts[faceVerts[f][1]])[i] -
423 (*m_verts[faceVerts[f][0]])[i];
424 elementBaxis[i] = (*m_verts[faceVerts[f][2]])[i] -
425 (*m_verts[faceVerts[f][1]])[i];
426 }
427 }
428 else if (baseVertex == m_verts[faceVerts[f][2]]->GetGlobalID())
429 {
430 for (i = 0; i < m_coordim; i++)
431 {
432 elementAaxis[i] = (*m_verts[faceVerts[f][1]])[i] -
433 (*m_verts[faceVerts[f][2]])[i];
434 elementBaxis[i] = (*m_verts[faceVerts[f][2]])[i] -
435 (*m_verts[faceVerts[f][0]])[i];
436 }
437 }
438 else
439 {
440 ASSERTL0(false, "Could not find matching vertex for the face");
441 }
442
443 // Now, construct the edge-vectors of the local coordinates of
444 // the Geometry2D-representation of the face
445 for (i = 0; i < m_coordim; i++)
446 {
447 faceAaxis[i] =
448 (*m_faces[f]->GetVertex(1))[i] - (*m_faces[f]->GetVertex(0))[i];
449 faceBaxis[i] =
450 (*m_faces[f]->GetVertex(2))[i] - (*m_faces[f]->GetVertex(0))[i];
451
452 elementAaxis_length += pow(elementAaxis[i], 2);
453 elementBaxis_length += pow(elementBaxis[i], 2);
454 faceAaxis_length += pow(faceAaxis[i], 2);
455 faceBaxis_length += pow(faceBaxis[i], 2);
456 }
457
458 elementAaxis_length = sqrt(elementAaxis_length);
459 elementBaxis_length = sqrt(elementBaxis_length);
460 faceAaxis_length = sqrt(faceAaxis_length);
461 faceBaxis_length = sqrt(faceBaxis_length);
462
463 // Calculate the inner product of both the A-axis
464 // (i.e. Elemental A axis and face A axis)
465 for (i = 0; i < m_coordim; i++)
466 {
467 dotproduct1 += elementAaxis[i] * faceAaxis[i];
468 }
469
470 NekDouble norm =
471 fabs(dotproduct1) / elementAaxis_length / faceAaxis_length;
472 orientation = 0;
473
474 // if the innerproduct is equal to the (absolute value of the ) products
475 // of the lengths of both vectors, then, the coordinate systems will NOT
476 // be transposed
477 if (fabs(norm - 1.0) < NekConstants::kNekZeroTol)
478 {
479 // if the inner product is negative, both A-axis point
480 // in reverse direction
481 if (dotproduct1 < 0.0)
482 {
483 orientation += 2;
484 }
485
486 // calculate the inner product of both B-axis
487 for (i = 0; i < m_coordim; i++)
488 {
489 dotproduct2 += elementBaxis[i] * faceBaxis[i];
490 }
491
492 norm = fabs(dotproduct2) / elementBaxis_length / faceBaxis_length;
493
494 // check that both these axis are indeed parallel
495 ASSERTL1(fabs(norm - 1.0) < NekConstants::kNekZeroTol,
496 "These vectors should be parallel");
497
498 // if the inner product is negative, both B-axis point
499 // in reverse direction
500 if (dotproduct2 < 0.0)
501 {
502 orientation++;
503 }
504 }
505 // The coordinate systems are transposed
506 else
507 {
508 orientation = 4;
509
510 // Calculate the inner product between the elemental A-axis
511 // and the B-axis of the face (which are now the corresponding axis)
512 dotproduct1 = 0.0;
513 for (i = 0; i < m_coordim; i++)
514 {
515 dotproduct1 += elementAaxis[i] * faceBaxis[i];
516 }
517
518 norm = fabs(dotproduct1) / elementAaxis_length / faceBaxis_length;
519
520 // check that both these axis are indeed parallel
521 ASSERTL1(fabs(norm - 1.0) < NekConstants::kNekZeroTol,
522 "These vectors should be parallel");
523
524 // if the result is negative, both axis point in reverse
525 // directions
526 if (dotproduct1 < 0.0)
527 {
528 orientation += 2;
529 }
530
531 // Do the same for the other two corresponding axis
532 dotproduct2 = 0.0;
533 for (i = 0; i < m_coordim; i++)
534 {
535 dotproduct2 += elementBaxis[i] * faceAaxis[i];
536 }
537
538 norm = fabs(dotproduct2) / elementBaxis_length / faceAaxis_length;
539
540 // check that both these axis are indeed parallel
541 ASSERTL1(fabs(norm - 1.0) < NekConstants::kNekZeroTol,
542 "These vectors should be parallel");
543
544 if (dotproduct2 < 0.0)
545 {
546 orientation++;
547 }
548 }
549
550 orientation = orientation + 5;
551
553 "Orientation of triangular face (id = " +
554 std::to_string(m_faces[f]->GetGlobalID()) +
555 ") is inconsistent with face " + std::to_string(f) +
556 " of tet element (id = " + std::to_string(m_globalID) +
557 ") since Dir2 is aligned with Dir1. Mesh setup "
558 "needs investigation");
559
560 // Fill the m_forient array
561 m_forient[f] = (StdRegions::Orientation)orientation;
562 }
563}
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
int m_coordim
Coordinate dimension of this geometry object.
Definition Geometry.h:183
std::array< StdRegions::Orientation, kNfaces > m_forient
Definition TetGeom.h:116
static const int kNfaces
Definition TetGeom.h:52
static const int kNqfaces
Definition TetGeom.h:50
static const int kNtfaces
Definition TetGeom.h:51
static const int kNverts
Definition TriGeom.h:58
static const NekDouble kNekZeroTol
scalarT< T > sqrt(scalarT< T > in)
Definition scalar.hpp:290

References ASSERTL0, ASSERTL1, Nektar::StdRegions::eDir1FwdDir2_Dir2FwdDir1, Nektar::SpatialDomains::Geometry::GetGlobalID(), Nektar::SpatialDomains::Geometry::GetVertex(), Nektar::NekConstants::kNekZeroTol, kNfaces, kNqfaces, kNtfaces, Nektar::SpatialDomains::TriGeom::kNverts, Nektar::SpatialDomains::Geometry::m_coordim, m_faces, m_forient, Nektar::SpatialDomains::Geometry::m_globalID, m_verts, and tinysimd::sqrt().

Referenced by TetGeom().

◆ SetUpLocalEdges()

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

Definition at line 113 of file TetGeom.cpp.

114{
115
116 // find edge 0
117 int i, j;
118 unsigned int check;
119
120 // First set up the 3 bottom edges
121
122 if (m_faces[0]->GetEid(0) != m_faces[1]->GetEid(0))
123 {
124 std::ostringstream errstrm;
125 errstrm << "Local edge 0 (eid=" << m_faces[0]->GetEid(0);
126 errstrm << ") on face " << m_faces[0]->GetGlobalID();
127 errstrm << " must be the same as local edge 0 (eid="
128 << m_faces[1]->GetEid(0);
129 errstrm << ") on face " << m_faces[1]->GetGlobalID();
130 ASSERTL0(false, errstrm.str());
131 }
132
133 int faceConnected;
134 for (faceConnected = 1; faceConnected < 4; faceConnected++)
135 {
136 check = 0;
137 for (i = 0; i < 3; i++)
138 {
139 if ((m_faces[0])->GetEid(i) == (m_faces[faceConnected])->GetEid(0))
140 {
141 m_edges[faceConnected - 1] =
142 static_cast<SegGeom *>((m_faces[0])->GetEdge(i));
143 check++;
144 }
145 }
146
147 if (check < 1)
148 {
149 std::ostringstream errstrm;
150 errstrm << "Face 0 does not share an edge with first edge of "
151 "adjacent face. Faces ";
152 errstrm << (m_faces[0])->GetGlobalID() << ", "
153 << (m_faces[faceConnected])->GetGlobalID();
154 ASSERTL0(false, errstrm.str());
155 }
156 else if (check > 1)
157 {
158 std::ostringstream errstrm;
159 errstrm << "Connected faces share more than one edge. Faces ";
160 errstrm << (m_faces[0])->GetGlobalID() << ", "
161 << (m_faces[faceConnected])->GetGlobalID();
162 ASSERTL0(false, errstrm.str());
163 }
164 }
165
166 // Then, set up the 3 vertical edges
167 check = 0;
168 for (i = 0; i < 3; i++) // Set up the vertical edge :face(1) and face(3)
169 {
170 for (j = 0; j < 3; j++)
171 {
172 if ((m_faces[1])->GetEid(i) == (m_faces[3])->GetEid(j))
173 {
174 m_edges[3] = static_cast<SegGeom *>((m_faces[1])->GetEdge(i));
175 check++;
176 }
177 }
178 }
179 if (check < 1)
180 {
181 std::ostringstream errstrm;
182 errstrm << "Connected faces do not share an edge. Faces ";
183 errstrm << (m_faces[1])->GetGlobalID() << ", "
184 << (m_faces[3])->GetGlobalID();
185 ASSERTL0(false, errstrm.str());
186 }
187 else if (check > 1)
188 {
189 std::ostringstream errstrm;
190 errstrm << "Connected faces share more than one edge. Faces ";
191 errstrm << (m_faces[1])->GetGlobalID() << ", "
192 << (m_faces[3])->GetGlobalID();
193 ASSERTL0(false, errstrm.str());
194 }
195 // Set up vertical edges: face(1) through face(3)
196 for (faceConnected = 1; faceConnected < 3; faceConnected++)
197 {
198 check = 0;
199 for (i = 0; i < 3; i++)
200 {
201 for (j = 0; j < 3; j++)
202 {
203 if ((m_faces[faceConnected])->GetEid(i) ==
204 (m_faces[faceConnected + 1])->GetEid(j))
205 {
206 m_edges[faceConnected + 3] = static_cast<SegGeom *>(
207 (m_faces[faceConnected])->GetEdge(i));
208 check++;
209 }
210 }
211 }
212
213 if (check < 1)
214 {
215 std::ostringstream errstrm;
216 errstrm << "Connected faces do not share an edge. Faces ";
217 errstrm << (m_faces[faceConnected])->GetGlobalID() << ", "
218 << (m_faces[faceConnected + 1])->GetGlobalID();
219 ASSERTL0(false, errstrm.str());
220 }
221 else if (check > 1)
222 {
223 std::ostringstream errstrm;
224 errstrm << "Connected faces share more than one edge. Faces ";
225 errstrm << (m_faces[faceConnected])->GetGlobalID() << ", "
226 << (m_faces[faceConnected + 1])->GetGlobalID();
227 ASSERTL0(false, errstrm.str());
228 }
229 }
230}
int GetEid(int i) const
Get the ID of edge i of this object.
Definition Geometry.cpp:110

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

Referenced by TetGeom().

◆ SetUpLocalVertices()

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

Definition at line 232 of file TetGeom.cpp.

233{
234
235 // Set up the first 2 vertices (i.e. vertex 0,1)
236 if ((m_edges[0]->GetVid(0) == m_edges[1]->GetVid(0)) ||
237 (m_edges[0]->GetVid(0) == m_edges[1]->GetVid(1)))
238 {
239 m_verts[0] = m_edges[0]->GetVertex(1);
240 m_verts[1] = m_edges[0]->GetVertex(0);
241 }
242 else if ((m_edges[0]->GetVid(1) == m_edges[1]->GetVid(0)) ||
243 (m_edges[0]->GetVid(1) == m_edges[1]->GetVid(1)))
244 {
245 m_verts[0] = m_edges[0]->GetVertex(0);
246 m_verts[1] = m_edges[0]->GetVertex(1);
247 }
248 else
249 {
250 std::ostringstream errstrm;
251 errstrm << "Connected edges do not share a vertex. Edges ";
252 errstrm << m_edges[0]->GetGlobalID() << ", "
253 << m_edges[1]->GetGlobalID();
254 ASSERTL0(false, errstrm.str());
255 }
256
257 // set up the other bottom vertices (i.e. vertex 2)
258 for (int i = 1; i < 2; i++)
259 {
260 if (m_edges[i]->GetVid(0) == m_verts[i]->GetGlobalID())
261 {
262 m_verts[2] = m_edges[i]->GetVertex(1);
263 }
264 else if (m_edges[i]->GetVid(1) == m_verts[i]->GetGlobalID())
265 {
266 m_verts[2] = m_edges[i]->GetVertex(0);
267 }
268 else
269 {
270 std::ostringstream errstrm;
271 errstrm << "Connected edges do not share a vertex. Edges ";
272 errstrm << m_edges[i]->GetGlobalID() << ", "
273 << m_edges[i - 1]->GetGlobalID();
274 ASSERTL0(false, errstrm.str());
275 }
276 }
277
278 // set up top vertex
279 if (m_edges[3]->GetVid(0) == m_verts[0]->GetGlobalID())
280 {
281 m_verts[3] = m_edges[3]->GetVertex(1);
282 }
283 else
284 {
285 m_verts[3] = m_edges[3]->GetVertex(0);
286 }
287
288 // Check the other edges match up.
289 int check = 0;
290 for (int i = 4; i < 6; ++i)
291 {
292 if ((m_edges[i]->GetVid(0) == m_verts[i - 3]->GetGlobalID() &&
293 m_edges[i]->GetVid(1) == m_verts[3]->GetGlobalID()) ||
294 (m_edges[i]->GetVid(1) == m_verts[i - 3]->GetGlobalID() &&
295 m_edges[i]->GetVid(0) == m_verts[3]->GetGlobalID()))
296 {
297 check++;
298 }
299 }
300 if (check != 2)
301 {
302 std::ostringstream errstrm;
303 errstrm << "Connected edges do not share a vertex. Edges ";
304 errstrm << m_edges[3]->GetGlobalID() << ", "
305 << m_edges[2]->GetGlobalID();
306 ASSERTL0(false, errstrm.str());
307 }
308}

References ASSERTL0, Nektar::SpatialDomains::Geometry::GetGlobalID(), Nektar::SpatialDomains::Geometry::GetVid(), m_edges, and m_verts.

Referenced by TetGeom().

◆ SetUpXmap()

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

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

Definition at line 644 of file TetGeom.cpp.

645{
646 std::vector<int> tmp;
647 tmp.push_back(m_faces[0]->GetXmap()->GetTraceNcoeffs(0));
648 int order0 = *std::max_element(tmp.begin(), tmp.end());
649
650 tmp.clear();
651 tmp.push_back(order0);
652 tmp.push_back(m_faces[0]->GetXmap()->GetTraceNcoeffs(1));
653 tmp.push_back(m_faces[0]->GetXmap()->GetTraceNcoeffs(2));
654 int order1 = *std::max_element(tmp.begin(), tmp.end());
655
656 tmp.clear();
657 tmp.push_back(order0);
658 tmp.push_back(order1);
659 tmp.push_back(m_faces[1]->GetXmap()->GetTraceNcoeffs(1));
660 tmp.push_back(m_faces[1]->GetXmap()->GetTraceNcoeffs(2));
661 tmp.push_back(m_faces[3]->GetXmap()->GetTraceNcoeffs(1));
662 int order2 = *std::max_element(tmp.begin(), tmp.end());
663
664 std::array<LibUtilities::BasisKey, 3> basis = {
665 LibUtilities::BasisKey(
667 LibUtilities::PointsKey(order0 + 1,
669 LibUtilities::BasisKey(
671 LibUtilities::PointsKey(order1,
672 LibUtilities::eGaussRadauMAlpha1Beta0)),
673 LibUtilities::BasisKey(
675 LibUtilities::PointsKey(order2,
676 LibUtilities::eGaussRadauMAlpha2Beta0))};
677
678 m_xmap = GetStdTetFactory().CreateInstance(basis);
679}
StdRegions::StdExpansionSharedPtr m_xmap
mapping containing isoparametric transformation.
Definition Geometry.h:189
StdRegions::StdExpansionSharedPtr GetXmap() const
Return the mapping object Geometry::m_xmap that represents the coordinate transformation from standar...
Definition Geometry.h:439
@ eGaussLobattoLegendre
1D Gauss-Lobatto-Legendre quadrature points
Definition PointsType.h:51
@ eModified_B
Principle Modified Functions .
Definition BasisType.h:49
@ eModified_C
Principle Modified Functions .
Definition BasisType.h:50
@ eModified_A
Principle Modified Functions .
Definition BasisType.h:48
XmapFactory< StdRegions::StdTetExp, 3 > & GetStdTetFactory()
Definition TetGeom.cpp:53

References Nektar::LibUtilities::eGaussLobattoLegendre, Nektar::LibUtilities::eModified_A, Nektar::LibUtilities::eModified_B, Nektar::LibUtilities::eModified_C, Nektar::SpatialDomains::GetStdTetFactory(), Nektar::SpatialDomains::Geometry::GetXmap(), m_faces, and Nektar::SpatialDomains::Geometry::m_xmap.

Referenced by v_Setup().

◆ v_FillGeom()

void Nektar::SpatialDomains::TetGeom::v_FillGeom ( )
overrideprotectedvirtual

Put all quadrature information into face/edge structure and backward transform.

Note verts, edges, and faces are listed according to anticlockwise convention but points in _coeffs have to be in array format from left to right.

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 689 of file TetGeom.cpp.

690{
691 if (m_state == ePtsFilled)
692 {
693 return;
694 }
695
696 int i, j, k;
697
698 for (i = 0; i < kNfaces; i++)
699 {
700 m_faces[i]->FillGeom();
701
702 int nFaceCoeffs = m_faces[i]->GetXmap()->GetNcoeffs();
703
704 Array<OneD, unsigned int> mapArray(nFaceCoeffs);
705 Array<OneD, int> signArray(nFaceCoeffs);
706
707 if (m_forient[i] < 9)
708 {
709 m_xmap->GetTraceToElementMap(
710 i, mapArray, signArray, m_forient[i],
711 m_faces[i]->GetXmap()->GetTraceNcoeffs(0),
712 m_faces[i]->GetXmap()->GetTraceNcoeffs(1));
713 }
714 else
715 {
716 m_xmap->GetTraceToElementMap(
717 i, mapArray, signArray, m_forient[i],
718 m_faces[i]->GetXmap()->GetTraceNcoeffs(1),
719 m_faces[i]->GetXmap()->GetTraceNcoeffs(0));
720 }
721
722 for (j = 0; j < m_coordim; j++)
723 {
724 const Array<OneD, const NekDouble> &coeffs =
725 m_faces[i]->GetCoeffs(j);
726
727 for (k = 0; k < nFaceCoeffs; k++)
728 {
729 NekDouble v = signArray[k] * coeffs[k];
730 m_coeffs[j][mapArray[k]] = v;
731 }
732 }
733 }
734
736}
GeomState m_state
Enumeration to dictate whether coefficients are filled.
Definition Geometry.h:191
Array< OneD, Array< OneD, NekDouble > > m_coeffs
Array containing expansion coefficients of m_xmap.
Definition Geometry.h:201
@ ePtsFilled
Geometric information has been generated.

References Nektar::SpatialDomains::ePtsFilled, Nektar::SpatialDomains::Geometry::GetXmap(), kNfaces, Nektar::SpatialDomains::Geometry::m_coeffs, Nektar::SpatialDomains::Geometry::m_coordim, m_faces, m_forient, Nektar::SpatialDomains::Geometry::m_state, and Nektar::SpatialDomains::Geometry::m_xmap.

Referenced by v_GenGeomFactors().

◆ v_GenGeomFactors()

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

Generate the geometry factors for this element.

Implements Nektar::SpatialDomains::Geometry.

Definition at line 592 of file TetGeom.cpp.

593{
594 if (!m_setupState)
595 {
597 }
598
600 {
601 GeomType Gtype = eRegular;
602
603 v_FillGeom();
604
605 // check to see if expansions are linear
606 m_straightEdge = 1;
607 if (m_xmap->GetBasisNumModes(0) != 2 ||
608 m_xmap->GetBasisNumModes(1) != 2 ||
609 m_xmap->GetBasisNumModes(2) != 2)
610 {
611 Gtype = eDeformed;
612 m_straightEdge = 0;
613 }
614
615 if (Gtype == eRegular)
616 {
617 m_isoParameter = Array<OneD, Array<OneD, NekDouble>>(3);
618 for (int i = 0; i < 3; ++i)
619 {
620 m_isoParameter[i] = Array<OneD, NekDouble>(4, 0.);
621 NekDouble A = (*m_verts[0])(i);
622 NekDouble B = (*m_verts[1])(i);
623 NekDouble C = (*m_verts[2])(i);
624 NekDouble D = (*m_verts[3])(i);
625 m_isoParameter[i][0] = 0.5 * (-A + B + C + D);
626
627 m_isoParameter[i][1] = 0.5 * (-A + B); // xi1
628 m_isoParameter[i][2] = 0.5 * (-A + C); // xi2
629 m_isoParameter[i][3] = 0.5 * (-A + D); // xi3
630 }
632 }
633
635 Gtype, m_coordim, m_xmap, m_coeffs);
637 }
638}
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
bool m_setupState
Wether or not the setup routines have been run.
Definition Geometry.h:193
Array< OneD, Array< OneD, NekDouble > > m_isoParameter
Definition Geometry.h:204
GeomState m_geomFactorsState
State of the geometric factors.
Definition Geometry.h:187
GeomFactorsSharedPtr m_geomFactors
Geometric factors.
Definition Geometry.h:185
void v_FillGeom() override
Put all quadrature information into face/edge structure and backward transform.
Definition TetGeom.cpp:689
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.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::SpatialDomains::eDeformed, Nektar::SpatialDomains::ePtsFilled, Nektar::SpatialDomains::eRegular, 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, m_verts, Nektar::SpatialDomains::Geometry::m_xmap, Nektar::SpatialDomains::Geometry3D::v_CalculateInverseIsoParam(), v_FillGeom(), and v_Setup().

◆ v_GetDir()

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

78{
79 if (faceidx == 0)
80 {
81 return facedir;
82 }
83 else if (faceidx == 1)
84 {
85 return 2 * facedir;
86 }
87 else
88 {
89 return 1 + facedir;
90 }
91}

◆ v_GetEdge()

Geometry1D * Nektar::SpatialDomains::TetGeom::v_GetEdge ( const int  i) const
inlinefinalprotectedvirtual

Returns edge i of this object.

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 92 of file TetGeom.h.

93 {
94 return static_cast<Geometry1D *>(m_edges[i]);
95 }

References m_edges.

◆ v_GetEdgeFaceMap()

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

104{
105 return EdgeFaceConnectivity[i][j];
106}
static const unsigned int EdgeFaceConnectivity[6][2]
Definition TetGeom.h:127

References EdgeFaceConnectivity.

◆ v_GetEdgeNormalToFaceVert()

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

109{
110 return EdgeNormalToFaceVert[i][j];
111}
static const unsigned int EdgeNormalToFaceVert[4][3]
Definition TetGeom.h:128

References EdgeNormalToFaceVert.

◆ v_GetEorient()

StdRegions::Orientation Nektar::SpatialDomains::TetGeom::v_GetEorient ( const int  i) const
inlinefinalprotectedvirtual

Returns the orientation of edge i with respect to the ordering of edges in the standard element.

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 102 of file TetGeom.h.

103 {
104 return m_eorient[i];
105 }

References m_eorient.

◆ v_GetFace()

Geometry2D * Nektar::SpatialDomains::TetGeom::v_GetFace ( const int  i) const
inlinefinalprotectedvirtual

Returns face i of this object.

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 97 of file TetGeom.h.

98 {
99 return static_cast<Geometry2D *>(m_faces[i]);
100 }

References m_faces.

◆ v_GetForient()

StdRegions::Orientation Nektar::SpatialDomains::TetGeom::v_GetForient ( const int  i) const
inlinefinalprotectedvirtual

Returns the orientation of face i with respect to the ordering of faces in the standard element.

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 107 of file TetGeom.h.

108 {
109 return m_forient[i];
110 }

References m_forient.

◆ v_GetNumEdges()

int Nektar::SpatialDomains::TetGeom::v_GetNumEdges ( ) const
inlinefinalprotectedvirtual

Get the number of edges of this object.

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 77 of file TetGeom.h.

78 {
79 return kNedges;
80 }

References kNedges.

◆ v_GetNumFaces()

int Nektar::SpatialDomains::TetGeom::v_GetNumFaces ( ) const
inlinefinalprotectedvirtual

Get the number of faces of this object.

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 82 of file TetGeom.h.

83 {
84 return kNfaces;
85 }

References kNfaces.

◆ v_GetNumVerts()

int Nektar::SpatialDomains::TetGeom::v_GetNumVerts ( ) const
inlinefinalprotectedvirtual

Get the number of vertices of this object.

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 72 of file TetGeom.h.

73 {
74 return kNverts;
75 }
static const int kNverts
Definition TetGeom.h:48

References kNverts.

◆ v_GetVertex()

PointGeom * Nektar::SpatialDomains::TetGeom::v_GetVertex ( const int  i) const
inlinefinalprotectedvirtual

Returns vertex i of this object.

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 87 of file TetGeom.h.

88 {
89 return m_verts[i];
90 }

References m_verts.

◆ v_GetVertexEdgeMap()

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

94{
95 return VertexEdgeConnectivity[i][j];
96}
static const unsigned int VertexEdgeConnectivity[4][3]
Definition TetGeom.h:125

References VertexEdgeConnectivity.

◆ v_GetVertexFaceMap()

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

99{
100 return VertexFaceConnectivity[i][j];
101}
static const unsigned int VertexFaceConnectivity[4][3]
Definition TetGeom.h:126

References VertexFaceConnectivity.

◆ v_Reset()

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

566{
567 Geometry::v_Reset(curvedEdges, curvedFaces);
568
569 for (int i = 0; i < 4; ++i)
570 {
571 m_faces[i]->Reset(curvedEdges, curvedFaces);
572 }
573}
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:372

References m_faces, and Nektar::SpatialDomains::Geometry::v_Reset().

◆ v_Setup()

void Nektar::SpatialDomains::TetGeom::v_Setup ( )
overrideprotectedvirtual

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 575 of file TetGeom.cpp.

576{
577 if (!m_setupState)
578 {
579 for (int i = 0; i < 4; ++i)
580 {
581 m_faces[i]->Setup();
582 }
583 SetUpXmap();
584 SetUpCoeffs(m_xmap->GetNcoeffs());
585 m_setupState = true;
586 }
587}
void SetUpCoeffs(const int nCoeffs)
Initialise the Geometry::m_coeffs array.
Definition Geometry.h:704
void SetUpXmap()
Set up the m_xmap object by determining the order of each direction from derived faces.
Definition TetGeom.cpp:644

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

Definition at line 127 of file TetGeom.h.

Referenced by v_GetEdgeFaceMap().

◆ EdgeNormalToFaceVert

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

Definition at line 128 of file TetGeom.h.

Referenced by v_GetEdgeNormalToFaceVert().

◆ kNedges

const int Nektar::SpatialDomains::TetGeom::kNedges = 6
static

Definition at line 49 of file TetGeom.h.

Referenced by SetUpEdgeOrientation(), and v_GetNumEdges().

◆ kNfaces

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

◆ kNfacets

const int Nektar::SpatialDomains::TetGeom::kNfacets = kNfaces
static

Definition at line 53 of file TetGeom.h.

◆ kNqfaces

const int Nektar::SpatialDomains::TetGeom::kNqfaces = 0
static

◆ kNtfaces

const int Nektar::SpatialDomains::TetGeom::kNtfaces = 4
static

◆ kNverts

const int Nektar::SpatialDomains::TetGeom::kNverts = 4
static

Definition at line 48 of file TetGeom.h.

Referenced by v_GetNumVerts().

◆ m_edges

std::array<SegGeom *, kNedges> Nektar::SpatialDomains::TetGeom::m_edges
protected

Definition at line 113 of file TetGeom.h.

Referenced by SetUpEdgeOrientation(), SetUpLocalEdges(), SetUpLocalVertices(), and v_GetEdge().

◆ m_eorient

std::array<StdRegions::Orientation, kNedges> Nektar::SpatialDomains::TetGeom::m_eorient
protected

Definition at line 115 of file TetGeom.h.

Referenced by SetUpEdgeOrientation(), and v_GetEorient().

◆ m_faces

std::array<TriGeom *, kNfaces> Nektar::SpatialDomains::TetGeom::m_faces
protected

◆ m_forient

std::array<StdRegions::Orientation, kNfaces> Nektar::SpatialDomains::TetGeom::m_forient
protected

Definition at line 116 of file TetGeom.h.

Referenced by SetUpFaceOrientation(), v_FillGeom(), and v_GetForient().

◆ m_verts

std::array<PointGeom *, kNverts> Nektar::SpatialDomains::TetGeom::m_verts
protected

◆ VertexEdgeConnectivity

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

Definition at line 125 of file TetGeom.h.

Referenced by v_GetVertexEdgeMap().

◆ VertexFaceConnectivity

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

Definition at line 126 of file TetGeom.h.

Referenced by v_GetVertexFaceMap().

◆ XMLElementType

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

Definition at line 54 of file TetGeom.h.