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