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

#include <SegGeom.h>

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

Public Member Functions

 SegGeom ()
 
 SegGeom (int id, const int coordim, const PointGeomSharedPtr vertex[], const CurveSharedPtr curve=CurveSharedPtr())
 
 SegGeom (const SegGeom &in)
 
SegGeomSharedPtr GenerateOneSpaceDimGeom (void)
 Generate a one dimensional space segment geometry where the vert[0] has the same x value and vert[1] is set to vert[0] plus the length of the original segment. More...
 
 ~SegGeom ()
 
CurveSharedPtr GetCurve ()
 
- Public Member Functions inherited from Nektar::SpatialDomains::Geometry1D
 Geometry1D ()
 
 Geometry1D (const int coordim)
 
virtual ~Geometry1D ()
 
- 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 GetEdgeOrientation (const SegGeom &edge1, const SegGeom &edge2)
 Get the orientation of edge1. More...
 

Static Public Attributes

static const int kNverts = 2
 
- Static Public Attributes inherited from Nektar::SpatialDomains::Geometry1D
static const int kDim = 1
 

Protected Member Functions

virtual PointGeomSharedPtr v_GetVertex (const int i) const override
 
virtual LibUtilities::ShapeType v_GetShapeType () const
 
virtual void v_GenGeomFactors () override
 
virtual void v_FillGeom () override
 Populate the coordinate mapping Geometry::m_coeffs information from any children geometry elements. 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
 
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 int v_GetNumVerts () const override
 Get the number of vertices of this object. More...
 
virtual NekDouble v_FindDistance (const Array< OneD, const NekDouble > &xs, Array< OneD, NekDouble > &xi) override
 
- Protected Member Functions inherited from Nektar::SpatialDomains::Geometry1D
virtual int v_GetShapeDim () const override
 Get the object's shape dimension. More...
 
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...
 
- 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 ()
 

Protected Attributes

SpatialDomains::PointGeomSharedPtr m_verts [kNverts]
 
StdRegions::Orientation m_porient [kNverts]
 
- 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
 

Private Member Functions

void SetUpXmap ()
 

Private Attributes

CurveSharedPtr m_curve
 Boolean indicating whether object owns the data. More...
 

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

Detailed Description

Definition at line 54 of file SegGeom.h.

Constructor & Destructor Documentation

◆ SegGeom() [1/3]

Nektar::SpatialDomains::SegGeom::SegGeom ( )

Definition at line 47 of file SegGeom.cpp.

48{
50}
LibUtilities::ShapeType m_shapeType
Type of shape.
Definition: Geometry.h:200

References Nektar::LibUtilities::eSegment, and Nektar::SpatialDomains::Geometry::m_shapeType.

◆ SegGeom() [2/3]

Nektar::SpatialDomains::SegGeom::SegGeom ( int  id,
const int  coordim,
const PointGeomSharedPtr  vertex[],
const CurveSharedPtr  curve = CurveSharedPtr() 
)

Definition at line 52 of file SegGeom.cpp.

54 : Geometry1D(coordim)
55{
57 m_globalID = id;
59 m_curve = curve;
60
61 m_verts[0] = vertex[0];
62 m_verts[1] = vertex[1];
63}
GeomState m_state
Enumeration to dictate whether coefficients are filled.
Definition: Geometry.h:194
CurveSharedPtr m_curve
Boolean indicating whether object owns the data.
Definition: SegGeom.h:96
SpatialDomains::PointGeomSharedPtr m_verts[kNverts]
Definition: SegGeom.h:79
@ eNotFilled
Geometric information has not been generated.

References Nektar::SpatialDomains::eNotFilled, Nektar::LibUtilities::eSegment, m_curve, Nektar::SpatialDomains::Geometry::m_globalID, Nektar::SpatialDomains::Geometry::m_shapeType, Nektar::SpatialDomains::Geometry::m_state, and m_verts.

◆ SegGeom() [3/3]

Nektar::SpatialDomains::SegGeom::SegGeom ( const SegGeom in)

Definition at line 65 of file SegGeom.cpp.

