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
40
41#include <iomanip>
42
43using namespace std;
44
46{
47
49{
50}
51
52Geometry2D::Geometry2D(const int coordim, CurveSharedPtr curve)
53 : Geometry(coordim), m_curve(curve)
54{
56 "Coordinate dimension should be at least 2 for a 2D geometry");
57}
58
60{
61}
62
64{
65 int nc = 1, d0 = m_manifold[0], d1 = m_manifold[1];
66 if (0 == m_edgeNormal.size())
67 {
70 x[0] = Array<OneD, NekDouble>(3);
71 x[1] = Array<OneD, NekDouble>(3);
72 m_verts[m_verts.size() - 1]->GetCoords(x[0]);
73 int i0 = 1, i1 = 0;
74 for (size_t i = 0; i < m_verts.size(); ++i)
75 {
76 i0 ^= 1;
77 i1 ^= 1;
78 m_verts[i]->GetCoords(x[i1]);
79 if (m_edges[i]->GetXmap()->GetBasis(0)->GetNumModes() > 2)
80 {
81 continue;
82 }
84 m_edgeNormal[i][0] = x[i0][d1] - x[i1][d1];
85 m_edgeNormal[i][1] = x[i1][d0] - x[i0][d0];
86 }
87 }
88
89 Array<OneD, NekDouble> vertex(3);
90 for (size_t i = 0; i < m_verts.size(); ++i)
91 {
92 m_verts[i]->GetCoords(vertex);
93 if (m_edgeNormal[i].size() == 0)
94 {
95 nc = 0; // not sure
96 continue;
97 }
98 NekDouble value = m_edgeNormal[i][0] * (gloCoord[d0] - vertex[d0]) +
99 m_edgeNormal[i][1] * (gloCoord[d1] - vertex[d1]);
100 if (value < 0)
101 {
102 return -1; // outside
103 }
104 }
105 // nc: 1 (side element), 0 (maybe inside), -1 (outside)
106 return nc;
107}
108
111{
112 int i0 = 0, i1 = 1, j1 = 0, j2 = 1;
113 if (m_straightEdge & 2)
114 {
115 i0 = 1;
116 i1 = 0;
117 }
118 if (m_straightEdge & 4)
119 {
120 j1 = 1;
121 j2 = 0;
122 }
124 NekDouble gamma = m_isoParameter[1][3];
125 NekDouble tty = (coords[i1] - gamma * coords[i0] - m_isoParameter[1][0]) *
126 m_isoParameter[1][1];
127 NekDouble denom = 1. / (m_isoParameter[0][2] + m_isoParameter[0][3] * tty);
128 NekDouble epsilon = -m_isoParameter[0][3] * beta * denom;
129 NekDouble h =
130 (m_isoParameter[0][0] + m_isoParameter[0][1] * tty - coords[i0]) *
131 denom;
132 Lcoords[j2] = -h / (0.5 + sqrt(0.25 - epsilon * h));
133 Lcoords[j1] = -beta * Lcoords[j2] + tty;
134}
135
137 const Array<OneD, const NekDouble> &coords,
140 NekDouble &dist)
141{
142 // Maximum iterations for convergence
143 const int MaxIterations = NekConstants::kNewtonIterations;
144 // |x-xp|^2 < EPSILON error tolerance
145 const NekDouble Tol = 1.e-8;
146 // |r,s| > LcoordDIV stop the search
147 const NekDouble LcoordDiv = 15.0;
148
150 m_geomFactors->GetJac(m_xmap->GetPointsKeys());
151
152 NekDouble ScaledTol =
153 Vmath::Vsum(Jac.size(), Jac, 1) / ((NekDouble)Jac.size());
154 ScaledTol *= Tol;
155
156 NekDouble xmap, ymap, F1, F2;
157 NekDouble derx_1, derx_2, dery_1, dery_2, jac;
158
159 // save intiial guess for later reference if required.
160 NekDouble init0 = Lcoords[0], init1 = Lcoords[1];
161
162 Array<OneD, NekDouble> DxD1(ptsx.size());
163 Array<OneD, NekDouble> DxD2(ptsx.size());
164 Array<OneD, NekDouble> DyD1(ptsx.size());
165 Array<OneD, NekDouble> DyD2(ptsx.size());
166
167 // Ideally this will be stored in m_geomfactors
168 m_xmap->PhysDeriv(ptsx, DxD1, DxD2);
169 m_xmap->PhysDeriv(ptsy, DyD1, DyD2);
170
171 int cnt = 0;
174
175 F1 = F2 = 2000; // Starting value of Function
176 NekDouble resid = sqrt(F1 * F1 + F2 * F2);
177 while (cnt++ < MaxIterations)
178 {
179 // evaluate lagrange interpolant at Lcoords
180 m_xmap->LocCoordToLocCollapsed(Lcoords, eta);
181 I[0] = m_xmap->GetBasis(0)->GetI(eta);
182 I[1] = m_xmap->GetBasis(1)->GetI(eta + 1);
183
184 // calculate the global point `corresponding to Lcoords
185 xmap = m_xmap->PhysEvaluate(I, ptsx);
186 ymap = m_xmap->PhysEvaluate(I, ptsy);
187
188 F1 = coords[0] - xmap;
189 F2 = coords[1] - ymap;
190
191 if (F1 * F1 + F2 * F2 < ScaledTol)
192 {
193 resid = sqrt(F1 * F1 + F2 * F2);
194 break;
195 }
196
197 // Interpolate derivative metric at Lcoords
198 derx_1 = m_xmap->PhysEvaluate(I, DxD1);
199 derx_2 = m_xmap->PhysEvaluate(I, DxD2);
200 dery_1 = m_xmap->PhysEvaluate(I, DyD1);
201 dery_2 = m_xmap->PhysEvaluate(I, DyD2);
202
203 jac = dery_2 * derx_1 - dery_1 * derx_2;
204
205 // use analytical inverse of derivitives which are
206 // also similar to those of metric factors.
207 Lcoords[0] =
208 Lcoords[0] +
209 (dery_2 * (coords[0] - xmap) - derx_2 * (coords[1] - ymap)) / jac;
210
211 Lcoords[1] =
212 Lcoords[1] +
213 (-dery_1 * (coords[0] - xmap) + derx_1 * (coords[1] - ymap)) / jac;
214
215 if (!(std::isfinite(Lcoords[0]) && std::isfinite(Lcoords[1])))
216 {
217 dist = 1e16;
218 std::ostringstream ss;
219 ss << "nan or inf found in NewtonIterationForLocCoord in element "
220 << GetGlobalID();
221 WARNINGL1(false, ss.str());
222 return;
223 }
224 if (fabs(Lcoords[0]) > LcoordDiv || fabs(Lcoords[1]) > LcoordDiv)
225 {
226 break; // lcoords have diverged so stop iteration
227 }
228 }
229
230 m_xmap->LocCoordToLocCollapsed(Lcoords, eta);
231 if (ClampLocCoords(eta, 0.))
232 {
233 I[0] = m_xmap->GetBasis(0)->GetI(eta);
234 I[1] = m_xmap->GetBasis(1)->GetI(eta + 1);
235 // calculate the global point corresponding to Lcoords
236 xmap = m_xmap->PhysEvaluate(I, ptsx);
237 ymap = m_xmap->PhysEvaluate(I, ptsy);
238 F1 = coords[0] - xmap;
239 F2 = coords[1] - ymap;
240 dist = sqrt(F1 * F1 + F2 * F2);
241 }
242 else
243 {
244 dist = 0.;
245 }
246
247 if (cnt >= MaxIterations)
248 {
249 Array<OneD, NekDouble> collCoords(2);
250 m_xmap->LocCoordToLocCollapsed(Lcoords, collCoords);
251
252 // if coordinate is inside element dump error!
253 if ((collCoords[0] >= -1.0 && collCoords[0] <= 1.0) &&
254 (collCoords[1] >= -1.0 && collCoords[1] <= 1.0))
255 {
256 std::ostringstream ss;
257
258 ss << "Reached MaxIterations (" << MaxIterations
259 << ") in Newton iteration ";
260 ss << "Init value (" << setprecision(4) << init0 << "," << init1
261 << ","
262 << ") ";
263 ss << "Fin value (" << Lcoords[0] << "," << Lcoords[1] << ","
264 << ") ";
265 ss << "Resid = " << resid << " Tolerance = " << sqrt(ScaledTol);
266
267 WARNINGL1(cnt < MaxIterations, ss.str());
268 }
269 }
270}
271
273 Array<OneD, NekDouble> &Lcoords)
274{
275 NekDouble dist = std::numeric_limits<double>::max();
276 Array<OneD, NekDouble> tmpcoords(2);
277 tmpcoords[0] = coords[m_manifold[0]];
278 tmpcoords[1] = coords[m_manifold[1]];
279 if (GetMetricInfo()->GetGtype() == eRegular)
280 {
281 tmpcoords[0] -= m_isoParameter[0][0];
282 tmpcoords[1] -= m_isoParameter[1][0];
283 Lcoords[0] = m_invIsoParam[0][0] * tmpcoords[0] +
284 m_invIsoParam[0][1] * tmpcoords[1];
285 Lcoords[1] = m_invIsoParam[1][0] * tmpcoords[0] +
286 m_invIsoParam[1][1] * tmpcoords[1];
287 }
288 else if (m_straightEdge)
289 {
290 SolveStraightEdgeQuad(tmpcoords, Lcoords);
291 }
292 else if (GetMetricInfo()->GetGtype() == eDeformed)
293 {
294 v_FillGeom();
295 // Determine nearest point of coords to values in m_xmap
296 int npts = m_xmap->GetTotPoints();
297 Array<OneD, NekDouble> ptsx(npts), ptsy(npts);
298 Array<OneD, NekDouble> tmpx(npts), tmpy(npts);
299
300 // Determine 3D manifold orientation
301 m_xmap->BwdTrans(m_coeffs[m_manifold[0]], ptsx);
302 m_xmap->BwdTrans(m_coeffs[m_manifold[1]], ptsy);
303
304 Array<OneD, NekDouble> eta(2, 0.);
305 m_xmap->LocCoordToLocCollapsed(Lcoords, eta);
306 ClampLocCoords(eta, 0.);
307
308 m_xmap->LocCollapsedToLocCoord(eta, Lcoords);
309
310 // Perform newton iteration to find local coordinates
311 NewtonIterationForLocCoord(tmpcoords, ptsx, ptsy, Lcoords, dist);
312 }
313 return dist;
314}
315
317{
318 return m_verts.size();
319}
320
322{
323 return m_edges.size();
324}
325
327{
328 ASSERTL2(i >= 0 && i < m_verts.size(), "Index out of range");
329 return m_verts[i];
330}
331
333{
334 ASSERTL2(i >= 0 && i < m_edges.size(), "Index out of range");
335 return m_edges[i];
336}
337
339{
340 ASSERTL2(i >= 0 && i < m_eorient.size(), "Index out of range");
341 return m_eorient[i];
342}
343
345{
346 return 2;
347}
348
351{
352 if (m_geomFactors->GetGtype() == eRegular)
353 {
354 xiOut = Array<OneD, NekDouble>(2, 0.0);
355
356 GetLocCoords(xs, xiOut);
357 ClampLocCoords(xiOut);
358
359 Array<OneD, NekDouble> gloCoord(3);
360 gloCoord[0] = GetCoord(0, xiOut);
361 gloCoord[1] = GetCoord(1, xiOut);
362 gloCoord[2] = GetCoord(2, xiOut);
363
364 return sqrt((xs[0] - gloCoord[0]) * (xs[0] - gloCoord[0]) +
365 (xs[1] - gloCoord[1]) * (xs[1] - gloCoord[1]) +
366 (xs[2] - gloCoord[2]) * (xs[2] - gloCoord[2]));
367 }
368 // If deformed edge then the inverse mapping is non-linear so need to
369 // numerically solve for the local coordinate
370 else if (m_geomFactors->GetGtype() == eDeformed)
371 {
372 // Choose starting based on closest quad
373 Array<OneD, NekDouble> xi(2, 0.0), eta(2, 0.0);
374 m_xmap->LocCollapsedToLocCoord(eta, xi);
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 nq = m_xmap->GetTotPoints();
381
382 Array<OneD, NekDouble> x(nq), y(nq), z(nq);
383 m_xmap->BwdTrans(m_coeffs[0], x);
384 m_xmap->BwdTrans(m_coeffs[1], y);
385 m_xmap->BwdTrans(m_coeffs[2], z);
386
387 Array<OneD, NekDouble> xderxi1(nq, 0.0), yderxi1(nq, 0.0),
388 zderxi1(nq, 0.0), xderxi2(nq, 0.0), yderxi2(nq, 0.0),
389 zderxi2(nq, 0.0), xderxi1xi1(nq, 0.0), yderxi1xi1(nq, 0.0),
390 zderxi1xi1(nq, 0.0), xderxi1xi2(nq, 0.0), yderxi1xi2(nq, 0.0),
391 zderxi1xi2(nq, 0.0), xderxi2xi1(nq, 0.0), yderxi2xi1(nq, 0.0),
392 zderxi2xi1(nq, 0.0), xderxi2xi2(nq, 0.0), yderxi2xi2(nq, 0.0),
393 zderxi2xi2(nq, 0.0);
394
395 // Get first & second derivatives & partial derivatives of x,y,z values
396 std::array<NekDouble, 3> xc_derxi, yc_derxi, zc_derxi;
397
398 m_xmap->PhysDeriv(x, xderxi1, xderxi2);
399 m_xmap->PhysDeriv(y, yderxi1, yderxi2);
400 m_xmap->PhysDeriv(z, zderxi1, zderxi2);
401
402 m_xmap->PhysDeriv(xderxi1, xderxi1xi1, xderxi1xi2);
403 m_xmap->PhysDeriv(yderxi1, yderxi1xi1, yderxi1xi2);
404 m_xmap->PhysDeriv(zderxi1, zderxi1xi1, zderxi1xi2);
405
406 m_xmap->PhysDeriv(yderxi2, yderxi2xi1, yderxi2xi2);
407 m_xmap->PhysDeriv(xderxi2, xderxi2xi1, xderxi2xi2);
408 m_xmap->PhysDeriv(zderxi2, zderxi2xi1, zderxi2xi2);
409
410 // Minimisation loop (Quasi-newton method)
411 NekDouble fx_prev = std::numeric_limits<NekDouble>::max();
412 for (int i = 0; i < NekConstants::kNewtonIterations; ++i)
413 {
414 // Compute the objective function, f(x_k) and its derivatives
415 NekDouble xc = m_xmap->PhysEvaluate(xi, x, xc_derxi);
416 NekDouble yc = m_xmap->PhysEvaluate(xi, y, yc_derxi);
417 NekDouble zc = m_xmap->PhysEvaluate(xi, z, zc_derxi);
418
419 NekDouble xc_derxi1xi1 = m_xmap->PhysEvaluate(xi, xderxi1xi1);
420 NekDouble yc_derxi1xi1 = m_xmap->PhysEvaluate(xi, yderxi1xi1);
421 NekDouble zc_derxi1xi1 = m_xmap->PhysEvaluate(xi, zderxi1xi1);
422
423 NekDouble xc_derxi1xi2 = m_xmap->PhysEvaluate(xi, xderxi1xi2);
424 NekDouble yc_derxi1xi2 = m_xmap->PhysEvaluate(xi, yderxi1xi2);
425 NekDouble zc_derxi1xi2 = m_xmap->PhysEvaluate(xi, zderxi1xi2);
426
427 NekDouble xc_derxi2xi2 = m_xmap->PhysEvaluate(xi, xderxi2xi2);
428 NekDouble yc_derxi2xi2 = m_xmap->PhysEvaluate(xi, yderxi2xi2);
429 NekDouble zc_derxi2xi2 = m_xmap->PhysEvaluate(xi, zderxi2xi2);
430
431 // Objective function is the distance to the search point
432 NekDouble xdiff = xc - xs[0];
433 NekDouble ydiff = yc - xs[1];
434 NekDouble zdiff = zc - xs[2];
435
436 NekDouble fx = xdiff * xdiff + ydiff * ydiff + zdiff * zdiff;
437
438 NekDouble fx_derxi1 = 2.0 * xdiff * xc_derxi[0] +
439 2.0 * ydiff * yc_derxi[0] +
440 2.0 * zdiff * zc_derxi[0];
441
442 NekDouble fx_derxi2 = 2.0 * xdiff * xc_derxi[1] +
443 2.0 * ydiff * yc_derxi[1] +
444 2.0 * zdiff * zc_derxi[1];
445
446 NekDouble fx_derxi1xi1 =
447 2.0 * xdiff * xc_derxi1xi1 + 2.0 * xc_derxi[0] * xc_derxi[0] +
448 2.0 * ydiff * yc_derxi1xi1 + 2.0 * yc_derxi[0] * yc_derxi[0] +
449 2.0 * zdiff * zc_derxi1xi1 + 2.0 * zc_derxi[0] * zc_derxi[0];
450
451 NekDouble fx_derxi1xi2 =
452 2.0 * xdiff * xc_derxi1xi2 + 2.0 * xc_derxi[1] * xc_derxi[0] +
453 2.0 * ydiff * yc_derxi1xi2 + 2.0 * yc_derxi[1] * yc_derxi[0] +
454 2.0 * zdiff * zc_derxi1xi2 + 2.0 * zc_derxi[1] * zc_derxi[0];
455
456 NekDouble fx_derxi2xi2 =
457 2.0 * xdiff * xc_derxi2xi2 + 2.0 * xc_derxi[1] * xc_derxi[1] +
458 2.0 * ydiff * yc_derxi2xi2 + 2.0 * yc_derxi[1] * yc_derxi[1] +
459 2.0 * zdiff * zc_derxi2xi2 + 2.0 * zc_derxi[1] * zc_derxi[1];
460
461 // Jacobian
462 NekDouble jac[2];
463 jac[0] = fx_derxi1;
464 jac[1] = fx_derxi2;
465
466 // Inverse of 2x2 hessian
467 NekDouble hessInv[2][2];
468
469 NekDouble det =
470 1 / (fx_derxi1xi1 * fx_derxi2xi2 - fx_derxi1xi2 * fx_derxi1xi2);
471 hessInv[0][0] = det * fx_derxi2xi2;
472 hessInv[0][1] = det * -fx_derxi1xi2;
473 hessInv[1][0] = det * -fx_derxi1xi2;
474 hessInv[1][1] = det * fx_derxi1xi1;
475
476 // Check for convergence
477 if (abs(fx - fx_prev) < 1e-12)
478 {
479 fx_prev = fx;
480 break;
481 }
482 else
483 {
484 fx_prev = fx;
485 }
486
487 NekDouble gamma = 1.0;
488 bool conv = false;
489
490 // Search direction: Newton's method
491 NekDouble pk[2];
492 pk[0] = -(hessInv[0][0] * jac[0] + hessInv[1][0] * jac[1]);
493 pk[1] = -(hessInv[0][1] * jac[0] + hessInv[1][1] * jac[1]);
494
495 // Backtracking line search
496 while (gamma > 1e-10)
497 {
498 Array<OneD, NekDouble> xi_pk(2);
499 xi_pk[0] = xi[0] + pk[0] * gamma;
500 xi_pk[1] = xi[1] + pk[1] * gamma;
501
502 Array<OneD, NekDouble> eta_pk(2, 0.0);
503 m_xmap->LocCoordToLocCollapsed(xi_pk, eta_pk);
504
505 if (eta_pk[0] <
506 (-1 - std::numeric_limits<NekDouble>::epsilon()) ||
507 eta_pk[0] >
508 (1 + std::numeric_limits<NekDouble>::epsilon()) ||
509 eta_pk[1] <
510 (-1 - std::numeric_limits<NekDouble>::epsilon()) ||
511 eta_pk[1] > (1 + std::numeric_limits<NekDouble>::epsilon()))
512 {
513 gamma /= 2.0;
514 continue;
515 }
516
517 std::array<NekDouble, 3> xc_pk_derxi, yc_pk_derxi, zc_pk_derxi;
518
519 NekDouble xc_pk = m_xmap->PhysEvaluate(xi_pk, x, xc_pk_derxi);
520 NekDouble yc_pk = m_xmap->PhysEvaluate(xi_pk, y, yc_pk_derxi);
521 NekDouble zc_pk = m_xmap->PhysEvaluate(xi_pk, z, zc_pk_derxi);
522
523 NekDouble xc_pk_diff = xc_pk - xs[0];
524 NekDouble yc_pk_diff = yc_pk - xs[1];
525 NekDouble zc_pk_diff = zc_pk - xs[2];
526
527 NekDouble fx_pk = xc_pk_diff * xc_pk_diff +
528 yc_pk_diff * yc_pk_diff +
529 zc_pk_diff * zc_pk_diff;
530
531 NekDouble fx_pk_derxi1 = 2.0 * xc_pk_diff * xc_pk_derxi[0] +
532 2.0 * yc_pk_diff * yc_pk_derxi[0] +
533 2.0 * zc_pk_diff * zc_pk_derxi[0];
534
535 NekDouble fx_pk_derxi2 = 2.0 * xc_pk_diff * xc_pk_derxi[1] +
536 2.0 * yc_pk_diff * yc_pk_derxi[1] +
537 2.0 * zc_pk_diff * zc_pk_derxi[1];
538
539 // Check Wolfe conditions using Armijo constants
540 // https://en.wikipedia.org/wiki/Wolfe_conditions
541 NekDouble tmp = pk[0] * fx_derxi1 + pk[1] * fx_derxi2;
542 NekDouble tmp2 = pk[0] * fx_pk_derxi1 + pk[1] * fx_pk_derxi2;
543 if ((fx_pk - (fx + c1 * gamma * tmp)) <
544 std::numeric_limits<NekDouble>::epsilon() &&
545 (-tmp2 - (-c2 * tmp)) <
546 std::numeric_limits<NekDouble>::epsilon())
547 {
548 conv = true;
549 break;
550 }
551
552 gamma /= 2.0;
553 }
554
555 if (!conv)
556 {
557 break;
558 }
559
560 xi[0] += gamma * pk[0];
561 xi[1] += gamma * pk[1];
562 }
563
564 xiOut = xi;
565 return sqrt(fx_prev);
566 }
567 else
568 {
569 ASSERTL0(false, "Geometry type unknown")
570 }
571
572 return -1.0;
573}
574
576{
577 NekDouble Jac = m_isoParameter[0][1] * m_isoParameter[1][2] -
578 m_isoParameter[1][1] * m_isoParameter[0][2];
579 Jac = 1. / Jac;
580 // a12, -a02, -a11, a01
584 m_invIsoParam[0][0] = m_isoParameter[1][2] * Jac;
585 m_invIsoParam[0][1] = -m_isoParameter[0][2] * Jac;
586 m_invIsoParam[1][0] = -m_isoParameter[1][1] * Jac;
587 m_invIsoParam[1][1] = m_isoParameter[0][1] * Jac;
588}
589
590} // namespace Nektar::SpatialDomains
#define WARNINGL1(condition, msg)
Definition: ErrorUtil.hpp:243
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed to...
Definition: ErrorUtil.hpp:265
int v_GetShapeDim() const override
Get the object's shape dimension.
Definition: Geometry2D.cpp:344
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:136
Array< OneD, int > m_manifold
Definition: Geometry2D.h:83
Array< OneD, Array< OneD, NekDouble > > m_edgeNormal
Definition: Geometry2D.h:84
void SolveStraightEdgeQuad(const Array< OneD, const NekDouble > &coords, Array< OneD, NekDouble > &Lcoords)
Definition: Geometry2D.cpp:109
PointGeomSharedPtr v_GetVertex(int i) const override
Definition: Geometry2D.cpp:326
Geometry1DSharedPtr v_GetEdge(int i) const override
Returns edge i of this object.
Definition: Geometry2D.cpp:332
std::vector< StdRegions::Orientation > m_eorient
Definition: Geometry2D.h:81
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 ge...
Definition: Geometry2D.cpp:272
int v_AllLeftCheck(const Array< OneD, const NekDouble > &gloCoord) override
Definition: Geometry2D.cpp:63
void v_CalculateInverseIsoParam() override
Definition: Geometry2D.cpp:575
StdRegions::Orientation v_GetEorient(const int i) const override
Returns the orientation of edge i with respect to the ordering of edges in the standard element.
Definition: Geometry2D.cpp:338
NekDouble v_FindDistance(const Array< OneD, const NekDouble > &xs, Array< OneD, NekDouble > &xi) override
Definition: Geometry2D.cpp:349
int v_GetNumEdges() const override
Get the number of edges of this object.
Definition: Geometry2D.cpp:321
int v_GetNumVerts() const override
Get the number of vertices of this object.
Definition: Geometry2D.cpp:316
Base class for shape geometry information.
Definition: Geometry.h:79
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:552
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:542
virtual void v_FillGeom()
Populate the coordinate mapping Geometry::m_coeffs information from any children geometry elements.
Definition: Geometry.cpp:375
Array< OneD, Array< OneD, NekDouble > > m_isoParameter
Definition: Geometry.h:204
Array< OneD, Array< OneD, NekDouble > > m_invIsoParam
Definition: Geometry.h:205
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:550
int GetGlobalID(void) const
Get the ID of this object.
Definition: Geometry.h:324
Array< OneD, Array< OneD, NekDouble > > m_coeffs
Array containing expansion coefficients of m_xmap.
Definition: Geometry.h:201
GeomFactorsSharedPtr GetMetricInfo()
Get the geometric factors for this object.
Definition: Geometry.h:308
StdRegions::StdExpansionSharedPtr m_xmap
mapping containing isoparametric transformation.
Definition: Geometry.h:189
StdRegions::StdExpansionSharedPtr GetXmap() const
Return the mapping object Geometry::m_xmap that represents the coordinate transformation from standar...
Definition: Geometry.h:433
GeomFactorsSharedPtr m_geomFactors
Geometric factors.
Definition: Geometry.h:185
int m_coordim
Coordinate dimension of this geometry object.
Definition: Geometry.h:183
@ beta
Gauss Radau pinned at x=-1,.
Definition: PointsType.h:59
static const unsigned int kNewtonIterations
std::shared_ptr< Curve > CurveSharedPtr
Definition: Curve.hpp:58
@ eRegular
Geometry is straight-sided with constant geometric factors.
@ eDeformed
Geometry is curved or has non-constant factors.
std::shared_ptr< PointGeom > PointGeomSharedPtr
Definition: Geometry.h:57
std::shared_ptr< Geometry1D > Geometry1DSharedPtr
Definition: Geometry.h:61
std::vector< double > z(NPUPPER)
double NekDouble
T Vsum(int n, const T *x, const int incx)
Subtract return sum(x)
Definition: Vmath.hpp:608
scalarT< T > abs(scalarT< T > in)
Definition: scalar.hpp:298
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294