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

Static Public Member Functions

static StdRegions::Orientation GetEdgeOrientation (const SegGeom &edge1, const SegGeom &edge2)
 Get the orientation of edge1.
 

Static Public Attributes

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

Protected Member Functions

PointGeomv_GetVertex (const int i) const override
 Returns vertex i of this object.
 
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.
 
void v_Reset (CurveMap &curvedEdges, CurveMap &curvedFaces) override
 Reset this geometry object: unset the current state, zero Geometry::m_coeffs and remove allocated GeomFactors.
 
void v_Setup () override
 
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.
 
int v_GetNumVerts () const override
 Get the number of vertices of this object.
 
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.
 
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.
 
- Protected Member Functions inherited from Nektar::SpatialDomains::Geometry
virtual int v_GetVid (int i) const
 Get the ID of vertex i of this object.
 
virtual Geometry1Dv_GetEdge (int i) const
 Returns edge i of this object.
 
virtual Geometry2Dv_GetFace (int i) const
 Returns face i of this object.
 
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.
 
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.
 
virtual int v_GetNumEdges () const
 Get the number of edges of this object.
 
virtual int v_GetNumFaces () const
 Get the number of faces of this object.
 
virtual StdRegions::StdExpansionSharedPtr v_GetXmap () const
 Return the mapping object Geometry::m_xmap that represents the coordinate transformation from standard element to physical element.
 
virtual bool v_ContainsPoint (const Array< OneD, const NekDouble > &gloCoord, Array< OneD, NekDouble > &locCoord, NekDouble tol, NekDouble &dist)
 Determine whether an element contains a particular Cartesian coordinate \(\vec{x} = (x,y,z)\).
 
virtual int v_AllLeftCheck (const Array< OneD, const NekDouble > &gloCoord)
 
virtual int v_GetVertexEdgeMap (int i, int j) const
 Returns the standard element edge IDs that are connected to a given vertex.
 
virtual int v_GetVertexFaceMap (int i, int j) const
 Returns the standard element face IDs that are connected to a given vertex.
 
virtual int v_GetEdgeFaceMap (int i, int j) const
 Returns the standard element edge IDs that are connected to a given face.
 
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.
 
virtual int v_GetDir (const int faceidx, const int facedir) const
 Returns the element coordinate direction corresponding to a given face coordinate direction.
 
void SetUpCoeffs (const int nCoeffs)
 Initialise the Geometry::m_coeffs array.
 
virtual void v_CalculateInverseIsoParam ()
 

Protected Attributes

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

Private Member Functions

void SetUpXmap ()
 

Private Attributes

CurveSharedPtr m_curve
 Boolean indicating whether object owns the data.
 

Additional Inherited Members

- Static Protected Member Functions inherited from Nektar::SpatialDomains::Geometry
static GeomFactorsSharedPtr ValidateRegGeomFactor (GeomFactorsSharedPtr geomFactor)
 Check to see if a geometric factor has already been created that contains the same regular information.
 
- Static Protected Attributes inherited from Nektar::SpatialDomains::Geometry
static GeomFactorsVector m_regGeomFactorsManager
 

Detailed Description

Definition at line 59 of file SegGeom.h.

Constructor & Destructor Documentation

◆ SegGeom() [1/3]

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

◆ SegGeom() [2/3]

Nektar::SpatialDomains::SegGeom::SegGeom ( int  id,
int  coordim,
std::array< PointGeom *, kNverts vertex,
const CurveSharedPtr  curve = CurveSharedPtr() 
)

Definition at line 69 of file SegGeom.cpp.

71 : Geometry1D(coordim)
72{
74 m_globalID = id;
76 m_curve = curve;
77 m_verts = vertex;
78}
GeomState m_state
Enumeration to dictate whether coefficients are filled.
Definition Geometry.h:191
std::array< SpatialDomains::PointGeom *, kNverts > m_verts
Definition SegGeom.h:86
CurveSharedPtr m_curve
Boolean indicating whether object owns the data.
Definition SegGeom.h:103
@ 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 80 of file SegGeom.cpp.

