Nektar++
Geometry2D.cpp
Go to the documentation of this file.
1 ////////////////////////////////////////////////////////////////////////////////
2 //
3 // File: Geometry2D.cpp
4 //
5 // For more information, please see: http://www.nektar.info/
6 //
7 // The MIT License
8 //
9 // Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10 // Department of Aeronautics, Imperial College London (UK), and Scientific
11 // Computing and Imaging Institute, University of Utah (USA).
12 //
13 // Permission is hereby granted, free of charge, to any person obtaining a
14 // copy of this software and associated documentation files (the "Software"),
15 // to deal in the Software without restriction, including without limitation
16 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
17 // and/or sell copies of the Software, and to permit persons to whom the
18 // Software is furnished to do so, subject to the following conditions:
19 //
20 // The above copyright notice and this permission notice shall be included
21 // in all copies or substantial portions of the Software.
22 //
23 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29 // DEALINGS IN THE SOFTWARE.
30 //
31 // Description: 2D geometry information
32 //
33 //
34 ////////////////////////////////////////////////////////////////////////////////
35 
36 #include <iomanip>
37 
39 #include <SpatialDomains/SegGeom.h>
40 
41 #include <iomanip>
42 
43 using namespace std;
44 
45 namespace Nektar
46 {
47 namespace SpatialDomains
48 {
49 
50 Geometry2D::Geometry2D()
51 {
52 }
53 
54 Geometry2D::Geometry2D(const int coordim, CurveSharedPtr curve)
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 }
60 
62 {
63 }
64 
66  const Array<OneD, const NekDouble> &coords,
67  const Array<OneD, const NekDouble> &ptsx,
69  NekDouble &dist)
70 {
71  // Maximum iterations for convergence
72  const int MaxIterations = 51;
73  // |x-xp|^2 < EPSILON error tolerance
74  const NekDouble Tol = 1.e-8;
75  // |r,s| > LcoordDIV stop the search
76  const NekDouble LcoordDiv = 15.0;
77 
79  m_geomFactors->GetJac(m_xmap->GetPointsKeys());
80 
81  NekDouble ScaledTol =
82  Vmath::Vsum(Jac.size(), Jac, 1) / ((NekDouble)Jac.size());
83  ScaledTol *= Tol;
84 
85  NekDouble xmap, ymap, F1, F2;
86  NekDouble derx_1, derx_2, dery_1, dery_2, jac;
87 
88  // save intiial guess for later reference if required.
89  NekDouble init0 = Lcoords[0], init1 = Lcoords[1];
90 
91  Array<OneD, NekDouble> DxD1(ptsx.size());
92  Array<OneD, NekDouble> DxD2(ptsx.size());
93  Array<OneD, NekDouble> DyD1(ptsx.size());
94  Array<OneD, NekDouble> DyD2(ptsx.size());
95 
96  // Ideally this will be stored in m_geomfactors
97  m_xmap->PhysDeriv(ptsx, DxD1, DxD2);
98  m_xmap->PhysDeriv(ptsy, DyD1, DyD2);
99 
100  int cnt = 0;
102  Array<OneD, NekDouble> eta(2);
103 
104  F1 = F2 = 2000; // Starting value of Function
105  NekDouble resid;
106  while (cnt++ < MaxIterations)
107  {
108  // evaluate lagrange interpolant at Lcoords
109  m_xmap->LocCoordToLocCollapsed(Lcoords, eta);
110  I[0] = m_xmap->GetBasis(0)->GetI(eta);
111  I[1] = m_xmap->GetBasis(1)->GetI(eta + 1);
112 
113  // calculate the global point `corresponding to Lcoords
114  xmap = m_xmap->PhysEvaluate(I, ptsx);
115  ymap = m_xmap->PhysEvaluate(I, ptsy);
116 
117  F1 = coords[0] - xmap;
118  F2 = coords[1] - ymap;
119 
120  if (F1 * F1 + F2 * F2 < ScaledTol)
121  {
122  resid = sqrt(F1 * F1 + F2 * F2);
123  break;
124  }
125 
126  // Interpolate derivative metric at Lcoords
127  derx_1 = m_xmap->PhysEvaluate(I, DxD1);
128  derx_2 = m_xmap->PhysEvaluate(I, DxD2);
129  dery_1 = m_xmap->PhysEvaluate(I, DyD1);
130  dery_2 = m_xmap->PhysEvaluate(I, DyD2);
131 
132  jac = dery_2 * derx_1 - dery_1 * derx_2;
133 
134  // use analytical inverse of derivitives which are
135  // also similar to those of metric factors.
136  Lcoords[0] =
137  Lcoords[0] +
138  (dery_2 * (coords[0] - xmap) - derx_2 * (coords[1] - ymap)) / jac;
139 
140  Lcoords[1] =
141  Lcoords[1] +
142  (-dery_1 * (coords[0] - xmap) + derx_1 * (coords[1] - ymap)) / jac;
143 
144  if (!(std::isfinite(Lcoords[0]) && std::isfinite(Lcoords[1])))
145  {
146  dist = 1e16;
147  std::ostringstream ss;
148  ss << "nan or inf found in NewtonIterationForLocCoord in element "
149  << GetGlobalID();
150  WARNINGL1(false, ss.str());
151  return;
152  }
153  if (fabs(Lcoords[0]) > LcoordDiv || fabs(Lcoords[1]) > LcoordDiv)
154  {
155  break; // lcoords have diverged so stop iteration
156  }
157  }
158 
159  m_xmap->LocCoordToLocCollapsed(Lcoords, eta);
160  if (ClampLocCoords(eta, 0.))
161  {
162  I[0] = m_xmap->GetBasis(0)->GetI(eta);
163  I[1] = m_xmap->GetBasis(1)->GetI(eta + 1);
164  // calculate the global point corresponding to Lcoords
165  xmap = m_xmap->PhysEvaluate(I, ptsx);
166  ymap = m_xmap->PhysEvaluate(I, ptsy);
167  F1 = coords[0] - xmap;
168  F2 = coords[1] - ymap;
169  dist = sqrt(F1 * F1 + F2 * F2);
170  }
171  else
172  {
173  dist = 0.;
174  }
175 
176  if (cnt >= MaxIterations)
177  {
178  Array<OneD, NekDouble> collCoords(2);
179  m_xmap->LocCoordToLocCollapsed(Lcoords, collCoords);
180 
181  // if coordinate is inside element dump error!
182  if ((collCoords[0] >= -1.0 && collCoords[0] <= 1.0) &&
183  (collCoords[1] >= -1.0 && collCoords[1] <= 1.0))
184  {
185  std::ostringstream ss;
186 
187  ss << "Reached MaxIterations (" << MaxIterations
188  << ") in Newton iteration ";
189  ss << "Init value (" << setprecision(4) << init0 << "," << init1
190  << ","
191  << ") ";
192  ss << "Fin value (" << Lcoords[0] << "," << Lcoords[1] << ","
193  << ") ";
194  ss << "Resid = " << resid << " Tolerance = " << sqrt(ScaledTol);
195 
196  WARNINGL1(cnt < MaxIterations, ss.str());
197  }
198  }
199 }
200 
202  Array<OneD, NekDouble> &Lcoords)
203 {
204  NekDouble dist = 0.;
205  if (GetMetricInfo()->GetGtype() == eRegular)
206  {
207  int v2;
209  {
210  v2 = 2;
211  }
213  {
214  v2 = 3;
215  }
216  else
217  {
218  v2 = 2;
219  ASSERTL0(false, "unrecognized 2D element type");
220  }
221 
222  NekDouble coords2 = (m_coordim == 3) ? coords[2] : 0.0;
223  PointGeom r(m_coordim, 0, coords[0], coords[1], coords2);
224 
225  // Edges
226  PointGeom er0, e10, e20;
227  PointGeom norm, orth1, orth2;
228  er0.Sub(r, *m_verts[0]);
229  e10.Sub(*m_verts[1], *m_verts[0]);
230  e20.Sub(*m_verts[v2], *m_verts[0]);
231 
232  // Obtain normal to element plane
233  norm.Mult(e10, e20);
234  // Obtain vector which are proportional to normal of e10 and e20.
235  orth1.Mult(norm, e10);
236  orth2.Mult(norm, e20);
237 
238  // Calculate length using L/|dv1| = (x-v0).n1/(dv1.n1) for coordiante 1
239  // Then rescale to [-1,1].
240  Lcoords[0] = er0.dot(orth2) / e10.dot(orth2);
241  Lcoords[0] = 2. * Lcoords[0] - 1.;
242  Lcoords[1] = er0.dot(orth1) / e20.dot(orth1);
243  Lcoords[1] = 2. * Lcoords[1] - 1.;
244 
245  // Set distance
246  Array<OneD, NekDouble> eta(2, 0.);
247  m_xmap->LocCoordToLocCollapsed(Lcoords, eta);
248  if (ClampLocCoords(eta, 0.))
249  {
250  Array<OneD, NekDouble> xi(2, 0.);
251  m_xmap->LocCollapsedToLocCoord(eta, xi);
252  xi[0] = (xi[0] + 1.) * 0.5; // re-scaled to ratio [0, 1]
253  xi[1] = (xi[1] + 1.) * 0.5;
254  for (int i = 0; i < m_coordim; ++i)
255  {
256  NekDouble tmp = xi[0] * e10[i] + xi[1] * e20[i] - er0[i];
257  dist += tmp * tmp;
258  }
259  dist = sqrt(dist);
260  }
261  }
262  else
263  {
264  v_FillGeom();
265  // Determine nearest point of coords to values in m_xmap
266  int npts = m_xmap->GetTotPoints();
267  Array<OneD, NekDouble> ptsx(npts), ptsy(npts);
268  Array<OneD, NekDouble> tmpx(npts), tmpy(npts);
269 
270  // Determine 3D manifold orientation
271  int idx = 0, idy = 1;
272  if (m_coordim == 3)
273  {
274  PointGeom e01, e21, norm;
275  e01.Sub(*m_verts[0], *m_verts[1]);
276  e21.Sub(*m_verts[2], *m_verts[1]);
277  norm.Mult(e01, e21);
278  int tmpi = 0;
279  double tmp = std::fabs(norm[0]);
280  if (tmp < fabs(norm[1]))
281  {
282  tmp = fabs(norm[1]);
283  tmpi = 1;
284  }
285  if (tmp < fabs(norm[2]))
286  {
287  tmpi = 2;
288  }
289  idx = (tmpi + 1) % 3;
290  idy = (tmpi + 2) % 3;
291  }
292  Array<OneD, NekDouble> tmpcoords(2);
293  tmpcoords[0] = coords[idx];
294  tmpcoords[1] = coords[idy];
295 
296  m_xmap->BwdTrans(m_coeffs[idx], ptsx);
297  m_xmap->BwdTrans(m_coeffs[idy], ptsy);
298 
299  const Array<OneD, const NekDouble> za = m_xmap->GetPoints(0);
300  const Array<OneD, const NekDouble> zb = m_xmap->GetPoints(1);
301 
302  // guess the first local coords based on nearest point
303  Vmath::Sadd(npts, -tmpcoords[0], ptsx, 1, tmpx, 1);
304  Vmath::Sadd(npts, -tmpcoords[1], ptsy, 1, tmpy, 1);
305  Vmath::Vmul(npts, tmpx, 1, tmpx, 1, tmpx, 1);
306  Vmath::Vvtvp(npts, tmpy, 1, tmpy, 1, tmpx, 1, tmpx, 1);
307 
308  int min_i = Vmath::Imin(npts, tmpx, 1);
309 
310  Array<OneD, NekDouble> eta(2, 0.);
311  eta[0] = za[min_i % za.size()];
312  eta[1] = zb[min_i / za.size()];
313  m_xmap->LocCollapsedToLocCoord(eta, Lcoords);
314 
315  // Perform newton iteration to find local coordinates
316  NewtonIterationForLocCoord(tmpcoords, ptsx, ptsy, Lcoords, dist);
317  }
318  return dist;
319 }
320 
322 {
323  return m_verts.size();
324 }
325 
327 {
328  return m_edges.size();
329 }
330 
332 {
333  ASSERTL2(i >= 0 && i < m_verts.size(), "Index out of range");
334  return m_verts[i];
335 }
336 
338 {
339  ASSERTL2(i >= 0 && i < m_edges.size(), "Index out of range");
340  return m_edges[i];
341 }
342 
344 {
345  ASSERTL2(i >= 0 && i < m_eorient.size(), "Index out of range");
346  return m_eorient[i];
347 }
348 
350 {
351  return 2;
352 }
353 
354 } // namespace SpatialDomains
355 } // namespace Nektar
#define WARNINGL1(condition, msg)
Definition: ErrorUtil.hpp:250
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed to...
Definition: ErrorUtil.hpp:272
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 Geometry1DSharedPtr v_GetEdge(int i) const
Returns edge i of this object.
Definition: Geometry2D.cpp:337
std::vector< StdRegions::Orientation > m_eorient
Definition: Geometry2D.h:80
virtual int v_GetShapeDim() const
Get the object's shape dimension.
Definition: Geometry2D.cpp:349
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.
Definition: Geometry2D.cpp:343
virtual int v_GetNumVerts() const
Get the number of vertices of this object.
Definition: Geometry2D.cpp:321
virtual int v_GetNumEdges() const
Get the number of edges of this object.
Definition: Geometry2D.cpp:326
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 ge...
Definition: Geometry2D.cpp:201
virtual PointGeomSharedPtr v_GetVertex(int i) const
Definition: Geometry2D.cpp:331
Base class for shape geometry information.
Definition: Geometry.h:82
bool ClampLocCoords(Array< OneD, NekDouble > &locCoord, NekDouble tol)
Clamp local coords to be within standard regions [-1, 1]^dim.
Definition: Geometry.cpp:496
virtual void v_FillGeom()
Populate the coordinate mapping Geometry::m_coeffs information from any children geometry elements.
Definition: Geometry.cpp:358
int GetGlobalID(void) const
Get the ID of this object.
Definition: Geometry.h:314
LibUtilities::ShapeType m_shapeType
Type of shape.
Definition: Geometry.h:194
Array< OneD, Array< OneD, NekDouble > > m_coeffs
Array containing expansion coefficients of m_xmap.
Definition: Geometry.h:198
GeomFactorsSharedPtr GetMetricInfo()
Get the geometric factors for this object.
Definition: Geometry.h:298
StdRegions::StdExpansionSharedPtr m_xmap
mapping containing isoparametric transformation.
Definition: Geometry.h:186
GeomFactorsSharedPtr m_geomFactors
Geometric factors.
Definition: Geometry.h:182
int m_coordim
Coordinate dimension of this geometry object.
Definition: Geometry.h:180
void Sub(PointGeom &a, PointGeom &b)
Definition: PointGeom.cpp:124
void Mult(PointGeom &a, PointGeom &b)
_this = a x b
Definition: PointGeom.cpp:133
NekDouble dot(PointGeom &a)
retun the dot product between this and input a
Definition: PointGeom.cpp:187
std::shared_ptr< Curve > CurveSharedPtr
Definition: Curve.hpp:60
@ eRegular
Geometry is straight-sided with constant geometric factors.
std::shared_ptr< PointGeom > PointGeomSharedPtr
Definition: Geometry.h:59
std::shared_ptr< Geometry1D > Geometry1DSharedPtr
Definition: Geometry.h:63
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:1
double NekDouble
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:209
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:574
T Vsum(int n, const T *x, const int incx)
Subtract return sum(x)
Definition: Vmath.cpp:895
int Imin(int n, const T *x, const int incx)
Return the index of the minimum element in x.
Definition: Vmath.cpp:1023
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:384
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:291