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 &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...
 
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)
 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...
 
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 Attributes

static const int kDim = 2
 

Protected Member Functions

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...
 
void NewtonIterationForLocCoord (const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &ptsx, const Array< OneD, const NekDouble > &ptsy, Array< OneD, NekDouble > &Lcoords, NekDouble &dist)
 
- 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 &dist)
 
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 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_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...
 

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...
 
Array< OneD, NekDoublem_boundingBox
 Array containing bounding box. 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.

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
Geometry()
Default constructor.
Definition: Geometry.cpp:54
int m_coordim
Coordinate dimension of this geometry object.
Definition: Geometry.h:185

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

◆ ~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 dist 
)
protected

Definition at line 65 of file Geometry2D.cpp.

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.size(), Jac, 1) /
83  ((NekDouble)Jac.size());
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.size());
93  Array<OneD, NekDouble> DxD2(ptsx.size());
94  Array<OneD, NekDouble> DyD1(ptsx.size());
95  Array<OneD, NekDouble> DyD2(ptsx.size());
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  NekDouble resid;
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( !(std::isfinite(Lcoords[0]) && std::isfinite(Lcoords[1])) )
146  {
147  dist = 1e16;
148  std::ostringstream ss;
149  ss << "nan or inf found in NewtonIterationForLocCoord in element "
150  << GetGlobalID();
151  WARNINGL1(false, ss.str());
152  return;
153  }
154  if (fabs(Lcoords[0]) > LcoordDiv || fabs(Lcoords[1]) > LcoordDiv)
155  {
156  break; // lcoords have diverged so stop iteration
157  }
158  }
159 
160  m_xmap->LocCoordToLocCollapsed(Lcoords, eta);
161  if(ClampLocCoords(eta, 0.))
162  {
163  I[0] = m_xmap->GetBasis(0)->GetI(eta);
164  I[1] = m_xmap->GetBasis(1)->GetI(eta + 1);
165  // calculate the global point corresponding to Lcoords
166  xmap = m_xmap->PhysEvaluate(I, ptsx);
167  ymap = m_xmap->PhysEvaluate(I, ptsy);
168  F1 = coords[0] - xmap;
169  F2 = coords[1] - ymap;
170  dist = sqrt(F1 * F1 + F2 * F2);
171  }
172  else
173  {
174  dist = 0.;
175  }
176 
177  if (cnt >= MaxIterations)
178  {
179  Array<OneD, NekDouble> collCoords(2);
180  m_xmap->LocCoordToLocCollapsed(Lcoords, collCoords);
181 
182  // if coordinate is inside element dump error!
183  if ((collCoords[0] >= -1.0 && collCoords[0] <= 1.0) &&
184  (collCoords[1] >= -1.0 && collCoords[1] <= 1.0))
185  {
186  std::ostringstream ss;
187 
188  ss << "Reached MaxIterations (" << MaxIterations
189  << ") in Newton iteration ";
190  ss << "Init value (" << setprecision(4) << init0 << "," << init1
191  << ","
192  << ") ";
193  ss << "Fin value (" << Lcoords[0] << "," << Lcoords[1] << ","
194  << ") ";
195  ss << "Resid = " << resid << " Tolerance = " << sqrt(ScaledTol);
196 
197  WARNINGL1(cnt < MaxIterations, ss.str());
198  }
199  }
200 }
#define WARNINGL1(condition, msg)
Definition: ErrorUtil.hpp:251
bool ClampLocCoords(Array< OneD, NekDouble > &locCoord, NekDouble tol)
Clamp local coords to be within standard regions [-1, 1]^dim.
Definition: Geometry.cpp:486
int GetGlobalID(void) const
Get the ID of this object.
Definition: Geometry.h:319
StdRegions::StdExpansionSharedPtr m_xmap
mapping containing isoparametric transformation.
Definition: Geometry.h:191
GeomFactorsSharedPtr m_geomFactors
Geometric factors.
Definition: Geometry.h:187
double NekDouble
T Vsum(int n, const T *x, const int incx)
Subtract return sum(x)
Definition: Vmath.cpp:846
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:267

References Nektar::SpatialDomains::Geometry::ClampLocCoords(), Nektar::SpatialDomains::Geometry::GetGlobalID(), Nektar::SpatialDomains::Geometry::m_geomFactors, Nektar::SpatialDomains::Geometry::m_xmap, tinysimd::sqrt(), Vmath::Vsum(), and WARNINGL1.

Referenced by v_GetLocCoords().

◆ 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 338 of file Geometry2D.cpp.

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

