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

2D geometry information More...

#include <Geometry2D.h>

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

Public Member Functions

 Geometry2D ()
 
 Geometry2D (const int coordim, CurveSharedPtr curve)
 
virtual ~Geometry2D ()
 
- 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...
 
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 &resid)
 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...
 
bool MinMaxCheck (const Array< OneD, const NekDouble > &gloCoord)
 Check if given global coord is within twice the min/max distance of the element. More...
 
void ClampLocCoords (Array< OneD, NekDouble > &locCoord, NekDouble tol)
 Clamp local coords to be within standard regions [-1, 1]^dim. More...
 
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...
 
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 Attributes

static const int kDim = 2
 

Protected Member Functions

void NewtonIterationForLocCoord (const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &ptsx, const Array< OneD, const NekDouble > &ptsy, Array< OneD, NekDouble > &Lcoords, NekDouble &resid)
 
- Protected Member Functions inherited from Nektar::SpatialDomains::Geometry
void GenGeomFactors ()
 Handles generation of geometry factors. More...
 
virtual Geometry2DSharedPtr v_GetFace (int i) const
 Returns face i of this object. 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_GetNumFaces () const
 Get the number of faces of this object. 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 &resid)
 
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 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 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...
 

Protected Attributes

PointGeomVector m_verts
 
SegGeomVector m_edges
 
std::vector< StdRegions::Orientationm_eorient
 
CurveSharedPtr m_curve
 
- 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...
 

Private Member Functions

virtual int v_GetShapeDim () const
 Get the object's shape dimension. More...
 
virtual PointGeomSharedPtr v_GetVertex (int i) const
 
virtual Geometry1DSharedPtr v_GetEdge (int i) const
 Returns edge i of this object. 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 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...
 

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

2D geometry information

Definition at line 68 of file Geometry2D.h.

Constructor & Destructor Documentation

◆ Geometry2D() [1/2]

Nektar::SpatialDomains::Geometry2D::Geometry2D ( )

Definition at line 50 of file Geometry2D.cpp.

51 {
52 }

◆ Geometry2D() [2/2]

Nektar::SpatialDomains::Geometry2D::Geometry2D ( const int  coordim,
CurveSharedPtr  curve 
)

Definition at line 54 of file Geometry2D.cpp.

References ASSERTL0, and Nektar::SpatialDomains::Geometry::m_coordim.