80 : Geometry1D(in)
81{
82 // From Geometry class
83 m_shapeType = in.m_shapeType;
84
85 // info from EdgeComponent class
86 m_globalID = in.m_globalID;
87 m_xmap = in.m_xmap;
88 SetUpCoeffs(m_xmap->GetNcoeffs());
89
90 // info from SegGeom class
91 m_coordim = in.m_coordim;
92 m_verts[0] = in.m_verts[0];
93 m_verts[1] = in.m_verts[1];
94
95 m_state = in.m_state;
96}
void SetUpCoeffs(const int nCoeffs)
Initialise the Geometry::m_coeffs array.
Definition Geometry.h:704
StdRegions::StdExpansionSharedPtr m_xmap
mapping containing isoparametric transformation.
Definition Geometry.h:189
int m_coordim
Coordinate dimension of this geometry object.
Definition Geometry.h:183

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

Member Function Documentation

◆ GenerateOneSpaceDimGeom()

SegGeomUniquePtr Nektar::SpatialDomains::SegGeom::GenerateOneSpaceDimGeom ( EntityHolder1D holder)

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

123{
125
126 // info about numbering
127 returnval->m_globalID = m_globalID;
128
129 // geometric information.
130 returnval->m_coordim = 1;
131 NekDouble x0 = (*m_verts[0])[0];
133 1, m_verts[0]->GetGlobalID(), x0, 0.0, 0.0);
134 vert0->SetGlobalID(vert0->GetGlobalID());
135 returnval->m_verts[0] = vert0.get();
136 holder.m_pointVec.push_back(std::move(vert0));
137
138 // Get information to calculate length.
139 const Array<OneD, const LibUtilities::BasisSharedPtr> base =
140 m_xmap->GetBase();
142 v.push_back(base[0]->GetPointsKey());
144
145 const Array<OneD, const NekDouble> jac = m_geomFactors->GetJac(v);
146
147 NekDouble len = 0.0;
148 if (jac.size() == 1)
149 {
150 len = jac[0] * 2.0;
151 }
152 else
153 {
154 Array<OneD, const NekDouble> w0 = base[0]->GetW();
155 len = 0.0;
156
157 for (int i = 0; i < jac.size(); ++i)
158 {
159 len += jac[i] * w0[i];
160 }
161 }
162 // Set up second vertex.
164 1, m_verts[1]->GetGlobalID(), x0 + len, 0.0, 0.0);
165 vert1->SetGlobalID(vert1->GetGlobalID());
166
167 returnval->m_verts[1] = vert1.get();
168 holder.m_pointVec.push_back(std::move(vert1));
169
170 // at present just use previous m_xmap[0];
171 returnval->m_xmap = m_xmap;
172 returnval->SetUpCoeffs(m_xmap->GetNcoeffs());
173 returnval->m_state = eNotFilled;
174
175 return returnval;
176}
static std::unique_ptr< DataType, UniquePtrDeleter > AllocateUniquePtr(const Args &...args)
int GetGlobalID(void) const
Get the ID of this object.
Definition Geometry.h:322
GeomFactorsSharedPtr m_geomFactors
Geometric factors.
Definition Geometry.h:185
void v_GenGeomFactors() override
Definition SegGeom.cpp:233
std::vector< PointsKey > PointsKeyVector
Definition Points.h:231
const Array< OneD, NekDouble > base(3, 0.0)
unique_ptr_objpool< SegGeom > SegGeomUniquePtr
Definition MeshGraph.h:98
unique_ptr_objpool< PointGeom > PointGeomUniquePtr
Definition MeshGraph.h:93

References Nektar::ObjPoolManager< DataType >::AllocateUniquePtr(), Nektar::SpatialDomains::eNotFilled, Nektar::SpatialDomains::Geometry::GetGlobalID(), Nektar::SpatialDomains::Geometry::m_geomFactors, Nektar::SpatialDomains::Geometry::m_globalID, Nektar::SpatialDomains::EntityHolder1D::m_pointVec, m_verts, Nektar::SpatialDomains::Geometry::m_xmap, and v_GenGeomFactors().

◆ GetCurve()

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

Definition at line 80 of file SegGeom.h.

81 {
82 return m_curve;
83 }

References m_curve.

