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 () override
 
- Public Member Functions inherited from Nektar::SpatialDomains::Geometry2D
 Geometry2D ()
 
 Geometry2D (const int coordim, CurveSharedPtr curve)
 
 ~Geometry2D () override
 
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

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_GenGeomFactors () override
 
void v_FillGeom () override
 
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::Geometry2D
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 SolveStraightEdgeQuad (const Array< OneD, const NekDouble > &coords, Array< OneD, NekDouble > &Lcoords)
 
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...
 
PointGeomSharedPtr v_GetVertex (int i) const override
 
Geometry1DSharedPtr v_GetEdge (int i) const override
 Returns edge i of this object. 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...
 
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...
 
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
 
int m_straightEdge
 
- Static Protected Attributes inherited from Nektar::SpatialDomains::Geometry
static GeomFactorsVector m_regGeomFactorsManager
 

Detailed Description

Definition at line 59 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 83 of file TriGeom.cpp.

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

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 54 of file TriGeom.cpp.

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

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

Definition at line 100 of file TriGeom.cpp.

101{
102}

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 125 of file TriGeom.cpp.

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

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

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 555 of file TriGeom.cpp.

556{
557 int order0 = m_edges[0]->GetXmap()->GetBasis(0)->GetNumModes();
558 int order1 =
559 max(order0, max(m_edges[1]->GetXmap()->GetBasis(0)->GetNumModes(),
560 m_edges[2]->GetXmap()->GetBasis(0)->GetNumModes()));
561
562 const LibUtilities::BasisKey B0(
564 LibUtilities::PointsKey(order0 + 1,
566 const LibUtilities::BasisKey B1(
568 LibUtilities::PointsKey(order1, LibUtilities::eGaussRadauMAlpha1Beta0));
569
571}
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:189
StdRegions::StdExpansionSharedPtr GetXmap() const
Return the mapping object Geometry::m_xmap that represents the coordinate transformation from standar...
Definition: Geometry.h:433
@ eGaussLobattoLegendre
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:51
@ eModified_B
Principle Modified Functions .
Definition: BasisType.h:49
@ eModified_A
Principle Modified Functions .
Definition: BasisType.h:48

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 307 of file TriGeom.cpp.

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

226{
227 if (!m_setupState)
228 {
230 }
231
233 {
234 GeomType Gtype = eRegular;
235
237
238 // check to see if expansions are linear
239 m_straightEdge = 1;
240 if (m_xmap->GetBasisNumModes(0) != 2 ||
241 m_xmap->GetBasisNumModes(1) != 2)
242 {
243 Gtype = eDeformed;
244 m_straightEdge = 0;
245 }
246
247 m_manifold = Array<OneD, int>(2);
248 m_manifold[0] = 0;
249 m_manifold[1] = 1;
250 if (m_coordim == 3)
251 {
252 PointGeom e01, e21, norm;
253 e01.Sub(*m_verts[0], *m_verts[1]);
254 e21.Sub(*m_verts[2], *m_verts[1]);
255 norm.Mult(e01, e21);
256 int tmpi = 0;
257 double tmp = std::fabs(norm[0]);
258 if (tmp < fabs(norm[1]))
259 {
260 tmp = fabs(norm[1]);
261 tmpi = 1;
262 }
263 if (tmp < fabs(norm[2]))
264 {
265 tmpi = 2;
266 }
267 m_manifold[0] = (tmpi + 1) % 3;
268 m_manifold[1] = (tmpi + 2) % 3;
269 }
270 if (Gtype == eRegular)
271 {
272 Array<OneD, Array<OneD, NekDouble>> verts(m_verts.size());
273 for (int i = 0; i < m_verts.size(); ++i)
274 {
275 verts[i] = Array<OneD, NekDouble>(3);
276 m_verts[i]->GetCoords(verts[i]);
277 }
278 // a00 + a01 xi1 + a02 xi2
279 // a10 + a11 xi1 + a12 xi2
280 m_isoParameter = Array<OneD, Array<OneD, NekDouble>>(2);
281 for (int i = 0; i < 2; ++i)
282 {
283 unsigned int d = m_manifold[i];
284 m_isoParameter[i] = Array<OneD, NekDouble>(3, 0.);
285 NekDouble A = verts[0][d];
286 NekDouble B = verts[1][d];
287 NekDouble C = verts[2][d];
288 m_isoParameter[i][0] = 0.5 * (B + C); // 1
289 m_isoParameter[i][1] = 0.5 * (-A + B); // xi1
290 m_isoParameter[i][2] = 0.5 * (-A + C); // xi2
291 }
293 }
294
296 Gtype, m_coordim, m_xmap, m_coeffs);
297
299 }
300}
Array< OneD, int > m_manifold
Definition: Geometry2D.h:83
void v_CalculateInverseIsoParam() override
Definition: Geometry2D.cpp:575
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
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 104 of file TriGeom.cpp.

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

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 517 of file TriGeom.cpp.

518{
519 return i == 0 ? 0 : 1;
520}

◆ 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 522 of file TriGeom.cpp.

523{
524 Geometry::v_Reset(curvedEdges, curvedFaces);
525 CurveMap::iterator it = curvedFaces.find(m_globalID);
526
527 if (it != curvedFaces.end())
528 {
529 m_curve = it->second;
530 }
531
532 for (int i = 0; i < 3; ++i)
533 {
534 m_edges[i]->Reset(curvedEdges, curvedFaces);
535 }
536
537 SetUpXmap();
538 SetUpCoeffs(m_xmap->GetNcoeffs());
539}
void SetUpCoeffs(const int nCoeffs)
Initialise the Geometry::m_coeffs array.
Definition: Geometry.h:687
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:384

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 541 of file TriGeom.cpp.

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

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 71 of file TriGeom.h.

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