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

Static Public 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

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

Constructor & Destructor Documentation

◆ SegGeom() [1/3]

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

Definition at line 45 of file SegGeom.cpp.

46{
48}
LibUtilities::ShapeType m_shapeType
Type of shape.
Definition: Geometry.h:203

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 50 of file SegGeom.cpp.

52 : Geometry1D(coordim)
53{
55 m_globalID = id;
57 m_curve = curve;
58
59 m_verts[0] = vertex[0];
60 m_verts[1] = vertex[1];
61}
GeomState m_state
Enumeration to dictate whether coefficients are filled.
Definition: Geometry.h:197
CurveSharedPtr m_curve
Boolean indicating whether object owns the data.
Definition: SegGeom.h:93
SpatialDomains::PointGeomSharedPtr m_verts[kNverts]
Definition: SegGeom.h:76
@ 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 63 of file SegGeom.cpp.

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

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

Definition at line 159 of file SegGeom.cpp.

160{
161}

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 105 of file SegGeom.cpp.

106{
108
109 // info about numbering
110 returnval->m_globalID = m_globalID;
111
112 // geometric information.
113 returnval->m_coordim = 1;
114 NekDouble x0 = (*m_verts[0])[0];
116 1, m_verts[0]->GetGlobalID(), x0, 0.0, 0.0);
117 vert0->SetGlobalID(vert0->GetGlobalID());
118 returnval->m_verts[0] = vert0;
119
120 // Get information to calculate length.
121 const Array<OneD, const LibUtilities::BasisSharedPtr> base =
122 m_xmap->GetBase();
124 v.push_back(base[0]->GetPointsKey());
126
127 const Array<OneD, const NekDouble> jac = m_geomFactors->GetJac(v);
128
129 NekDouble len = 0.0;
130 if (jac.size() == 1)
131 {
132 len = jac[0] * 2.0;
133 }
134 else
135 {
136 Array<OneD, const NekDouble> w0 = base[0]->GetW();
137 len = 0.0;
138
139 for (int i = 0; i < jac.size(); ++i)
140 {
141 len += jac[i] * w0[i];
142 }
143 }
144 // Set up second vertex.
146 1, m_verts[1]->GetGlobalID(), x0 + len, 0.0, 0.0);
147 vert1->SetGlobalID(vert1->GetGlobalID());
148
149 returnval->m_verts[1] = vert1;
150
151 // at present just use previous m_xmap[0];
152 returnval->m_xmap = m_xmap;
153 returnval->SetUpCoeffs(m_xmap->GetNcoeffs());
154 returnval->m_state = eNotFilled;
155
156 return returnval;
157}
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:191
void v_GenGeomFactors() override
Definition: SegGeom.cpp:218
std::vector< PointsKey > PointsKeyVector
Definition: Points.h:231
std::shared_ptr< SegGeom > SegGeomSharedPtr
Definition: Geometry2D.h:59
std::shared_ptr< PointGeom > PointGeomSharedPtr
Definition: Geometry.h:57
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 68 of file SegGeom.h.

69 {
70 return m_curve;
71 }

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 194 of file SegGeom.cpp.

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

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 81 of file SegGeom.cpp.

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

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::LibUtilities::eGaussLobattoLegendre, Nektar::LibUtilities::eModified_A, 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 241 of file SegGeom.cpp.

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

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

219{
220 if (!m_setupState)
221 {
223 }
224
226 {
229
230 if (m_xmap->GetBasisNumModes(0) != 2)
231 {
232 gType = eDeformed;
233 }
234
236 gType, m_coordim, m_xmap, m_coeffs);
238 }
239}
bool m_setupState
Wether or not the setup routines have been run.
Definition: Geometry.h:199
GeomState m_geomFactorsState
State of the geometric factors.
Definition: Geometry.h:193
void v_FillGeom() override
Populate the coordinate mapping Geometry::m_coeffs information from any children geometry elements.
Definition: SegGeom.cpp:241
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 168 of file SegGeom.cpp.

170{
171 ASSERTL1(m_state == ePtsFilled, "Geometry is not in physical space");
172
173 Array<OneD, NekDouble> tmp(m_xmap->GetTotPoints());
174 m_xmap->BwdTrans(m_coeffs[i], tmp);
175
176 return m_xmap->PhysEvaluate(Lcoord, tmp);
177}
#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_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 345 of file SegGeom.cpp.

346{
347 return kNverts;
348}
static const int kNverts
Definition: SegGeom.h:73

References kNverts.

◆ v_GetShapeType()

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

Definition at line 163 of file SegGeom.cpp.

164{
166}

References Nektar::LibUtilities::eSegment.

◆ v_GetVertex()

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

Implements Nektar::SpatialDomains::Geometry.

Definition at line 333 of file SegGeom.cpp.

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

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 309 of file SegGeom.cpp.

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

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 73 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 93 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 77 of file SegGeom.h.

◆ m_verts

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

Definition at line 76 of file SegGeom.h.

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