Referenced by export_GeomElements(), and Nektar::SpatialDomains::ZoneBase::ZoneBase().

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

211{
213
214 if ((*edge1.GetVertex(0) == *edge2.GetVertex(0)) ||
215 (*edge1.GetVertex(0) == *edge2.GetVertex(1)))
216 {
217 // Backward direction. Vertex 0 is connected to edge 2.
218 returnval = StdRegions::eBackwards;
219 }
220 else if ((*edge1.GetVertex(1) != *edge2.GetVertex(0)) &&
221 (*edge1.GetVertex(1) != *edge2.GetVertex(1)))
222 {
223 // Not forward either, then we have a problem.
224 std::ostringstream errstrm;
225 errstrm << "Connected edges do not share a vertex. Edges ";
226 errstrm << edge1.GetGlobalID() << ", " << edge2.GetGlobalID();
227 ASSERTL0(false, errstrm.str());
228 }
229
230 return returnval;
231}
#define ASSERTL0(condition, msg)

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

Referenced by Nektar::SpatialDomains::TriGeom::TriGeom().

◆ SetUpXmap()

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

Definition at line 98 of file SegGeom.cpp.

99{
100 if (m_curve)
101 {
102 int npts = m_curve->m_points.size();
103 LibUtilities::PointsKey pkey(npts + 1,
105 const LibUtilities::BasisKey B(LibUtilities::eModified_A, npts, pkey);
107 }
108 else
109 {
110 const LibUtilities::BasisKey B(
112 LibUtilities::PointsKey(2, LibUtilities::eGaussLobattoLegendre));
114 }
115}
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
@ 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 256 of file SegGeom.cpp.

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

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

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

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

234{
235 if (!m_setupState)
236 {
238 }
239
241 {
244
245 if (m_xmap->GetBasisNumModes(0) != 2)
246 {
247 gType = eDeformed;
248 }
249
251 gType, m_coordim, m_xmap, m_coeffs);
253 }
254}
bool m_setupState
Wether or not the setup routines have been run.
Definition Geometry.h:193
GeomState m_geomFactorsState
State of the geometric factors.
Definition Geometry.h:187
void v_FillGeom() override
Populate the coordinate mapping Geometry::m_coeffs information from any children geometry elements.
Definition SegGeom.cpp:256
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 183 of file SegGeom.cpp.

185{
186 ASSERTL1(m_state == ePtsFilled, "Geometry is not in physical space");
187
188 Array<OneD, NekDouble> tmp(m_xmap->GetTotPoints());
189 m_xmap->BwdTrans(m_coeffs[i], tmp);
190
191 return m_xmap->PhysEvaluate(Lcoord, tmp);
192}
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....

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

361{
362 return kNverts;
363}
static const int kNverts
Definition SegGeom.h:62

References kNverts.

◆ v_GetShapeType()

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

Definition at line 178 of file SegGeom.cpp.

179{
181}

References Nektar::LibUtilities::eSegment.

◆ v_GetVertex()

PointGeom * Nektar::SpatialDomains::SegGeom::v_GetVertex ( const int  i) const
overrideprotectedvirtual

Returns vertex i of this object.

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 348 of file SegGeom.cpp.

349{
350 PointGeom *returnval = nullptr;
351
352 if (i >= 0 && i < kNverts)
353 {
354 returnval = m_verts[i];
355 }
356
357 return returnval;
358}

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

325{
326 Geometry::v_Reset(curvedEdges, curvedFaces);
327 CurveMap::iterator it = curvedEdges.find(m_globalID);
328
329 if (it != curvedEdges.end())
330 {
331 m_curve = it->second;
332 }
333
334 SetUpXmap();
335 SetUpCoeffs(m_xmap->GetNcoeffs());
336}
virtual void v_Reset(CurveMap &curvedEdges, CurveMap &curvedFaces)
Reset this geometry object: unset the current state, zero Geometry::m_coeffs and remove allocated Geo...
Definition Geometry.cpp:372

References m_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

◆ kNfacets

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

Definition at line 63 of file SegGeom.h.

◆ kNverts

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

Definition at line 62 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 103 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 87 of file SegGeom.h.

◆ m_verts

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

Definition at line 86 of file SegGeom.h.

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