References ASSERTL2, and m_edges.

◆ 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 344 of file Geometry2D.cpp.

345 {
346  ASSERTL2(i >= 0 && i < m_eorient.size(), "Index out of range");
347  return m_eorient[i];
348 }
std::vector< StdRegions::Orientation > m_eorient
Definition: Geometry2D.h:80

References ASSERTL2, and m_eorient.

◆ v_GetLocCoords()

NekDouble Nektar::SpatialDomains::Geometry2D::v_GetLocCoords ( const Array< OneD, const NekDouble > &  coords,
Array< OneD, NekDouble > &  Lcoords 
)
protectedvirtual

Determine the local collapsed coordinates that correspond to a given Cartesian coordinate for this geometry object.

For curvilinear and non-affine elements (i.e. where the Jacobian varies as a function of the standard element coordinates), this is a non-linear optimisation problem that requires the use of a Newton iteration. Note therefore that this can be an expensive operation.

Note that, clearly, the provided Cartesian coordinate lie outside the element. The function therefore returns the minimum distance from some position in the element to . Lcoords will also be constrained to fit within the range \([-1,1]^d\) where \( d \) is the dimension of the element.

Parameters
coordsInput Cartesian global coordinates
LcoordsCorresponding local coordinates
Returns
Distance between obtained coordinates and provided ones.

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 202 of file Geometry2D.cpp.