65 : Geometry1D(in)
66{
67 // From Geometry class
68 m_shapeType = in.m_shapeType;
69
70 // info from EdgeComponent class
71 m_globalID = in.m_globalID;
72 m_xmap = in.m_xmap;
73 SetUpCoeffs(m_xmap->GetNcoeffs());
74
75 // info from SegGeom class
76 m_coordim = in.m_coordim;
77 m_verts[0] = in.m_verts[0];
78 m_verts[1] = in.m_verts[1];
79
80 m_state = in.m_state;
81}
void SetUpCoeffs(const int nCoeffs)
Initialise the Geometry::m_coeffs array.
Definition: Geometry.h:690
StdRegions::StdExpansionSharedPtr m_xmap
mapping containing isoparametric transformation.
Definition: Geometry.h:192
int m_coordim
Coordinate dimension of this geometry object.
Definition: Geometry.h:186

References Nektar::SpatialDomains::Geometry::m_coordim, Nektar::SpatialDomains::Geometry::m_globalID, Nektar::SpatialDomains::Geometry::m_shapeType, Nektar::SpatialDomains::Geometry::m_state, m_verts, Nektar::SpatialDomains::Geometry::m_xmap, and Nektar::SpatialDomains::Geometry::SetUpCoeffs().

◆ ~SegGeom()

Nektar::SpatialDomains::SegGeom::~SegGeom ( )

Definition at line 161 of file SegGeom.cpp.

162{
163}

Member Function Documentation

◆ GenerateOneSpaceDimGeom()

SegGeomSharedPtr Nektar::SpatialDomains::SegGeom::GenerateOneSpaceDimGeom ( void  )

Generate a one dimensional space segment geometry where the vert[0] has the same x value and vert[1] is set to vert[0] plus the length of the original segment.

Definition at line 107 of file SegGeom.cpp.

108{
110
111 // info about numbering
112 returnval->m_globalID = m_globalID;
113
114 // geometric information.
115 returnval->m_coordim = 1;
116 NekDouble x0 = (*m_verts[0])[0];
118 1, m_verts[0]->GetGlobalID(), x0, 0.0, 0.0);
119 vert0->SetGlobalID(vert0->GetGlobalID());
120 returnval->m_verts[0] = vert0;
121
122 // Get information to calculate length.
123 const Array<OneD, const LibUtilities::BasisSharedPtr> base =
124 m_xmap->GetBase();
126 v.push_back(base[0]->GetPointsKey());
128
129 const Array<OneD, const NekDouble> jac = m_geomFactors->GetJac(v);
130
131 NekDouble len = 0.0;
132 if (jac.size() == 1)
133 {
134 len = jac[0] * 2.0;
135 }
136 else
137 {
138 Array<OneD, const NekDouble> w0 = base[0]->GetW();
139 len = 0.0;
140
141 for (int i = 0; i < jac.size(); ++i)
142 {
143 len += jac[i] * w0[i];
144 }
145 }
146 // Set up second vertex.
148 1, m_verts[1]->GetGlobalID(), x0 + len, 0.0, 0.0);
149 vert1->SetGlobalID(vert1->GetGlobalID());
150
151 returnval->m_verts[1] = vert1;
152
153 // at present just use previous m_xmap[0];
154 returnval->m_xmap = m_xmap;
155 returnval->SetUpCoeffs(m_xmap->GetNcoeffs());
156 returnval->m_state = eNotFilled;
157
158 return returnval;
159}
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
int GetGlobalID(void) const
Get the ID of this object.
Definition: Geometry.h:327
GeomFactorsSharedPtr m_geomFactors
Geometric factors.
Definition: Geometry.h:188
virtual void v_GenGeomFactors() override
Definition: SegGeom.cpp:220
std::vector< PointsKey > PointsKeyVector
Definition: Points.h:236
std::shared_ptr< SegGeom > SegGeomSharedPtr
Definition: Geometry2D.h:62
std::shared_ptr< PointGeom > PointGeomSharedPtr
Definition: Geometry.h:60
double NekDouble

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::SpatialDomains::eNotFilled, Nektar::SpatialDomains::Geometry::GetGlobalID(), Nektar::SpatialDomains::Geometry::m_geomFactors, Nektar::SpatialDomains::Geometry::m_globalID, m_verts, Nektar::SpatialDomains::Geometry::m_xmap, and v_GenGeomFactors().

◆ GetCurve()

CurveSharedPtr Nektar::SpatialDomains::SegGeom::GetCurve ( )
inline

Definition at line 71 of file SegGeom.h.

72 {
73 return m_curve;
74 }

References m_curve.