55  : Geometry(coordim), m_curve(curve)
56 {
57  ASSERTL0(m_coordim > 1,
58  "Coordinate dimension should be at least 2 for a 2D geometry");
59 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
int m_coordim
Coordinate dimension of this geometry object.
Definition: Geometry.h:183
Geometry()
Default constructor.
Definition: Geometry.cpp:54

◆ ~Geometry2D()

Nektar::SpatialDomains::Geometry2D::~Geometry2D ( )
virtual

Definition at line 61 of file Geometry2D.cpp.

62 {
63 }

Member Function Documentation

◆ NewtonIterationForLocCoord()

void Nektar::SpatialDomains::Geometry2D::NewtonIterationForLocCoord ( const Array< OneD, const NekDouble > &  coords,
const Array< OneD, const NekDouble > &  ptsx,
const Array< OneD, const NekDouble > &  ptsy,
Array< OneD, NekDouble > &  Lcoords,
NekDouble resid 
)
protected

Definition at line 65 of file Geometry2D.cpp.

References Nektar::SpatialDomains::Geometry::m_geomFactors, Nektar::SpatialDomains::Geometry::m_xmap, Vmath::Vsum(), and WARNINGL1.

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

71 {
72  // Maximum iterations for convergence
73  const int MaxIterations = 51;
74  // |x-xp|^2 < EPSILON error tolerance
75  const NekDouble Tol = 1.e-8;
76  // |r,s| > LcoordDIV stop the search
77  const NekDouble LcoordDiv = 15.0;
78 
79  Array<OneD, const NekDouble> Jac =
80  m_geomFactors->GetJac(m_xmap->GetPointsKeys());
81 
82  NekDouble ScaledTol = Vmath::Vsum(Jac.num_elements(), Jac, 1) /
83  ((NekDouble)Jac.num_elements());
84  ScaledTol *= Tol;
85 
86  NekDouble xmap, ymap, F1, F2;
87  NekDouble derx_1, derx_2, dery_1, dery_2, jac;
88 
89  // save intiial guess for later reference if required.
90  NekDouble init0 = Lcoords[0], init1 = Lcoords[1];
91 
92  Array<OneD, NekDouble> DxD1(ptsx.num_elements());
93  Array<OneD, NekDouble> DxD2(ptsx.num_elements());
94  Array<OneD, NekDouble> DyD1(ptsx.num_elements());
95  Array<OneD, NekDouble> DyD2(ptsx.num_elements());
96 
97  // Ideally this will be stored in m_geomfactors
98  m_xmap->PhysDeriv(ptsx, DxD1, DxD2);
99  m_xmap->PhysDeriv(ptsy, DyD1, DyD2);
100 
101  int cnt = 0;
102  Array<OneD, DNekMatSharedPtr> I(2);
103  Array<OneD, NekDouble> eta(2);
104 
105  F1 = F2 = 2000; // Starting value of Function
106 
107  while (cnt++ < MaxIterations)
108  {
109  // evaluate lagrange interpolant at Lcoords
110  m_xmap->LocCoordToLocCollapsed(Lcoords, eta);
111  I[0] = m_xmap->GetBasis(0)->GetI(eta);
112  I[1] = m_xmap->GetBasis(1)->GetI(eta + 1);
113 
114  // calculate the global point `corresponding to Lcoords
115  xmap = m_xmap->PhysEvaluate(I, ptsx);
116  ymap = m_xmap->PhysEvaluate(I, ptsy);
117 
118  F1 = coords[0] - xmap;
119  F2 = coords[1] - ymap;
120 
121  if (F1 * F1 + F2 * F2 < ScaledTol)
122  {
123  resid = sqrt(F1 * F1 + F2 * F2);
124  break;
125  }
126 
127  // Interpolate derivative metric at Lcoords
128  derx_1 = m_xmap->PhysEvaluate(I, DxD1);
129  derx_2 = m_xmap->PhysEvaluate(I, DxD2);
130  dery_1 = m_xmap->PhysEvaluate(I, DyD1);
131  dery_2 = m_xmap->PhysEvaluate(I, DyD2);
132 
133  jac = dery_2 * derx_1 - dery_1 * derx_2;
134 
135  // use analytical inverse of derivitives which are
136  // also similar to those of metric factors.
137  Lcoords[0] =
138  Lcoords[0] +
139  (dery_2 * (coords[0] - xmap) - derx_2 * (coords[1] - ymap)) / jac;
140 
141  Lcoords[1] =
142  Lcoords[1] +
143  (-dery_1 * (coords[0] - xmap) + derx_1 * (coords[1] - ymap)) / jac;
144 
145  if (fabs(Lcoords[0]) > LcoordDiv || fabs(Lcoords[1]) > LcoordDiv)
146  {
147  break; // lcoords have diverged so stop iteration
148  }
149  }
150 
151  resid = sqrt(F1 * F1 + F2 * F2);
152 
153  if (cnt >= MaxIterations)
154  {
155  Array<OneD, NekDouble> collCoords(2);
156  m_xmap->LocCoordToLocCollapsed(Lcoords, collCoords);
157 
158  // if coordinate is inside element dump error!
159  if ((collCoords[0] >= -1.0 && collCoords[0] <= 1.0) &&
160  (collCoords[1] >= -1.0 && collCoords[1] <= 1.0))
161  {
162  std::ostringstream ss;
163 
164  ss << "Reached MaxIterations (" << MaxIterations
165  << ") in Newton iteration ";
166  ss << "Init value (" << setprecision(4) << init0 << "," << init1
167  << ","
168  << ") ";
169  ss << "Fin value (" << Lcoords[0] << "," << Lcoords[1] << ","
170  << ") ";
171  ss << "Resid = " << resid << " Tolerance = " << sqrt(ScaledTol);
172 
173  WARNINGL1(cnt < MaxIterations, ss.str());
174  }
175  }
176 }
StdRegions::StdExpansionSharedPtr m_xmap
mapping containing isoparametric transformation.
Definition: Geometry.h:189
GeomFactorsSharedPtr m_geomFactors
Geometric factors.
Definition: Geometry.h:185
double NekDouble
#define WARNINGL1(condition, msg)
Definition: ErrorUtil.hpp:251
T Vsum(int n, const T *x, const int incx)
Subtract return sum(x)
Definition: Vmath.cpp:740

◆ v_GetEdge()

Geometry1DSharedPtr Nektar::SpatialDomains::Geometry2D::v_GetEdge ( int  i) const
privatevirtual

Returns edge i of this object.

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 194 of file Geometry2D.cpp.

References ASSERTL2, and m_edges.

195 {
196  ASSERTL2(i >= 0 && i < m_edges.size(), "Index out of range");
197  return m_edges[i];
198 }
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed t...
Definition: ErrorUtil.hpp:274

◆ v_GetEorient()

StdRegions::Orientation Nektar::SpatialDomains::Geometry2D::v_GetEorient ( const int  i) const
privatevirtual

Returns the orientation of edge i with respect to the ordering of edges in the standard element.

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 200 of file Geometry2D.cpp.

References ASSERTL2, and m_eorient.

201 {
202  ASSERTL2(i >= 0 && i < m_eorient.size(), "Index out of range");
203  return m_eorient[i];
204 }
std::vector< StdRegions::Orientation > m_eorient
Definition: Geometry2D.h:80
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed t...
Definition: ErrorUtil.hpp:274

◆ v_GetNumEdges()

int Nektar::SpatialDomains::Geometry2D::v_GetNumEdges ( ) const
privatevirtual

Get the number of edges of this object.

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 183 of file Geometry2D.cpp.

References m_edges.

184 {
185  return m_edges.size();
186 }

◆ v_GetNumVerts()

int Nektar::SpatialDomains::Geometry2D::v_GetNumVerts ( ) const
privatevirtual

Get the number of vertices of this object.

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 178 of file Geometry2D.cpp.

References m_verts.

179 {
180  return m_verts.size();
181 }

◆ v_GetShapeDim()

int Nektar::SpatialDomains::Geometry2D::v_GetShapeDim ( ) const
privatevirtual

Get the object's shape dimension.

For example, a segment is one dimensional and quadrilateral is two dimensional.

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 206 of file Geometry2D.cpp.

207 {
208  return 2;
209 }

◆ v_GetVertex()

PointGeomSharedPtr Nektar::SpatialDomains::Geometry2D::v_GetVertex ( int  i) const
privatevirtual

Implements Nektar::SpatialDomains::Geometry.

Definition at line 188 of file Geometry2D.cpp.

References ASSERTL2, and m_verts.

189 {
190  ASSERTL2(i >= 0 && i < m_verts.size(), "Index out of range");
191  return m_verts[i];
192 }
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed t...
Definition: ErrorUtil.hpp:274

Member Data Documentation

◆ kDim

const int Nektar::SpatialDomains::Geometry2D::kDim = 2
static

Definition at line 75 of file Geometry2D.h.

◆ m_curve

CurveSharedPtr Nektar::SpatialDomains::Geometry2D::m_curve
protected

◆ m_edges

SegGeomVector Nektar::SpatialDomains::Geometry2D::m_edges
protected

◆ m_eorient

std::vector<StdRegions::Orientation> Nektar::SpatialDomains::Geometry2D::m_eorient
protected

◆ m_verts

PointGeomVector Nektar::SpatialDomains::Geometry2D::m_verts
protected