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

#include <TriGeom.h>

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

Public Member Functions

 TriGeom ()
 
 TriGeom (const TriGeom &in)
 
 TriGeom (const int id, const SegGeomSharedPtr edges[], const CurveSharedPtr curve=CurveSharedPtr())
 
 ~TriGeom ()
 
- Public Member Functions inherited from Nektar::SpatialDomains::Geometry2D
 Geometry2D ()
 
 Geometry2D (const int coordim, CurveSharedPtr curve)
 
virtual ~Geometry2D ()
 
CurveSharedPtr GetCurve ()
 
- Public Member Functions inherited from Nektar::SpatialDomains::Geometry
 Geometry ()
 Default constructor. More...
 
 Geometry (int coordim)
 Constructor when supplied a coordinate dimension. More...
 
virtual ~Geometry ()
 Default destructor. More...
 
int GetCoordim () const
 Return the coordinate dimension of this object (i.e. the dimension of the space in which this object is embedded). More...
 
void SetCoordim (int coordim)
 Sets the coordinate dimension of this object (i.e. the dimension of the space in which this object is embedded). More...
 
GeomFactorsSharedPtr GetGeomFactors ()
 Get the geometric factors for this object, generating them if required. More...
 
GeomFactorsSharedPtr GetRefGeomFactors (const Array< OneD, const LibUtilities::BasisSharedPtr > &tbasis)
 
GeomFactorsSharedPtr GetMetricInfo ()
 Get the geometric factors for this object. More...
 
LibUtilities::ShapeType GetShapeType (void)
 Get the geometric shape type of this object. More...
 
int GetGlobalID (void) const
 Get the ID of this object. More...
 
void SetGlobalID (int globalid)
 Set the ID of this object. More...
 
int GetVid (int i) const
 Get the ID of vertex i of this object. More...
 
int GetEid (int i) const
 Get the ID of edge i of this object. More...
 
int GetFid (int i) const
 Get the ID of face i of this object. More...
 
int GetTid (int i) const
 Get the ID of trace i of this object. More...
 
PointGeomSharedPtr GetVertex (int i) const
 Returns vertex i of this object. More...
 
Geometry1DSharedPtr GetEdge (int i) const
 Returns edge i of this object. More...
 
Geometry2DSharedPtr GetFace (int i) const
 Returns face i of this object. More...
 
StdRegions::Orientation GetEorient (const int i) const
 Returns the orientation of edge i with respect to the ordering of edges in the standard element. More...
 
StdRegions::Orientation GetForient (const int i) const
 Returns the orientation of face i with respect to the ordering of faces in the standard element. More...
 
int GetNumVerts () const
 Get the number of vertices of this object. More...
 
int GetNumEdges () const
 Get the number of edges of this object. More...
 
int GetNumFaces () const
 Get the number of faces of this object. More...
 
int GetShapeDim () const
 Get the object's shape dimension. More...
 
StdRegions::StdExpansionSharedPtr GetXmap () const
 Return the mapping object Geometry::m_xmap that represents the coordinate transformation from standard element to physical element. More...
 
const Array< OneD, const NekDouble > & GetCoeffs (const int i) const
 Return the coefficients of the transformation Geometry::m_xmap in coordinate direction i. More...
 
void FillGeom ()
 Populate the coordinate mapping Geometry::m_coeffs information from any children geometry elements. More...
 
std::array< NekDouble, 6 > GetBoundingBox ()
 Generates the bounding box for the element. More...
 
void ClearBoundingBox ()
 
bool ContainsPoint (const Array< OneD, const NekDouble > &gloCoord, NekDouble tol=0.0)
 Determine whether an element contains a particular Cartesian coordinate \((x,y,z)\). More...
 
bool ContainsPoint (const Array< OneD, const NekDouble > &gloCoord, Array< OneD, NekDouble > &locCoord, NekDouble tol)
 Determine whether an element contains a particular Cartesian coordinate \((x,y,z)\). More...
 
bool ContainsPoint (const Array< OneD, const NekDouble > &gloCoord, Array< OneD, NekDouble > &locCoord, NekDouble tol, NekDouble &dist)
 Determine whether an element contains a particular Cartesian coordinate \(\vec{x} = (x,y,z)\). More...
 
NekDouble GetLocCoords (const Array< OneD, const NekDouble > &coords, Array< OneD, NekDouble > &Lcoords)
 Determine the local collapsed coordinates that correspond to a given Cartesian coordinate for this geometry object. More...
 