◆ GetEdgeOrientation()

StdRegions::Orientation Nektar::SpatialDomains::SegGeom::GetEdgeOrientation ( const SegGeom edge1,
const SegGeom edge2 
)
static

Get the orientation of edge1.

If edge1 is connected to edge2 in the same direction as the points comprising edge1 then it is forward, otherwise it is backward.

For example, assume edge1 is comprised of points 1 and 2, and edge2 is comprised of points 2 and 3, then edge1 is forward.

If edge1 is comprised of points 2 and 1 and edge2 is comprised of points 3 and 2, then edge1 is backward.

Since both edges are passed, it does not need any information from the EdgeComponent instance.

Definition at line 196 of file SegGeom.cpp.

198{
200
201 if ((*edge1.GetVertex(0) == *edge2.GetVertex(0)) ||
202 (*edge1.GetVertex(0) == *edge2.GetVertex(1)))
203 {
204 // Backward direction. Vertex 0 is connected to edge 2.
205 returnval = StdRegions::eBackwards;
206 }
207 else if ((*edge1.GetVertex(1) != *edge2.GetVertex(0)) &&
208 (*edge1.GetVertex(1) != *edge2.GetVertex(1)))
209 {
210 // Not forward either, then we have a problem.
211 std::ostringstream errstrm;
212 errstrm << "Connected edges do not share a vertex. Edges ";
213 errstrm << edge1.GetGlobalID() << ", " << edge2.GetGlobalID();
214 ASSERTL0(false, errstrm.str());
215 }
216
217 return returnval;
218}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215

References ASSERTL0, Nektar::StdRegions::eBackwards, Nektar::StdRegions::eForwards, Nektar::SpatialDomains::Geometry::GetGlobalID(), and Nektar::SpatialDomains::Geometry::GetVertex().

Referenced by Nektar::SpatialDomains::QuadGeom::QuadGeom(), and Nektar::SpatialDomains::TriGeom::TriGeom().

◆ SetUpXmap()

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

Definition at line 83 of file SegGeom.cpp.