204 {
205  NekDouble dist = 0.;
206  if (GetMetricInfo()->GetGtype() == eRegular)
207  {
208  int v2;
210  {
211  v2 = 2;
212  }
214  {
215  v2 = 3;
216  }
217  else
218  {
219  v2 = 2;
220  ASSERTL0(false, "unrecognized 2D element type");
221  }
222 
223  NekDouble coords2 = (m_coordim == 3) ? coords[2] : 0.0;
224  PointGeom r(m_coordim, 0, coords[0], coords[1], coords2);
225 
226  // Edges
227  PointGeom er0, e10, e20;
228  PointGeom norm, orth1, orth2;
229  er0.Sub(r, *m_verts[0]);
230  e10.Sub(*m_verts[1], *m_verts[0]);
231  e20.Sub(*m_verts[v2], *m_verts[0]);
232 
233  // Obtain normal to element plane
234  norm.Mult(e10, e20);
235  // Obtain vector which are proportional to normal of e10 and e20.
236  orth1.Mult(norm, e10);
237  orth2.Mult(norm, e20);
238 
239  // Calculate length using L/|dv1| = (x-v0).n1/(dv1.n1) for coordiante 1
240  // Then rescale to [-1,1].
241  Lcoords[0] = er0.dot(orth2) / e10.dot(orth2);
242  Lcoords[0] = 2. * Lcoords[0] - 1.;
243  Lcoords[1] = er0.dot(orth1) / e20.dot(orth1);
244  Lcoords[1] = 2. * Lcoords[1] - 1.;
245 
246  // Set distance
247  Array<OneD, NekDouble> eta(2, 0.);
248  m_xmap->LocCoordToLocCollapsed(Lcoords, eta);
249  if(ClampLocCoords(eta, 0.))
250  {
251  Array<OneD, NekDouble> xi(2, 0.);
252  m_xmap->LocCollapsedToLocCoord(eta, xi);
253  xi[0] = (xi[0] + 1.) * 0.5; //re-scaled to ratio [0, 1]
254  xi[1] = (xi[1] + 1.) * 0.5;
255  for (int i = 0; i < m_coordim; ++i)
256  {
257  NekDouble tmp = xi[0]*e10[i] + xi[1]*e20[i] - er0[i];
258  dist += tmp * tmp;
259  }
260  dist = sqrt(dist);
261  }
262  }
263  else
264  {
265  v_FillGeom();
266  // Determine nearest point of coords to values in m_xmap
267  int npts = m_xmap->GetTotPoints();
268  Array<OneD, NekDouble> ptsx(npts), ptsy(npts);
269  Array<OneD, NekDouble> tmpx(npts), tmpy(npts);
270 
271  // Determine 3D manifold orientation
272  int idx = 0, idy = 1;
273  if(m_coordim == 3)
274  {
275  PointGeom e01, e21, norm;
276  e01.Sub(*m_verts[0], *m_verts[1]);
277  e21.Sub(*m_verts[2], *m_verts[1]);
278  norm.Mult(e01, e21);
279  int tmpi = 0;
280  double tmp = std::fabs(norm[0]);
281  if(tmp < fabs(norm[1]))
282  {
283  tmp = fabs(norm[1]);
284  tmpi = 1;
285  }
286  if(tmp < fabs(norm[2]))
287  {
288  tmpi = 2;
289  }
290  idx = (tmpi + 1) % 3;
291  idy = (tmpi + 2) % 3;
292  }
293  Array<OneD, NekDouble> tmpcoords(2);
294  tmpcoords[0] = coords[idx];
295  tmpcoords[1] = coords[idy];
296 
297  m_xmap->BwdTrans(m_coeffs[idx], ptsx);
298  m_xmap->BwdTrans(m_coeffs[idy], ptsy);
299 
300  const Array<OneD, const NekDouble> za = m_xmap->GetPoints(0);
301  const Array<OneD, const NekDouble> zb = m_xmap->GetPoints(1);
302 
303  // guess the first local coords based on nearest point
304  Vmath::Sadd(npts, -tmpcoords[0], ptsx, 1, tmpx, 1);
305  Vmath::Sadd(npts, -tmpcoords[1], ptsy, 1, tmpy, 1);
306  Vmath::Vmul(npts, tmpx, 1, tmpx, 1, tmpx, 1);
307  Vmath::Vvtvp(npts, tmpy, 1, tmpy, 1, tmpx, 1, tmpx, 1);
308 
309  int min_i = Vmath::Imin(npts, tmpx, 1);
310 
311  Array<OneD, NekDouble> eta(2, 0.);
312  eta[0] = za[min_i % za.size()];
313  eta[1] = zb[min_i / za.size()];
314  m_xmap->LocCollapsedToLocCoord(eta, Lcoords);
315 
316  // Perform newton iteration to find local coordinates
317  NewtonIterationForLocCoord(tmpcoords, ptsx, ptsy, Lcoords, dist);
318  }
319  return dist;
320 }
void NewtonIterationForLocCoord(const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &ptsx, const Array< OneD, const NekDouble > &ptsy, Array< OneD, NekDouble > &Lcoords, NekDouble &dist)
Definition: Geometry2D.cpp:65
virtual void v_FillGeom()
Populate the coordinate mapping Geometry::m_coeffs information from any children geometry elements.
Definition: Geometry.cpp:346
LibUtilities::ShapeType m_shapeType
Type of shape.
Definition: Geometry.h:199
Array< OneD, Array< OneD, NekDouble > > m_coeffs
Array containing expansion coefficients of m_xmap.
Definition: Geometry.h:203
GeomFactorsSharedPtr GetMetricInfo()
Get the geometric factors for this object.
Definition: Geometry.h:303
@ eRegular
Geometry is straight-sided with constant geometric factors.
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:192
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.cpp:513
int Imin(int n, const T *x, const int incx)
Return the index of the minimum element in x.
Definition: Vmath.cpp:966
void Sadd(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Add vector y = alpha - x.
Definition: Vmath.cpp:341

References ASSERTL0, Nektar::SpatialDomains::Geometry::ClampLocCoords(), Nektar::SpatialDomains::PointGeom::dot(), Nektar::LibUtilities::eQuadrilateral, Nektar::SpatialDomains::eRegular, Nektar::LibUtilities::eTriangle, Nektar::SpatialDomains::Geometry::GetMetricInfo(), Vmath::Imin(), Nektar::SpatialDomains::Geometry::m_coeffs, Nektar::SpatialDomains::Geometry::m_coordim, Nektar::SpatialDomains::Geometry::m_shapeType, m_verts, Nektar::SpatialDomains::Geometry::m_xmap, Nektar::SpatialDomains::PointGeom::Mult(), NewtonIterationForLocCoord(), Vmath::Sadd(), tinysimd::sqrt(), Nektar::SpatialDomains::PointGeom::Sub(), Nektar::SpatialDomains::Geometry::v_FillGeom(), Vmath::Vmul(), and Vmath::Vvtvp().

◆ 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 327 of file Geometry2D.cpp.

328 {
329  return m_edges.size();
330 }

References m_edges.

◆ 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 322 of file Geometry2D.cpp.

323 {
324  return m_verts.size();
325 }

References m_verts.

◆ 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 350 of file Geometry2D.cpp.

351 {
352  return 2;
353 }

◆ v_GetVertex()

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

Implements Nektar::SpatialDomains::Geometry.

Definition at line 332 of file Geometry2D.cpp.

333 {
334  ASSERTL2(i >= 0 && i < m_verts.size(), "Index out of range");
335  return m_verts[i];
336 }

References ASSERTL2, and m_verts.

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