NekDouble GetCoord (const int i, const Array< OneD, const NekDouble > &Lcoord)
 Given local collapsed coordinate Lcoord, return the value of physical coordinate in direction i. More...
 
int PreliminaryCheck (const Array< OneD, const NekDouble > &gloCoord)
 A fast and robust check if a given global coord is outside of a deformed element. For regular elements, this check is unnecessary. More...
 
bool MinMaxCheck (const Array< OneD, const NekDouble > &gloCoord)
 Check if given global coord is within the BoundingBox of the element. More...
 
bool ClampLocCoords (Array< OneD, NekDouble > &locCoord, NekDouble tol=std::numeric_limits< NekDouble >::epsilon())
 Clamp local coords to be within standard regions [-1, 1]^dim. More...
 
NekDouble FindDistance (const Array< OneD, const NekDouble > &xs, Array< OneD, NekDouble > &xi)
 
int GetVertexEdgeMap (int i, int j) const
 Returns the standard element edge IDs that are connected to a given vertex. More...
 
int GetVertexFaceMap (int i, int j) const
 Returns the standard element face IDs that are connected to a given vertex. More...
 
int GetEdgeFaceMap (int i, int j) const
 Returns the standard element edge IDs that are connected to a given face. More...
 
int GetEdgeNormalToFaceVert (int i, int j) const
 Returns the standard lement edge IDs that are normal to a given face vertex. More...
 
int GetDir (const int i, const int j=0) const
 Returns the element coordinate direction corresponding to a given face coordinate direction. More...
 
void Reset (CurveMap &curvedEdges, CurveMap &curvedFaces)
 Reset this geometry object: unset the current state, zero Geometry::m_coeffs and remove allocated GeomFactors. More...
 
void Setup ()
 

Static Public Member Functions

static StdRegions::Orientation GetFaceOrientation (const TriGeom &face1, const TriGeom &face2, bool doRot, int dir, NekDouble angle, NekDouble tol)
 
static StdRegions::Orientation GetFaceOrientation (const PointGeomVector &face1, const PointGeomVector &face2, bool doRot, int dir, NekDouble angle, NekDouble tol)
 

Static Public Attributes

static const int kNedges = 3
 Get the orientation of face1. More...
 
static const int kNverts = 3
 
- Static Public Attributes inherited from Nektar::SpatialDomains::Geometry2D
static const int kDim = 2
 

Protected Member Functions

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

Private Member Functions

void SetUpXmap ()
 

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::Geometry2D
PointGeomVector m_verts
 
SegGeomVector m_edges
 
std::vector< StdRegions::Orientationm_eorient
 
CurveSharedPtr m_curve
 
Array< OneD, int > m_manifold
 
Array< OneD, Array< OneD, NekDouble > > m_edgeNormal
 
- Protected Attributes inherited from Nektar::SpatialDomains::Geometry
int m_coordim
 Coordinate dimension of this geometry object. More...
 
GeomFactorsSharedPtr m_geomFactors
 Geometric factors. More...
 
GeomState m_geomFactorsState
 State of the geometric factors. More...
 
StdRegions::StdExpansionSharedPtr m_xmap
 \(\chi\) mapping containing isoparametric transformation. More...
 
GeomState m_state
 Enumeration to dictate whether coefficients are filled. More...
 
bool m_setupState
 Wether or not the setup routines have been run. More...
 
GeomType m_geomType
 Type of geometry. More...
 
LibUtilities::ShapeType m_shapeType
 Type of shape. More...
 
int m_globalID
 Global ID. More...
 
Array< OneD, Array< OneD, NekDouble > > m_coeffs
 Array containing expansion coefficients of m_xmap. More...
 
Array< OneD, NekDoublem_boundingBox
 Array containing bounding box. More...
 
Array< OneD, Array< OneD, NekDouble > > m_isoParameter
 
Array< OneD, Array< OneD, NekDouble > > m_invIsoParam
 
bool m_straightEdge
 
- Static Protected Attributes inherited from Nektar::SpatialDomains::Geometry
static GeomFactorsVector m_regGeomFactorsManager
 

Detailed Description

Definition at line 61 of file TriGeom.h.

Constructor & Destructor Documentation

◆ TriGeom() [1/3]

Nektar::SpatialDomains::TriGeom::TriGeom ( )

◆ TriGeom() [2/3]

Nektar::SpatialDomains::TriGeom::TriGeom ( const TriGeom in)

Definition at line 85 of file TriGeom.cpp.