84{
85 if (m_curve)
86 {
87 int npts = m_curve->m_points.size();
88 LibUtilities::PointsKey pkey(npts + 1,
90 const LibUtilities::BasisKey B(LibUtilities::eModified_A, npts, pkey);
92 }
93 else
94 {
95 const LibUtilities::BasisKey B(
97 LibUtilities::PointsKey(2, LibUtilities::eGaussLobattoLegendre));
99 }
100}
@ eGaussLobattoLegendre
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:53
@ eModified_A
Principle Modified Functions .
Definition: BasisType.h:50

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::LibUtilities::eGaussLobattoLegendre, Nektar::LibUtilities::eModified_A, m_curve, and Nektar::SpatialDomains::Geometry::m_xmap.

Referenced by v_Reset(), and v_Setup().

◆ v_FillGeom()

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

Populate the coordinate mapping Geometry::m_coeffs information from any children geometry elements.

See also
v_FillGeom()

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 243 of file SegGeom.cpp.

244{
245 if (m_state != ePtsFilled)
246 {
247 int i;
248
249 if (m_coordim > 0 && m_curve)
250 {
251 int npts = m_curve->m_points.size();
252 LibUtilities::PointsKey pkey(npts + 1,
254 Array<OneD, NekDouble> tmp(npts);
255
256 if (m_verts[0]->dist(*(m_curve->m_points[0])) >
258 {
259 std::string err =
260 "Vertex 0 is separated from first point by more than ";
261 std::stringstream strstrm;
262 strstrm << NekConstants::kVertexTheSameDouble << " in edge "
263 << m_globalID;
264 err += strstrm.str();
265 NEKERROR(ErrorUtil::ewarning, err.c_str());
266 }
267
268 if (m_verts[1]->dist(*(m_curve->m_points[npts - 1])) >
270 {
271 std::string err =
272 "Vertex 1 is separated from last point by more than ";
273 std::stringstream strstrm;
274 strstrm << NekConstants::kVertexTheSameDouble << " in edge "
275 << m_globalID;
276 err += strstrm.str();
277 NEKERROR(ErrorUtil::ewarning, err.c_str());
278 }
279
280 LibUtilities::PointsKey fkey(npts, m_curve->m_ptype);
282 LibUtilities::PointsManager()[fkey]->GetI(pkey);
283 NekVector<NekDouble> out(npts + 1);
284
285 for (int i = 0; i < m_coordim; ++i)
286 {
287 // Load up coordinate values into tmp
288 for (int j = 0; j < npts; ++j)
289 {
290 tmp[j] = (m_curve->m_points[j]->GetPtr())[i];
291 }
292
293 // Interpolate to GLL points
294 NekVector<NekDouble> in(npts, tmp, eWrapper);
295 out = (*I0) * in;
296
297 m_xmap->FwdTrans(out.GetPtr(), m_coeffs[i]);
298 }
299 }
300
301 for (i = 0; i < m_coordim; ++i)
302 {
303 m_coeffs[i][0] = (*m_verts[0])[i];
304 m_coeffs[i][1] = (*m_verts[1])[i];
305 }
306
308 }
309}
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:209
Array< OneD, Array< OneD, NekDouble > > m_coeffs
Array containing expansion coefficients of m_xmap.
Definition: Geometry.h:204
PointsManagerT & PointsManager(void)
static const NekDouble kVertexTheSameDouble
@ ePtsFilled
Geometric information has been generated.
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:75

References Nektar::LibUtilities::eGaussLobattoLegendre, Nektar::SpatialDomains::ePtsFilled, Nektar::ErrorUtil::ewarning, Nektar::eWrapper, Nektar::NekVector< DataType >::GetPtr(), Nektar::NekConstants::kVertexTheSameDouble, Nektar::SpatialDomains::Geometry::m_coeffs, Nektar::SpatialDomains::Geometry::m_coordim, m_curve, Nektar::SpatialDomains::Geometry::m_globalID, Nektar::SpatialDomains::Geometry::m_state, m_verts, Nektar::SpatialDomains::Geometry::m_xmap, NEKERROR, and Nektar::LibUtilities::PointsManager().

Referenced by v_GenGeomFactors().

◆ v_FindDistance()

NekDouble Nektar::SpatialDomains::SegGeom::v_FindDistance ( const Array< OneD, const NekDouble > &  xs,
Array< OneD, NekDouble > &  xi 
)
overrideprotectedvirtual

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 352 of file SegGeom.cpp.

354{
355 if (m_geomFactors->GetGtype() == eRegular)
356 {
357 xiOut = Array<OneD, NekDouble>(1, 0.0);
358
359 GetLocCoords(xs, xiOut);
360 ClampLocCoords(xiOut);
361
362 Array<OneD, NekDouble> gloCoord(m_coordim);
363 NekDouble tmp = 0;
364 for (int i = 0; i < m_coordim; ++i)
365 {
366 gloCoord[i] = GetCoord(i, xiOut);
367 tmp += (xs[i] - gloCoord[i]) * (xs[i] - gloCoord[i]);
368 }
369
370 return sqrt(tmp);
371 }
372 // If deformed edge then the inverse mapping is non-linear so need to
373 // numerically solve for the local coordinate
374 else if (m_geomFactors->GetGtype() == eDeformed)
375 {
376 Array<OneD, NekDouble> xi(1, 0.0);
377
378 // Armijo constants:
379 // https://en.wikipedia.org/wiki/Backtracking_line_search
380 const NekDouble c1 = 1e-4, c2 = 0.9;
381
382 int dim = GetCoordim();
383 int nq = m_xmap->GetTotPoints();
384
385 Array<OneD, Array<OneD, NekDouble>> x(dim), xder(dim), xder2(dim);
386 // Get x,y,z phys values from coefficients
387 for (int i = 0; i < dim; ++i)
388 {
389 x[i] = Array<OneD, NekDouble>(nq);
390 xder[i] = Array<OneD, NekDouble>(nq);
391 xder2[i] = Array<OneD, NekDouble>(nq);
392
393 m_xmap->BwdTrans(m_coeffs[i], x[i]);
394 }
395
396 NekDouble fx_prev = std::numeric_limits<NekDouble>::max();
397
398 // Minimisation loop (Quasi-newton method)
399 for (int i = 0; i < NekConstants::kNewtonIterations; ++i)
400 {
401 // Compute the objective function, f(x_k) and its derivatives
402 Array<OneD, NekDouble> xc(dim);
403 Array<OneD, std::array<NekDouble, 3>> xc_der(dim);
404 Array<OneD, std::array<NekDouble, 6>> xc_der2(dim);
405 NekDouble fx = 0, fxp = 0, fxp2 = 0, xcDiff = 0;
406 for (int j = 0; j < dim; ++j)
407 {
408 xc[j] = m_xmap->PhysEvaluate(xi, x[j], xc_der[j], xc_der2[j]);
409
410 xcDiff = xc[j] - xs[j];
411 // Objective function is the distance to the search point
412 fx += xcDiff * xcDiff;
413 fxp += xc_der[j][0] * xcDiff;
414 fxp2 += xc_der2[j][0] * xcDiff + xc_der[j][0] * xc_der[j][0];
415 }
416
417 fxp *= 2;
418 fxp2 *= 2;
419
420 // Check for convergence
421 if (std::abs(fx - fx_prev) < 1e-12)
422 {
423 fx_prev = fx;
424 break;
425 }
426 else
427 {
428 fx_prev = fx;
429 }
430
431 NekDouble gamma = 1.0;
432 bool conv = false;
433
434 // Search direction: Newton's method
435 NekDouble pk = -fxp / fxp2;
436
437 // Perform backtracking line search
438 while (gamma > 1e-10)
439 {
440 Array<OneD, NekDouble> xi_pk(1);
441 xi_pk[0] = xi[0] + pk * gamma;
442
443 if (xi_pk[0] < -1.0 || xi_pk[0] > 1.0)
444 {
445 gamma /= 2.0;
446 continue;
447 }
448
449 Array<OneD, NekDouble> xc_pk(dim);
450 Array<OneD, std::array<NekDouble, 3>> xc_der_pk(dim);
451 NekDouble fx_pk = 0, fxp_pk = 0, xc_pkDiff = 0;
452 for (int j = 0; j < dim; ++j)
453 {
454 xc_pk[j] = m_xmap->PhysEvaluate(xi_pk, x[j], xc_der_pk[j]);
455
456 xc_pkDiff = xc_pk[j] - xs[j];
457 fx_pk += xc_pkDiff * xc_pkDiff;
458 fxp_pk += xc_der_pk[j][0] * xc_pkDiff;
459 }
460
461 fxp_pk *= 2;
462
463 // Check Wolfe conditions using Armijo constants
464 // https://en.wikipedia.org/wiki/Wolfe_conditions
465 if ((fx_pk - (fx + c1 * gamma * pk * fxp)) <
466 std::numeric_limits<NekDouble>::epsilon() &&
467 (-pk * fxp_pk + c2 * pk * fxp) <
468 std::numeric_limits<NekDouble>::epsilon())
469 {
470 conv = true;
471 break;
472 }
473
474 gamma /= 2.0;
475 }
476
477 if (!conv)
478 {
479 break;
480 }
481
482 xi[0] += gamma * pk;
483 }
484
485 xiOut = xi;
486 return sqrt(fx_prev);
487 }
488 else
489 {
490 ASSERTL0(false, "Geometry type unknown")
491 }
492
493 return -1.0;
494}
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.
Definition: Geometry.h:555
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 ge...
Definition: Geometry.h:545
bool ClampLocCoords(Array< OneD, NekDouble > &locCoord, NekDouble tol=std::numeric_limits< NekDouble >::epsilon())
Clamp local coords to be within standard regions [-1, 1]^dim.
Definition: Geometry.cpp:556
int GetCoordim() const
Return the coordinate dimension of this object (i.e. the dimension of the space in which this object ...
Definition: Geometry.h:284
static const unsigned int kNewtonIterations
@ eRegular
Geometry is straight-sided with constant geometric factors.
@ eDeformed
Geometry is curved or has non-constant factors.
scalarT< T > abs(scalarT< T > in)
Definition: scalar.hpp:298
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294

References tinysimd::abs(), ASSERTL0, Nektar::SpatialDomains::Geometry::ClampLocCoords(), Nektar::SpatialDomains::eDeformed, Nektar::SpatialDomains::eRegular, Nektar::SpatialDomains::Geometry::GetCoord(), Nektar::SpatialDomains::Geometry::GetCoordim(), Nektar::SpatialDomains::Geometry::GetLocCoords(), Nektar::NekConstants::kNewtonIterations, Nektar::SpatialDomains::Geometry::m_coeffs, Nektar::SpatialDomains::Geometry::m_coordim, Nektar::SpatialDomains::Geometry::m_geomFactors, Nektar::SpatialDomains::Geometry::m_xmap, and tinysimd::sqrt().

◆ v_GenGeomFactors()

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

Implements Nektar::SpatialDomains::Geometry.

Definition at line 220 of file SegGeom.cpp.

221{
222 if (!m_setupState)
223 {
225 }
226
228 {
231
232 if (m_xmap->GetBasisNumModes(0) != 2)
233 {
234 gType = eDeformed;
235 }
236
238 gType, m_coordim, m_xmap, m_coeffs);
240 }
241}
bool m_setupState
Wether or not the setup routines have been run.
Definition: Geometry.h:196
GeomState m_geomFactorsState
State of the geometric factors.
Definition: Geometry.h:190
virtual void v_Setup() override
Definition: SegGeom.cpp:325
virtual void v_FillGeom() override
Populate the coordinate mapping Geometry::m_coeffs information from any children geometry elements.
Definition: SegGeom.cpp:243
GeomType
Indicates the type of element geometry.

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_setupState, Nektar::SpatialDomains::Geometry::m_xmap, v_FillGeom(), and v_Setup().

Referenced by GenerateOneSpaceDimGeom().

◆ v_GetCoord()

NekDouble Nektar::SpatialDomains::SegGeom::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 170 of file SegGeom.cpp.

172{
173 ASSERTL1(m_state == ePtsFilled, "Geometry is not in physical space");
174
175 Array<OneD, NekDouble> tmp(m_xmap->GetTotPoints());
176 m_xmap->BwdTrans(m_coeffs[i], tmp);
177
178 return m_xmap->PhysEvaluate(Lcoord, tmp);
179}
#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_GetNumVerts()

int Nektar::SpatialDomains::SegGeom::v_GetNumVerts ( ) const
overrideprotectedvirtual

Get the number of vertices of this object.

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 347 of file SegGeom.cpp.

348{
349 return kNverts;
350}
static const int kNverts
Definition: SegGeom.h:76

References kNverts.

◆ v_GetShapeType()

LibUtilities::ShapeType Nektar::SpatialDomains::SegGeom::v_GetShapeType ( ) const
protectedvirtual

Definition at line 165 of file SegGeom.cpp.

166{
168}

References Nektar::LibUtilities::eSegment.

◆ v_GetVertex()

PointGeomSharedPtr Nektar::SpatialDomains::SegGeom::v_GetVertex ( const int  i) const
overrideprotectedvirtual

Implements Nektar::SpatialDomains::Geometry.

Definition at line 335 of file SegGeom.cpp.

336{
337 PointGeomSharedPtr returnval;
338
339 if (i >= 0 && i < kNverts)
340 {
341 returnval = m_verts[i];
342 }
343
344 return returnval;
345}

References kNverts, and m_verts.

◆ v_Reset()

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

312{
313 Geometry::v_Reset(curvedEdges, curvedFaces);
314 CurveMap::iterator it = curvedEdges.find(m_globalID);
315
316 if (it != curvedEdges.end())
317 {
318 m_curve = it->second;
319 }
320
321 SetUpXmap();
322 SetUpCoeffs(m_xmap->GetNcoeffs());
323}
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 m_curve, 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::SegGeom::v_Setup ( )
overrideprotectedvirtual

Member Data Documentation

◆ kNverts

const int Nektar::SpatialDomains::SegGeom::kNverts = 2
static

Definition at line 76 of file SegGeom.h.

Referenced by v_GetNumVerts(), and v_GetVertex().

◆ m_curve

CurveSharedPtr Nektar::SpatialDomains::SegGeom::m_curve
private

Boolean indicating whether object owns the data.

Definition at line 96 of file SegGeom.h.

Referenced by GetCurve(), SegGeom(), SetUpXmap(), v_FillGeom(), and v_Reset().

◆ m_porient

StdRegions::Orientation Nektar::SpatialDomains::SegGeom::m_porient[kNverts]
protected

Definition at line 80 of file SegGeom.h.

◆ m_verts

SpatialDomains::PointGeomSharedPtr Nektar::SpatialDomains::SegGeom::m_verts[kNverts]
protected

Definition at line 79 of file SegGeom.h.

Referenced by GenerateOneSpaceDimGeom(), SegGeom(), v_FillGeom(), and v_GetVertex().