85 : Geometry2D(in)
86{
87 // From Geometry
88 m_shapeType = in.m_shapeType;
89
90 // From TriFaceComponent
91 m_globalID = in.m_globalID;
92
93 // From TriGeom
94 m_verts = in.m_verts;
95 m_edges = in.m_edges;
96 for (int i = 0; i < kNedges; i++)
97 {
98 m_eorient[i] = in.m_eorient[i];
99 }
100}
std::vector< StdRegions::Orientation > m_eorient
Definition: Geometry2D.h:84
static const int kNedges
Get the orientation of face1.
Definition: TriGeom.h:72

References kNedges, Nektar::SpatialDomains::Geometry2D::m_edges, Nektar::SpatialDomains::Geometry2D::m_eorient, Nektar::SpatialDomains::Geometry::m_globalID, Nektar::SpatialDomains::Geometry::m_shapeType, and Nektar::SpatialDomains::Geometry2D::m_verts.

◆ TriGeom() [3/3]

Nektar::SpatialDomains::TriGeom::TriGeom ( const int  id,
const SegGeomSharedPtr  edges[],
const CurveSharedPtr  curve = CurveSharedPtr() 
)

Definition at line 56 of file TriGeom.cpp.

58 : Geometry2D(edges[0]->GetVertex(0)->GetCoordim(), curve)
59{
60 int j;
61
63 m_globalID = id;
64
65 // Copy the edge shared pointers.
66 m_edges.insert(m_edges.begin(), edges, edges + TriGeom::kNedges);
67 m_eorient.resize(kNedges);
68
69 for (j = 0; j < kNedges; ++j)
70 {
71 m_eorient[j] =
72 SegGeom::GetEdgeOrientation(*edges[j], *edges[(j + 1) % kNedges]);
73 m_verts.push_back(
74 edges[j]->GetVertex(m_eorient[j] == StdRegions::eForwards ? 0 : 1));
75 }
76
80
81 m_coordim = edges[0]->GetVertex(0)->GetCoordim();
82 ASSERTL0(m_coordim > 1, "Cannot call function with dim == 1");
83}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
PointGeomSharedPtr GetVertex(int i) const
Returns vertex i of this object.
Definition: Geometry.h:358
int GetCoordim() const
Return the coordinate dimension of this object (i.e. the dimension of the space in which this object ...
Definition: Geometry.h:284
int m_coordim
Coordinate dimension of this geometry object.
Definition: Geometry.h:186
static StdRegions::Orientation GetEdgeOrientation(const SegGeom &edge1, const SegGeom &edge2)
Get the orientation of edge1.
Definition: SegGeom.cpp:196

References ASSERTL0, Nektar::StdRegions::eBackwards, Nektar::StdRegions::eForwards, Nektar::LibUtilities::eTriangle, Nektar::SpatialDomains::SegGeom::GetEdgeOrientation(), Nektar::SpatialDomains::Geometry::GetVertex(), kNedges, Nektar::SpatialDomains::Geometry::m_coordim, Nektar::SpatialDomains::Geometry2D::m_edges, Nektar::SpatialDomains::Geometry2D::m_eorient, Nektar::SpatialDomains::Geometry::m_globalID, Nektar::SpatialDomains::Geometry::m_shapeType, and Nektar::SpatialDomains::Geometry2D::m_verts.

◆ ~TriGeom()

Nektar::SpatialDomains::TriGeom::~TriGeom ( )

Definition at line 102 of file TriGeom.cpp.

103{
104}

Member Function Documentation

◆ GetFaceOrientation() [1/2]

StdRegions::Orientation Nektar::SpatialDomains::TriGeom::GetFaceOrientation ( const PointGeomVector face1,
const PointGeomVector face2,
bool  doRot,
int  dir,
NekDouble  angle,
NekDouble  tol 
)
static

Definition at line 127 of file TriGeom.cpp.

130{
131 int i, j, vmap[3] = {-1, -1, -1};
132
133 if (doRot)
134 {
135 PointGeom rotPt;
136
137 for (i = 0; i < 3; ++i)
138 {
139 rotPt.Rotate((*face1[i]), dir, angle);
140 for (j = 0; j < 3; ++j)
141 {
142 if (rotPt.dist(*face2[j]) < tol)
143 {
144 vmap[j] = i;
145 break;
146 }
147 }
148 }
149 }
150 else
151 {
152
153 NekDouble x, y, z, x1, y1, z1, cx = 0.0, cy = 0.0, cz = 0.0;
154
155 // For periodic faces, we calculate the vector between the centre
156 // points of the two faces. (For connected faces this will be
157 // zero). We can then use this to determine alignment later in the
158 // algorithm.
159 for (i = 0; i < 3; ++i)
160 {
161 cx += (*face2[i])(0) - (*face1[i])(0);
162 cy += (*face2[i])(1) - (*face1[i])(1);
163 cz += (*face2[i])(2) - (*face1[i])(2);
164 }
165 cx /= 3;
166 cy /= 3;
167 cz /= 3;
168
169 // Now construct a mapping which takes us from the vertices of one
170 // face to the other. That is, vertex j of face2 corresponds to
171 // vertex vmap[j] of face1.
172 for (i = 0; i < 3; ++i)
173 {
174 x = (*face1[i])(0);
175 y = (*face1[i])(1);
176 z = (*face1[i])(2);
177 for (j = 0; j < 3; ++j)
178 {
179 x1 = (*face2[j])(0) - cx;
180 y1 = (*face2[j])(1) - cy;
181 z1 = (*face2[j])(2) - cz;
182 if (sqrt((x1 - x) * (x1 - x) + (y1 - y) * (y1 - y) +
183 (z1 - z) * (z1 - z)) < 1e-8)
184 {
185 vmap[j] = i;
186 break;
187 }
188 }
189 }
190 }
191
192 if (vmap[1] == (vmap[0] + 1) % 3)
193 {
194 switch (vmap[0])
195 {
196 case 0:
198 break;
199 case 1:
201 break;
202 case 2:
204 break;
205 }
206 }
207 else
208 {
209 switch (vmap[0])
210 {
211 case 0:
213 break;
214 case 1:
216 break;
217 case 2:
219 break;
220 }
221 }
222
223 ASSERTL0(false, "Unable to determine triangle orientation");
225}
std::vector< double > z(NPUPPER)
double NekDouble
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294

References ASSERTL0, Nektar::SpatialDomains::PointGeom::dist(), Nektar::StdRegions::eDir1BwdDir1_Dir2BwdDir2, Nektar::StdRegions::eDir1BwdDir1_Dir2FwdDir2, Nektar::StdRegions::eDir1BwdDir2_Dir2BwdDir1, Nektar::StdRegions::eDir1FwdDir1_Dir2FwdDir2, Nektar::StdRegions::eDir1FwdDir2_Dir2BwdDir1, Nektar::StdRegions::eDir1FwdDir2_Dir2FwdDir1, Nektar::StdRegions::eNoOrientation, Nektar::SpatialDomains::PointGeom::Rotate(), tinysimd::sqrt(), and Nektar::UnitTests::z().

◆ GetFaceOrientation() [2/2]

StdRegions::Orientation Nektar::SpatialDomains::TriGeom::GetFaceOrientation ( const TriGeom face1,
const TriGeom face2,
bool  doRot,
int  dir,
NekDouble  angle,
NekDouble  tol 
)
static

Definition at line 117 of file TriGeom.cpp.

122{
123 return GetFaceOrientation(face1.m_verts, face2.m_verts, doRot, dir, angle,
124 tol);
125}
static StdRegions::Orientation GetFaceOrientation(const TriGeom &face1, const TriGeom &face2, bool doRot, int dir, NekDouble angle, NekDouble tol)
Definition: TriGeom.cpp:117

References GetFaceOrientation(), and Nektar::SpatialDomains::Geometry2D::m_verts.

Referenced by Nektar::MultiRegions::DisContField::FindPeriodicTraces(), and GetFaceOrientation().

◆ SetUpXmap()

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

Definition at line 559 of file TriGeom.cpp.

560{
561 int order0 = m_edges[0]->GetXmap()->GetBasis(0)->GetNumModes();
562 int order1 =
563 max(order0, max(m_edges[1]->GetXmap()->GetBasis(0)->GetNumModes(),
564 m_edges[2]->GetXmap()->GetBasis(0)->GetNumModes()));
565
566 const LibUtilities::BasisKey B0(
568 LibUtilities::PointsKey(order0 + 1,
570 const LibUtilities::BasisKey B1(
572 LibUtilities::PointsKey(order1, LibUtilities::eGaussRadauMAlpha1Beta0));
573
575}
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
StdRegions::StdExpansionSharedPtr m_xmap
mapping containing isoparametric transformation.
Definition: Geometry.h:192
StdRegions::StdExpansionSharedPtr GetXmap() const
Return the mapping object Geometry::m_xmap that represents the coordinate transformation from standar...
Definition: Geometry.h:436
@ eGaussLobattoLegendre
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:53
@ eModified_B
Principle Modified Functions .
Definition: BasisType.h:51
@ eModified_A
Principle Modified Functions .
Definition: BasisType.h:50

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::LibUtilities::eGaussLobattoLegendre, Nektar::LibUtilities::eModified_A, Nektar::LibUtilities::eModified_B, Nektar::SpatialDomains::Geometry::GetXmap(), Nektar::SpatialDomains::Geometry2D::m_edges, and Nektar::SpatialDomains::Geometry::m_xmap.

Referenced by v_Reset(), and v_Setup().

◆ v_FillGeom()

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

Note verts and edges 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 309 of file TriGeom.cpp.

310{
311 // check to see if geometry structure is already filled
312 if (m_state == ePtsFilled)
313 {
314 return;
315 }
316
317 int i, j, k;
318 int nEdgeCoeffs = m_xmap->GetTraceNcoeffs(0);
319
320 if (m_curve)
321 {
322 int pdim = LibUtilities::PointsManager()[LibUtilities::PointsKey(
323 2, m_curve->m_ptype)]
324 ->GetPointsDim();
325
326 // Deal with 2D points type separately
327 // (e.g. electrostatic or Fekete points) to 1D tensor
328 // product.
329 if (pdim == 2)
330 {
331 int N = m_curve->m_points.size();
332 int nEdgePts =
333 (-1 + (int)sqrt(static_cast<NekDouble>(8 * N + 1))) / 2;
334
335 ASSERTL0(nEdgePts * (nEdgePts + 1) / 2 == N,
336 "NUMPOINTS should be a triangle number for"
337 " triangle curved face " +
338 boost::lexical_cast<string>(m_globalID));
339
340 // Sanity check 1: are curved vertices consistent with
341 // triangle vertices?
342 for (i = 0; i < 3; ++i)
343 {
344 NekDouble dist = m_verts[i]->dist(*(m_curve->m_points[i]));
346 {
347 std::stringstream ss;
348 ss << "Curved vertex " << i << " of triangle " << m_globalID
349 << " is separated from expansion vertex by"
350 << " more than " << NekConstants::kVertexTheSameDouble
351 << " (dist = " << dist << ")";
352 NEKERROR(ErrorUtil::ewarning, ss.str().c_str());
353 }
354 }
355
356 // Sanity check 2: are curved edges from the face curvature
357 // consistent with curved edges?
358 for (i = 0; i < kNedges; ++i)
359 {
360 CurveSharedPtr edgeCurve = m_edges[i]->GetCurve();
361
362 ASSERTL0(edgeCurve->m_points.size() == nEdgePts,
363 "Number of edge points does not correspond "
364 "to number of face points in triangle " +
365 boost::lexical_cast<string>(m_globalID));
366
367 const int offset = 3 + i * (nEdgePts - 2);
368 NekDouble maxDist = 0.0;
369
370 // Account for different ordering of nodal coordinates
371 // vs. Cartesian ordering of element.
373
374 if (i == 2)
375 {
376 orient = orient == StdRegions::eForwards
379 }
380
381 if (orient == StdRegions::eForwards)
382 {
383 for (j = 0; j < nEdgePts - 2; ++j)
384 {
385 NekDouble dist = m_curve->m_points[offset + j]->dist(
386 *(edgeCurve->m_points[j + 1]));
387 maxDist = dist > maxDist ? dist : maxDist;
388 }
389 }
390 else
391 {
392 for (j = 0; j < nEdgePts - 2; ++j)
393 {
394 NekDouble dist = m_curve->m_points[offset + j]->dist(
395 *(edgeCurve->m_points[nEdgePts - 2 - j]));
396 maxDist = dist > maxDist ? dist : maxDist;
397 }
398 }
399
401 {
402 std::stringstream ss;
403 ss << "Curved edge " << i << " of triangle " << m_globalID
404 << " has a point separated from edge interior"
405 << " points by more than "
407 << " (maxdist = " << maxDist << ")";
408 NEKERROR(ErrorUtil::ewarning, ss.str().c_str());
409 }
410 }
411
412 const LibUtilities::PointsKey P0(
414 const LibUtilities::PointsKey P1(
415 nEdgePts, LibUtilities::eGaussRadauMAlpha1Beta0);
416 const LibUtilities::BasisKey T0(LibUtilities::eOrtho_A, nEdgePts,
417 P0);
418 const LibUtilities::BasisKey T1(LibUtilities::eOrtho_B, nEdgePts,
419 P1);
420 Array<OneD, NekDouble> phys(
421 max(nEdgePts * nEdgePts, m_xmap->GetTotPoints()));
422 Array<OneD, NekDouble> tmp(nEdgePts * nEdgePts);
423
424 for (i = 0; i < m_coordim; ++i)
425 {
426 // Create a StdNodalTriExp.
429 AllocateSharedPtr(T0, T1, m_curve->m_ptype);
430
431 for (j = 0; j < N; ++j)
432 {
433 phys[j] = (m_curve->m_points[j]->GetPtr())[i];
434 }
435
436 t->BwdTrans(phys, tmp);
437
438 // Interpolate points to standard region.
440 P0, P1, tmp, m_xmap->GetBasis(0)->GetPointsKey(),
441 m_xmap->GetBasis(1)->GetPointsKey(), phys);
442
443 // Forwards transform to get coefficient space.
444 m_xmap->FwdTrans(phys, m_coeffs[i]);
445 }
446 }
447 else if (pdim == 1)
448 {
449 int npts = m_curve->m_points.size();
450 int nEdgePts = (int)sqrt(static_cast<NekDouble>(npts));
451 Array<OneD, NekDouble> tmp(npts);
452 Array<OneD, NekDouble> phys(m_xmap->GetTotPoints());
453 LibUtilities::PointsKey curveKey(nEdgePts, m_curve->m_ptype);
454
455 // Sanity checks:
456 // - Curved faces should have square number of points;
457 // - Each edge should have sqrt(npts) points.
458 ASSERTL0(nEdgePts * nEdgePts == npts,
459 "NUMPOINTS should be a square number for"
460 " triangle " +
461 boost::lexical_cast<string>(m_globalID));
462
463 for (i = 0; i < kNedges; ++i)
464 {
465 ASSERTL0(m_edges[i]->GetXmap()->GetNcoeffs() == nEdgePts,
466 "Number of edge points does not correspond to "
467 "number of face points in triangle " +
468 boost::lexical_cast<string>(m_globalID));
469 }
470
471 for (i = 0; i < m_coordim; ++i)
472 {
473 for (j = 0; j < npts; ++j)
474 {
475 tmp[j] = (m_curve->m_points[j]->GetPtr())[i];
476 }
477
478 // Interpolate curve points to standard triangle
479 // points.
480 LibUtilities::Interp2D(curveKey, curveKey, tmp,
481 m_xmap->GetBasis(0)->GetPointsKey(),
482 m_xmap->GetBasis(1)->GetPointsKey(),
483 phys);
484
485 // Forwards transform to get coefficient space.
486 m_xmap->FwdTrans(phys, m_coeffs[i]);
487 }
488 }
489 else
490 {
491 ASSERTL0(false, "Only 1D/2D points distributions "
492 "supported.");
493 }
494 }
495
496 Array<OneD, unsigned int> mapArray(nEdgeCoeffs);
497 Array<OneD, int> signArray(nEdgeCoeffs);
498
499 for (i = 0; i < kNedges; i++)
500 {
501 m_edges[i]->FillGeom();
502 m_xmap->GetTraceToElementMap(i, mapArray, signArray, m_eorient[i]);
503
504 nEdgeCoeffs = m_edges[i]->GetXmap()->GetNcoeffs();
505
506 for (j = 0; j < m_coordim; j++)
507 {
508 for (k = 0; k < nEdgeCoeffs; k++)
509 {
510 m_coeffs[j][mapArray[k]] =
511 signArray[k] * m_edges[i]->GetCoeffs(j)[k];
512 }
513 }
514 }
515
517}
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:209
GeomState m_state
Enumeration to dictate whether coefficients are filled.
Definition: Geometry.h:194
Array< OneD, Array< OneD, NekDouble > > m_coeffs
Array containing expansion coefficients of m_xmap.
Definition: Geometry.h:204
PointsManagerT & PointsManager(void)
void Interp2D(const BasisKey &fbasis0, const BasisKey &fbasis1, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, const BasisKey &tbasis1, Array< OneD, NekDouble > &to)
this function interpolates a 2D function evaluated at the quadrature points of the 2D basis,...
Definition: Interp.cpp:103
@ eOrtho_A
Principle Orthogonal Functions .
Definition: BasisType.h:44
@ eOrtho_B
Principle Orthogonal Functions .
Definition: BasisType.h:46
static const NekDouble kVertexTheSameDouble
std::shared_ptr< Curve > CurveSharedPtr
Definition: Curve.hpp:60
@ ePtsFilled
Geometric information has been generated.
std::shared_ptr< StdNodalTriExp > StdNodalTriExpSharedPtr

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), ASSERTL0, Nektar::StdRegions::eBackwards, Nektar::StdRegions::eForwards, Nektar::LibUtilities::eGaussLobattoLegendre, Nektar::LibUtilities::eOrtho_A, Nektar::LibUtilities::eOrtho_B, Nektar::SpatialDomains::ePtsFilled, Nektar::ErrorUtil::ewarning, Nektar::SpatialDomains::Geometry::GetXmap(), Nektar::LibUtilities::Interp2D(), kNedges, Nektar::NekConstants::kVertexTheSameDouble, Nektar::SpatialDomains::Geometry::m_coeffs, Nektar::SpatialDomains::Geometry::m_coordim, Nektar::SpatialDomains::Geometry2D::m_curve, Nektar::SpatialDomains::Geometry2D::m_edges, Nektar::SpatialDomains::Geometry2D::m_eorient, Nektar::SpatialDomains::Geometry::m_globalID, Nektar::SpatialDomains::Geometry::m_state, Nektar::SpatialDomains::Geometry2D::m_verts, Nektar::SpatialDomains::Geometry::m_xmap, NEKERROR, Nektar::LibUtilities::PointsManager(), and tinysimd::sqrt().

Referenced by v_GenGeomFactors().

◆ v_GenGeomFactors()

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

Implements Nektar::SpatialDomains::Geometry.

Definition at line 227 of file TriGeom.cpp.

228{
229 if (!m_setupState)
230 {
232 }
233
235 {
236 GeomType Gtype = eRegular;
237
239
240 // check to see if expansions are linear
241 m_straightEdge = true;
242 if (m_xmap->GetBasisNumModes(0) != 2 ||
243 m_xmap->GetBasisNumModes(1) != 2)
244 {
245 Gtype = eDeformed;
246 m_straightEdge = false;
247 }
248
249 m_manifold = Array<OneD, int>(2);
250 m_manifold[0] = 0;
251 m_manifold[1] = 1;
252 if (m_coordim == 3)
253 {
254 PointGeom e01, e21, norm;
255 e01.Sub(*m_verts[0], *m_verts[1]);
256 e21.Sub(*m_verts[2], *m_verts[1]);
257 norm.Mult(e01, e21);
258 int tmpi = 0;
259 double tmp = std::fabs(norm[0]);
260 if (tmp < fabs(norm[1]))
261 {
262 tmp = fabs(norm[1]);
263 tmpi = 1;
264 }
265 if (tmp < fabs(norm[2]))
266 {
267 tmpi = 2;
268 }
269 m_manifold[0] = (tmpi + 1) % 3;
270 m_manifold[1] = (tmpi + 2) % 3;
271 }
272 if (Gtype == eRegular)
273 {
274 Array<OneD, Array<OneD, NekDouble>> verts(m_verts.size());
275 for (int i = 0; i < m_verts.size(); ++i)
276 {
277 verts[i] = Array<OneD, NekDouble>(3);
278 m_verts[i]->GetCoords(verts[i]);
279 }
280 // a00 + a01 xi1 + a02 xi2
281 // a10 + a11 xi1 + a12 xi2
282 m_isoParameter = Array<OneD, Array<OneD, NekDouble>>(2);
283 for (int i = 0; i < 2; ++i)
284 {
285 unsigned int d = m_manifold[i];
286 m_isoParameter[i] = Array<OneD, NekDouble>(3, 0.);
287 NekDouble A = verts[0][d];
288 NekDouble B = verts[1][d];
289 NekDouble C = verts[2][d];
290 m_isoParameter[i][0] = 0.5 * (B + C); // 1
291 m_isoParameter[i][1] = 0.5 * (-A + B); // xi1
292 m_isoParameter[i][2] = 0.5 * (-A + C); // xi2
293 }
295 }
296
298 Gtype, m_coordim, m_xmap, m_coeffs);
299
301 }
302}
Array< OneD, int > m_manifold
Definition: Geometry2D.h:86
virtual void v_CalculateInverseIsoParam() override
Definition: Geometry2D.cpp:625
bool m_setupState
Wether or not the setup routines have been run.
Definition: Geometry.h:196
Array< OneD, Array< OneD, NekDouble > > m_isoParameter
Definition: Geometry.h:207
GeomState m_geomFactorsState
State of the geometric factors.
Definition: Geometry.h:190
GeomFactorsSharedPtr m_geomFactors
Geometric factors.
Definition: Geometry.h:188
virtual void v_FillGeom() override
Definition: TriGeom.cpp:309
virtual void v_Setup() override
Definition: TriGeom.cpp:545
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.
std::vector< double > d(NPUPPER *NPUPPER)

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::UnitTests::d(), 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::Geometry2D::m_manifold, Nektar::SpatialDomains::Geometry::m_setupState, Nektar::SpatialDomains::Geometry::m_straightEdge, Nektar::SpatialDomains::Geometry2D::m_verts, Nektar::SpatialDomains::Geometry::m_xmap, Nektar::SpatialDomains::PointGeom::Mult(), Nektar::SpatialDomains::PointGeom::Sub(), Nektar::SpatialDomains::Geometry2D::v_CalculateInverseIsoParam(), v_FillGeom(), and v_Setup().

◆ v_GetCoord()

NekDouble Nektar::SpatialDomains::TriGeom::v_GetCoord ( const int  i,
const Array< OneD, const NekDouble > &  Lcoord 
)
overrideprotectedvirtual

Given local collapsed coordinate Lcoord, return the value of physical coordinate in direction i.

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 106 of file TriGeom.cpp.

108{
109 ASSERTL1(m_state == ePtsFilled, "Geometry is not in physical space");
110
111 Array<OneD, NekDouble> tmp(m_xmap->GetTotPoints());
112 m_xmap->BwdTrans(m_coeffs[i], tmp);
113
114 return m_xmap->PhysEvaluate(Lcoord, tmp);
115}
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:249

References ASSERTL1, Nektar::SpatialDomains::ePtsFilled, Nektar::SpatialDomains::Geometry::m_coeffs, Nektar::SpatialDomains::Geometry::m_state, and Nektar::SpatialDomains::Geometry::m_xmap.

◆ v_GetDir()

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

520{
521 boost::ignore_unused(j); // required in 3D shapes
522
523 return i == 0 ? 0 : 1;
524}

◆ v_Reset()

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

527{
528 Geometry::v_Reset(curvedEdges, curvedFaces);
529 CurveMap::iterator it = curvedFaces.find(m_globalID);
530
531 if (it != curvedFaces.end())
532 {
533 m_curve = it->second;
534 }
535
536 for (int i = 0; i < 3; ++i)
537 {
538 m_edges[i]->Reset(curvedEdges, curvedFaces);
539 }
540
541 SetUpXmap();
542 SetUpCoeffs(m_xmap->GetNcoeffs());
543}
void SetUpCoeffs(const int nCoeffs)
Initialise the Geometry::m_coeffs array.
Definition: Geometry.h:690
virtual void v_Reset(CurveMap &curvedEdges, CurveMap &curvedFaces)
Reset this geometry object: unset the current state, zero Geometry::m_coeffs and remove allocated Geo...
Definition: Geometry.cpp:390

References Nektar::SpatialDomains::Geometry2D::m_curve, Nektar::SpatialDomains::Geometry2D::m_edges, Nektar::SpatialDomains::Geometry::m_globalID, Nektar::SpatialDomains::Geometry::m_xmap, Nektar::SpatialDomains::Geometry::SetUpCoeffs(), SetUpXmap(), and Nektar::SpatialDomains::Geometry::v_Reset().

◆ v_Setup()

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

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 545 of file TriGeom.cpp.

546{
547 if (!m_setupState)
548 {
549 for (int i = 0; i < 3; ++i)
550 {
551 m_edges[i]->Setup();
552 }
553 SetUpXmap();
554 SetUpCoeffs(m_xmap->GetNcoeffs());
555 m_setupState = true;
556 }
557}

References Nektar::SpatialDomains::Geometry2D::m_edges, Nektar::SpatialDomains::Geometry::m_setupState, Nektar::SpatialDomains::Geometry::m_xmap, Nektar::SpatialDomains::Geometry::SetUpCoeffs(), and SetUpXmap().

Referenced by v_GenGeomFactors().

Member Data Documentation

◆ kNedges

const int Nektar::SpatialDomains::TriGeom::kNedges = 3
static

◆ kNverts

const int Nektar::SpatialDomains::TriGeom::kNverts = 3
static

Definition at line 73 of file TriGeom.h.

Referenced by Nektar::SpatialDomains::TetGeom::SetUpFaceOrientation().