Nektar++
ProcessWallNormalData.cpp
Go to the documentation of this file.
1////////////////////////////////////////////////////////////////////////////////
2//
3// File: ProcessWallNormalData.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: Get the wall-normal data at a given origin.
32//
33////////////////////////////////////////////////////////////////////////////////
34
35#include <iostream>
36#include <string>
37
42
44
45using namespace std;
46
47namespace Nektar
48{
49namespace FieldUtils
50{
51
54 ModuleKey(eProcessModule, "wallNormalData"),
56 "Export data in wall-normal direction from a single point on the "
57 "wall.");
58
61{
62 f->m_writeBndFld = false; // turned on in upstream ProcessBoundaryExtract
63
64 m_config["xorig"] =
65 ConfigOption(false, "0.5,0.0,0.0",
66 "The point to be projected onto the wall to get the \
67 sampling origin. default=[0.5,0,0]");
68 m_config["projDir"] =
69 ConfigOption(false, "0.0,1.0,0.0",
70 "The direction to project the point onto the wall to \
71 get the sampling origin. default=[0,1,0]");
72 m_config["maxDist"] = ConfigOption(
73 false, "1.0", "Distance to limit projection distance to find the \
74 desired sampling origin. defalut=1.0");
75 m_config["distH"] = ConfigOption(
76 false, "0.01", "Sampling distance along the wall normal at the \
77 sampling origin. default=0.1");
78 m_config["nptsH"] = ConfigOption(
79 false, "5", "Number of sampling points along the wall normal. \
80 default=5");
82 false, "0.1", "The parameter that controls the sampling points' \
83 distribution. d should be in the range (0,inf). d \
84 in (0,0.95] gives controled points; d in (0.95,inf) \
85 gives evenly spaced points");
86
87 // To allow some functions called somewhere else outside Process()
88 m_spacedim = m_f->m_exp[0]->GetCoordim(0) + m_f->m_numHomogeneousDir;
89}
90
92{
93}
94
95/*Note
96 * This module is used to get field data in the wall-normal direction.
97 * The input cases can be 2D, 2.5D and 3D.
98 * The data will be exported with .pts extension.
99 *
100 * The user defined parameters are: bnd, xorig, searchDir, maxDist, distH,
101 * nptsH, and d.
102 * - bnd=0 is the boundary id. This boundary should contain the desired origin.
103 * - xorig="x,y,z" are the coordinates of the input sampling origin. They are
104 * used as the references to get exact origin on the wall.
105 * - projDir="0,1,0" is the projection direction to find the point on the wall
106 * - maxDist=1.0 is the distance that constrains the projection on the wall.
107 * The projection distance should be smaller than this distance. This
108 * parameter is set in case the boundary is wavey so that multiple
109 * projections exist at the same time.
110 * - distH=0.01 is the sampling depth in the wall-normal direction.
111 * - nptsH=11 is the number of sampling points along wall-normal direction.
112 * - d=0.1 is a destribution control parameter of the sampling points. It
113 * should be in the range (0,inf). d in range (0,0.95] for controlled array.
114 * d>0.95 gives evenly spaced array.
115 */
116void ProcessWallNormalData::v_Process(po::variables_map &vm)
117{
119
120 // Get dim to store data
121 const int nfields = m_f->m_variables.size();
122 const int nCoordDim = m_f->m_exp[0]->GetCoordim(0);
123 const int nBndLcoordDim = nCoordDim - 1;
124 const int totVars = m_spacedim + nfields;
125
126 // Initialize the sampling parameters
127 std::vector<NekDouble> xorig, searchDir;
128 ASSERTL0(ParseUtils::GenerateVector(m_config["xorig"].as<string>(), xorig),
129 "Failed to interpret origin coordinates string");
130 ASSERTL0(
131 ParseUtils::GenerateVector(m_config["projDir"].as<string>(), searchDir),
132 "Failed to interpret search direction string");
133 const NekDouble maxDist = m_config["maxDist"].as<NekDouble>();
134 const NekDouble distH = m_config["distH"].as<NekDouble>();
135 const int nptsH = m_config["nptsH"].as<int>();
136 const NekDouble delta = m_config["d"].as<NekDouble>();
137
138 Array<OneD, NekDouble> orig(3), projDir(3); // gloCoord of the origin
139 for (int i = 0; i < 3; ++i)
140 {
141 orig[i] = (xorig.size() >= (i + 1)) ? xorig[i] : 0.0;
142 projDir[i] = (searchDir.size() >= (i + 1)) ? searchDir[i] : 0.0;
143 }
144 if (nCoordDim == 2)
145 {
146 projDir[2] = 0.0;
147 }
148 Vmath::Smul(3, 1.0 / sqrt(Vmath::Dot(3, projDir, 1, projDir, 1)), projDir,
149 1, projDir, 1);
150
151 // Update z according to the closest plane in 2.5D cases
152 if (m_f->m_numHomogeneousDir == 1)
153 {
154 int nPlanes = m_f->m_exp[0]->GetHomogeneousBasis()->GetZ().size();
155 NekDouble lHom = m_f->m_exp[0]->GetHomoLen();
156 if (orig[2] < 0.0 || orig[2] > lHom)
157 {
158 orig[2] = 0.0;
159 }
160 else
161 {
162 NekDouble dZ = lHom / nPlanes;
163 NekDouble zTmp, zCur = 0.0, distTmp, distCur = 999.0;
164 for (int i = 0; i <= nPlanes; ++i)
165 {
166 zTmp = dZ * i;
167 distTmp = fabs(orig[2] - zTmp);
168 if (distTmp < distCur)
169 {
170 distCur = distTmp;
171 zCur = zTmp;
172 if (i == nPlanes)
173 {
174 zCur = 0.0;
175 }
176 }
177 }
178 orig[2] = zCur;
179 }
180 }
181
182 // Get the bnd id
184 m_f->m_exp[0]->GetGraph());
186 bcs.GetBoundaryRegions();
187 map<int, int> BndRegionMap;
188 int cnt = 0;
189 for (auto &breg_it : bregions)
190 {
191 BndRegionMap[breg_it.first] = cnt++;
192 }
193 int bnd = BndRegionMap[m_f->m_bndRegionsToWrite[0]];
194
195 // Get expansion list for boundary
197 for (int i = 0; i < nfields; ++i)
198 {
199 BndExp[i] = m_f->m_exp[i]->UpdateBndCondExpansion(bnd);
200 }
201
202 //-------------------------------------------------------------------------
203 // Find the element that contains the origin, and give the precise origin
205 Array<OneD, NekDouble> locCoord(nBndLcoordDim, -999.0);
206 NekDouble projDist;
207 int elmtid;
208 bool isInside = false;
209
210 // Search and get precise locCoord
211 const int nElmts = BndExp[0]->GetNumElmts();
212 for (elmtid = 0; elmtid < nElmts; ++elmtid) // nElmts
213 {
214 bndGeom = BndExp[0]->GetExp(elmtid)->GetGeom();
215 isInside = BndElmtContainsPoint(bndGeom, orig, projDir, locCoord,
216 projDist, maxDist);
217 if (isInside)
218 {
219 break;
220 }
221 }
222 ASSERTL0(isInside, "Failed to find the sampling origin on the boundary.");
223
224 // Then Update the precise sampling position
225 StdRegions::StdExpansionSharedPtr bndXmap = bndGeom->GetXmap();
226 const int npts = bndXmap->GetTotPoints();
228 Array<OneD, Array<OneD, const NekDouble>> bndCoeffs(nCoordDim);
229
230 for (int i = 0; i < nCoordDim; ++i)
231 {
232 pts[i] = Array<OneD, NekDouble>(npts);
233 bndCoeffs[i] = bndGeom->GetCoeffs(i); // 0/1/2 for x/y/z
234 bndXmap->BwdTrans(bndCoeffs[i], pts[i]);
235 orig[i] = bndXmap->PhysEvaluate(locCoord, pts[i]); // Update wall point
236 }
237
238 // Get outward-pointing normal vectors for all quadrature points on bnd
239 // Use these normals as diretion references
241 m_f->m_exp[0]->GetBoundaryNormals(bnd, normalsQ);
242
243 // Get the normals in the element
244 // Set point key to get point id
245 Array<OneD, LibUtilities::PointsKey> from_key(nBndLcoordDim);
246 int from_nPtsPerElmt = 1;
247 for (int i = 0; i < nBndLcoordDim; ++i)
248 {
249 from_key[i] = BndExp[0]->GetExp(elmtid)->GetBasis(i)->GetPointsKey();
250 from_nPtsPerElmt *= from_key[i].GetNumPoints();
251 }
252
253 // Get ref normal
254 Array<OneD, NekDouble> normalsRef(3, 0.0);
255 int refId = elmtid * from_nPtsPerElmt + 0; // the 1st normal in the bnd elmt
256 for (int i = 0; i < m_spacedim; ++i)
257 {
258 normalsRef[i] = normalsQ[i][refId];
259 }
260
261 // Get the precise normals and correct the direction to be inward
262 Array<OneD, NekDouble> normals(3, 0.0);
263 GetNormals(bndGeom, locCoord, normals);
264
265 if (Vmath::Dot(3, normals, normalsRef) > 0.0)
266 {
267 Vmath::Neg(3, normals, 1);
268 }
269
270 // Output the info if -v is set
271 if (m_f->m_verbose)
272 {
273 cout << "------ wallNormalData module ------\n";
274 cout << "Input point:\n";
275 cout << " - [Px,Py,Pz] = [" << xorig[0] << ", " << xorig[1];
276 if (xorig.size() >= 3)
277 {
278 cout << ", " << xorig[2] << "]\n";
279 }
280 else
281 {
282 cout << ", " << 0.0 << "]\n";
283 }
284 cout << "Projection direction:\n";
285 cout << " - [vx,vy,vz] = [" << projDir[0] << ", " << projDir[1] << ", "
286 << projDir[2] << "]\n";
287 cout << "Sampling origin on the wall:\n";
288 cout << " - [Ox,Oy,Oz] = [" << orig[0] << ", " << orig[1] << ", "
289 << orig[2] << "]\n";
290 cout << "Normals at the origin:\n";
291 cout << " - [nx,ny,nz] = [" << normals[0] << ", " << normals[1] << ", "
292 << normals[2] << "]\n";
293 cout
294 << "Ref normals (at quadrature points in the projected element):\n";
295 for (int i = 0; i < from_nPtsPerElmt; ++i)
296 {
297 cout << " - " << i << ": [nx,ny,nz] = ["
298 << -normalsQ[0][elmtid * from_nPtsPerElmt + i] << ", "
299 << -normalsQ[1][elmtid * from_nPtsPerElmt + i];
300 if (m_spacedim == 3)
301 {
302 cout << ", " << -normalsQ[2][elmtid * from_nPtsPerElmt + i]
303 << "]\n";
304 }
305 else
306 {
307 cout << "]\n";
308 }
309 }
310 cout << endl;
311 }
312
313 //-------------------------------------------------------------------------
314 // Set the depth of sampling array
315 Array<OneD, NekDouble> h(nptsH, 0.0);
316 if (delta > 0 && delta <= 0.95)
317 {
318 // Use expression in Agrawal's paper:
319 // h = 1-tanh((1-ksi)*atanh(sqrt(1-delta)))/sqrt(1-delta), ksi in [0,1]
320 NekDouble tmp1;
321 const NekDouble tmp2 = 1.0 / (static_cast<NekDouble>(nptsH) - 1.0);
322 const NekDouble tmp3 = sqrt(1.0 - delta);
323 const NekDouble tmp4 = atanh(tmp3);
324 const NekDouble tmp5 = 1.0 / tmp3;
325 for (int i = 1; i < nptsH; ++i)
326 {
327 tmp1 = 1.0 - i * tmp2; // tmp1 = 1-ksi
328 h[i] = 1 - tanh(tmp1 * tmp4) * tmp5;
329 h[i] *= distH; // physical distance in normal direction
330 }
331 }
332 else if (delta > 0.95)
333 {
334 // Evenly spaced array
335 const NekDouble tmp1 = 1.0 / (static_cast<NekDouble>(nptsH) - 1.0);
336 for (int i = 1; i < nptsH; ++i)
337 {
338 h[i] = i * tmp1;
339 h[i] *= distH; // physical distance in normal direction
340 }
341 }
342 else
343 {
344 ASSERTL0(false, "Input error. Delta needs to be greater than 0.0.");
345 }
346
347 // Set pts coordinates and interpoate the data
349 for (int i = 0; i < totVars; ++i)
350 {
351 ptsH[i] = Array<OneD, NekDouble>(nptsH, 0.0);
352 }
353
354 for (int i = 0; i < m_spacedim; ++i)
355 {
356 for (int j = 0; j < nptsH; ++j)
357 {
358 ptsH[i][j] = orig[i] + h[j] * normals[i]; // x0+dist*nx
359 }
360 }
361
363 m_spacedim, m_f->m_variables, ptsH, LibUtilities::NullPtsInfoMap);
364
366 interp.Interpolate(m_f->m_exp, m_f->m_fieldPts,
368}
369
370/**
371 * @brief Project a single point along the given direction to a plane
372 * @param gloCoord Global coordinate of the point. size=3.
373 * @param projDir Projection direction, which is also the normal vector of
374 * the target plane. size=3, norm=1.
375 * @param distToOrig The distance from the origin (0,0,0) to the target plane
376 * @param projGloCoord The global coordinate of the projecion result.
377 */
379 const Array<OneD, const NekDouble> &gloCoord,
380 const Array<OneD, const NekDouble> &projDir, const NekDouble distToOrig,
381 Array<OneD, NekDouble> &projGloCoord)
382{
383 NekDouble tmp1;
385
386 tmp1 = Vmath::Dot(3, gloCoord, 1, projDir, 1); // |xn| = x dot n
387 Vmath::Smul(3, distToOrig - tmp1, projDir, 1, tmp2,
388 1); // d_xn = (dist-|xn|) n
389 Vmath::Vadd(3, gloCoord, 1, tmp2, 1, projGloCoord, 1); // x' = x + d_xn
390}
391
392/**
393 * @brief Project the vertices of the elmt along the given direction to a plane
394 * @param pts Global coordinate of the vertices of the elmt. size=2/3.
395 * @param projDir Projection direction, which is also the normal vector of
396 * the target plane. size=3, norm=1.
397 * @param distToOrig The distance from the origin (0,0,0) to the target plane.
398 * @param projPts The global coordinate of the projecion result.
399 */
401 const Array<OneD, const Array<OneD, NekDouble>> &pts,
402 const Array<OneD, const NekDouble> &projDir, const NekDouble distToOrig,
404{
405 const int nCoordDim = pts.size(); // 2 for 2.5D cases
406 const int npts = pts[0].size();
407
408 NekDouble tmp1;
409 Array<OneD, NekDouble> singlePnt(nCoordDim), tmp2(nCoordDim);
410
411 for (int i = 0; i < npts; ++i)
412 {
413 // Get a point
414 for (int j = 0; j < nCoordDim; ++j)
415 {
416 singlePnt[j] = pts[j][i];
417 }
418
419 // Projection
420 tmp1 = Vmath::Dot(nCoordDim, singlePnt, 1, projDir, 1);
421 Vmath::Smul(nCoordDim, distToOrig - tmp1, projDir, 1, tmp2, 1);
422 Vmath::Vadd(nCoordDim, singlePnt, 1, tmp2, 1, singlePnt, 1);
423
424 // Save to projPts
425 for (int j = 0; j < nCoordDim; ++j)
426 {
427 projPts[j][i] = singlePnt[j];
428 }
429 }
430}
431
432/**
433 * @brief Determine if the projected point is inside the projected element.
434 * @param projGloCoord The global coordinate of the projected single point.
435 * @param projPts The global coordinate of the projected vertices,size=2/3
436 * @param projDir Projection direction, which is used as the reference
437 * direction in the 3D routine. size=3, norm=1.
438 * @param paralTol Tolerence to check if two vectors are parallel.
439 * @param angleTol Tolerence to check if the total angle is 2*pi.
440 * @return Inside (true) or not (false)
441 */
443 const Array<OneD, const NekDouble> &projGloCoord,
444 const Array<OneD, const Array<OneD, NekDouble>> &projPts,
445 const NekDouble paralTol)
446{
447 const NekDouble npts = projPts[0].size();
448
449 Array<OneD, NekDouble> vec1(2, 0.0), vec2(2, 0.0);
450
451 vec1[0] = projPts[0][0] - projGloCoord[0];
452 vec1[1] = projPts[1][0] - projGloCoord[1];
453 vec2[0] = projPts[0][npts - 1] - projGloCoord[0];
454 vec2[1] = projPts[1][npts - 1] - projGloCoord[1];
455
456 Vmath::Smul(2, 1.0 / sqrt(Vmath::Dot(2, vec1, 1, vec1, 1)), vec1, 1, vec1,
457 1);
458 Vmath::Smul(2, 1.0 / sqrt(Vmath::Dot(2, vec2, 1, vec2, 1)), vec2, 1, vec2,
459 1);
460
461 if ((fabs(vec1[0] + vec2[0]) + fabs(vec1[1] + vec2[1])) < paralTol)
462 {
463 // On the line segement, true
464 return true;
465 }
466 else
467 {
468 return false;
469 }
470}
471
473 const Array<OneD, const NekDouble> &projGloCoord,
474 const Array<OneD, const Array<OneD, NekDouble>> &projPts,
475 const Array<OneD, const NekDouble> &projDir, const NekDouble paralTol,
476 const NekDouble angleTol)
477{
478
479 // Generate a polygen (quad), re-order the points on the edge
480 const NekDouble npts = projPts[0].size();
481 const int nptsEdge =
482 sqrt(npts); // num of points on edge for geo representation
483 const int nptsPolygon = 4 * nptsEdge - 4;
485 for (int i = 0; i < 3; ++i)
486 {
487 ptsPolygon[i] = Array<OneD, NekDouble>(nptsPolygon);
488
489 for (int j = 0; j < nptsEdge; ++j)
490 {
491 ptsPolygon[i][j] = projPts[i][j];
492 }
493 for (int j = 0; j < nptsEdge - 2; ++j)
494 {
495 ptsPolygon[i][nptsEdge + j] = projPts[i][(j + 2) * nptsEdge - 1];
496 }
497 for (int j = 0; j < nptsEdge; ++j)
498 {
499 ptsPolygon[i][2 * nptsEdge - 2 + j] = projPts[i][npts - 1 - j];
500 }
501 for (int j = 0; j < nptsEdge - 2; ++j)
502 {
503 ptsPolygon[i][3 * nptsEdge - 2 + j] =
504 projPts[i][nptsEdge * (nptsEdge - j - 2)];
505 }
506 }
507
508 // Determine relation using angle method (sum=2*pi)
509 NekDouble angleCos, angleAbs, angleSign, angleSum = 0.0;
510 Array<OneD, NekDouble> vec1(3, 0.0), vec2(3, 0.0), vec3(3, 0.0);
511 int id1, id2;
512
513 for (int i = 0; i < nptsPolygon; ++i)
514 {
515 id1 = i;
516 id2 = (id1 == (nptsPolygon - 1)) ? 0 : (id1 + 1);
517
518 for (int j = 0; j < 3; ++j)
519 {
520 vec1[j] = ptsPolygon[j][id1] - projGloCoord[j];
521 vec2[j] = ptsPolygon[j][id2] - projGloCoord[j];
522 }
523
524 Vmath::Smul(3, 1.0 / sqrt(Vmath::Dot(3, vec1, 1, vec1, 1)), vec1, 1,
525 vec1, 1);
526 Vmath::Smul(3, 1.0 / sqrt(Vmath::Dot(3, vec2, 1, vec2, 1)), vec2, 1,
527 vec2, 1);
528
529 if ((fabs(vec1[0] + vec2[0]) + fabs(vec1[1] + vec2[1]) +
530 fabs(vec1[2] + vec2[2])) < paralTol)
531 {
532 // On the line segement, true
533 // This branch is used to aviod angle=pi but angleSign=-1 (not 1)
534 return true;
535 }
536 else
537 {
538 // Off the line segement, compute angle
539 // vec3 = vec1 x vec2, not scaled
540 vec3[0] = vec1[1] * vec2[2] - vec1[2] * vec2[1];
541 vec3[1] = vec1[2] * vec2[0] - vec1[0] * vec2[2];
542 vec3[2] = vec1[0] * vec2[1] - vec1[1] * vec2[0];
543
544 // Use the projDir as reference direction (positive)
545 angleSign = Vmath::Dot(3, vec3, 1, projDir, 1) > 0.0 ? 1.0 : -1.0;
546 angleCos = Vmath::Dot(3, vec1, 1, vec2, 1);
547
548 // Avoid error in using acos()
549 if (angleCos > 1.0)
550 {
551 angleCos = 1.0;
552 }
553 else if (angleCos < -1.0)
554 {
555 angleCos = -1.0;
556 }
557
558 angleAbs = acos(angleCos);
559 angleSum += angleSign * angleAbs;
560 }
561 }
562
563 // Check
564 angleSum = fabs(angleSum);
565 if (fabs(angleSum - 2.0 * M_PI) < angleTol)
566 {
567 return true;
568 }
569 else
570 {
571 return false;
572 }
573}
574
575/**
576 * @brief Use iteration to get the locCoord. This routine should be used after
577 * we have checked the projected point is inside the projected element.
578 * @param bndGeom Geometry to get the xmap.
579 * @param gloCoord Global coordinate of the point. size=3.
580 * @param pts Global coordinate of the vertices of the elmt. size=2/3.
581 * @param dieUse The main direction(s) used to compute local coordinate
582 * @param locCoord Iteration results for local coordinate(s)
583 * @param dist Returned distance in physical space if the collapsed
584 * locCoord is out of range [-1,1].
585 * @param iterTol Tolerence for iteration.
586 * @param iterMax Maximum iteration steps
587 * @return Converged (true) or not (false)
588 */
591 const Array<OneD, const NekDouble> &gloCoord,
592 const Array<OneD, const Array<OneD, NekDouble>> &pts,
593 const Array<OneD, const int> &dirUse, Array<OneD, NekDouble> &locCoord,
594 const NekDouble iterTol, const int iterMax)
595{
596 // Initial settings
597 Array<OneD, NekDouble> etaLR(2); // range [-1,1]
598 etaLR[0] = -1.0; // left
599 etaLR[1] = 1.0; // right
600 NekDouble tmpL, tmpR; // tmp values for L/R
601
602 StdRegions::StdExpansionSharedPtr bndXmap = bndGeom->GetXmap();
603 locCoord[0] = -2.0;
604 int cnt = 0;
605 bool isConverge = false;
606 while (cnt < iterMax)
607 {
608 tmpL = bndXmap->PhysEvaluate(etaLR, pts[dirUse[0]]);
609 tmpR = bndXmap->PhysEvaluate(etaLR + 1, pts[dirUse[0]]);
610
611 if (fabs(gloCoord[dirUse[0]] - tmpL) >=
612 fabs(gloCoord[dirUse[0]] - tmpR))
613 {
614 etaLR[0] = 0.5 * (etaLR[0] + etaLR[1]);
615 }
616 else
617 {
618 etaLR[1] = 0.5 * (etaLR[0] + etaLR[1]);
619 }
620
621 if ((etaLR[1] - etaLR[0]) < iterTol)
622 {
623 locCoord[0] = 0.5 * (etaLR[0] + etaLR[1]);
624 isConverge = true;
625 break;
626 }
627
628 ++cnt;
629 }
630
631 // Warning if failed
632 if (cnt >= iterMax)
633 {
634 WARNINGL1(false, "Bisection iteration is not converged");
635 }
636
637 return isConverge;
638}
639
642 const Array<OneD, const NekDouble> &gloCoord,
643 const Array<OneD, const Array<OneD, NekDouble>> &pts,
644 const Array<OneD, const int> &dirUse, Array<OneD, NekDouble> &locCoord,
645 NekDouble &dist, const NekDouble iterTol, const int iterMax)
646{
647
648 const NekDouble LcoordDiv = 15.0;
649
650 StdRegions::StdExpansionSharedPtr bndXmap = bndGeom->GetXmap();
651
653 bndGeom->GetMetricInfo()->GetJac(bndXmap->GetPointsKeys());
654 NekDouble scaledTol =
655 Vmath::Vsum(Jac.size(), Jac, 1) / ((NekDouble)Jac.size());
656 scaledTol *= iterTol;
657
658 // Set the gloCoord used to compute locCoord
659 const int dir1 = dirUse[0];
660 const int dir2 = dirUse[1];
661
662 Array<OneD, NekDouble> Dx1D1(pts[dir1].size());
663 Array<OneD, NekDouble> Dx1D2(pts[dir1].size());
664 Array<OneD, NekDouble> Dx2D1(pts[dir1].size());
665 Array<OneD, NekDouble> Dx2D2(pts[dir1].size());
666
667 // Ideally this will be stored in m_geomfactors
668 bndXmap->PhysDeriv(pts[dir1], Dx1D1, Dx1D2);
669 bndXmap->PhysDeriv(pts[dir2], Dx2D1, Dx2D2);
670
671 // Initial the locCoord, in [-1,1]
672 locCoord[0] = 0.0;
673 locCoord[1] = 0.0;
674
675 NekDouble x1map, x2map, F1, F2;
676 NekDouble derx1_1, derx1_2, derx2_1, derx2_2, jac;
677 NekDouble resid;
678 int cnt = 0;
679 bool isConverge = false;
680
681 F1 = F2 = 2000; // Starting value of Function
682 resid = sqrt(F1 * F1 + F2 * F2);
683 while (cnt++ < iterMax)
684 {
685 x1map = bndXmap->PhysEvaluate(locCoord, pts[dir1]);
686 x2map = bndXmap->PhysEvaluate(locCoord, pts[dir2]);
687
688 F1 = gloCoord[dir1] - x1map;
689 F2 = gloCoord[dir2] - x2map;
690
691 if (F1 * F1 + F2 * F2 < scaledTol)
692 {
693 resid = sqrt(F1 * F1 + F2 * F2);
694 isConverge = true;
695 break;
696 }
697
698 // Interpolate derivative metric at locCoord
699 derx1_1 = bndXmap->PhysEvaluate(locCoord, Dx1D1);
700 derx1_2 = bndXmap->PhysEvaluate(locCoord, Dx1D2);
701 derx2_1 = bndXmap->PhysEvaluate(locCoord, Dx2D1);
702 derx2_2 = bndXmap->PhysEvaluate(locCoord, Dx2D2);
703
704 jac = derx2_2 * derx1_1 - derx2_1 * derx1_2;
705
706 // Use analytical inverse of derivitives which are
707 // also similar to those of metric factors.
708 locCoord[0] = locCoord[0] + (derx2_2 * (gloCoord[dir1] - x1map) -
709 derx1_2 * (gloCoord[dir2] - x2map)) /
710 jac;
711
712 locCoord[1] = locCoord[1] + (-derx2_1 * (gloCoord[dir1] - x1map) +
713 derx1_1 * (gloCoord[dir2] - x2map)) /
714 jac;
715
716 // locCoord have diverged so stop iteration
717 if (!(std::isfinite(locCoord[0]) && std::isfinite(locCoord[1])))
718 {
719 dist = 1e16;
720 std::ostringstream ss;
721 ss << "nan or inf found in NewtonIterForLocCoordOnProjBndElmt in "
722 "element "
723 << bndGeom->GetGlobalID();
724 WARNINGL1(false, ss.str());
725 return false;
726 }
727 if (fabs(locCoord[0]) > LcoordDiv || fabs(locCoord[1]) > LcoordDiv)
728 {
729 break;
730 }
731 }
732
733 // Check distance for collapsed coordinate
735 bndXmap->LocCoordToLocCollapsed(locCoord, eta);
736
737 if (bndGeom->ClampLocCoords(eta, 0.0))
738 {
739 // calculate the global point corresponding to locCoord
740 x1map = bndXmap->PhysEvaluate(eta, pts[dir1]);
741 x2map = bndXmap->PhysEvaluate(eta, pts[dir2]);
742
743 F1 = gloCoord[dir1] - x1map;
744 F2 = gloCoord[dir2] - x2map;
745
746 dist = sqrt(F1 * F1 + F2 * F2);
747 }
748 else
749 {
750 dist = 0.0;
751 }
752
753 // Warning if failed
754 if (cnt >= iterMax)
755 {
756 Array<OneD, NekDouble> collCoords(2);
757 bndXmap->LocCoordToLocCollapsed(locCoord, collCoords);
758
759 // if coordinate is inside element dump error!
760 if ((collCoords[0] >= -1.0 && collCoords[0] <= 1.0) &&
761 (collCoords[1] >= -1.0 && collCoords[1] <= 1.0))
762 {
763 std::ostringstream ss;
764
765 ss << "Reached MaxIterations (" << iterMax
766 << ") in Newton iteration ";
767 ss << "Init value (" << setprecision(4) << 0 << "," << 0 << ","
768 << ") ";
769 ss << "Fin value (" << locCoord[0] << "," << locCoord[1] << ","
770 << ") ";
771 ss << "Resid = " << resid << " Tolerance = " << sqrt(scaledTol);
772
773 WARNINGL1(cnt < iterMax, ss.str());
774 }
775 }
776
777 return isConverge;
778}
779
780/**
781 * @brief Check if a point can be projected onto an oundary element in a given
782 * direction. If yes, give the local coordinates of the projected point.
783 * we have checked the projected point is inside the projected element.
784 * @param bndGeom Pointer to the geometry of the boundary element.
785 * @param gloCoord Global coordinate of the point. size=3.
786 * @param projDir Projection direction, which is used as the reference
787 * direction in the 3D routine. size=3, norm=1.
788 * @param locCoord Iteration results for local coordinates (if inside).
789 * @param projDist Projection distance betweem the point to the wall point.
790 * @param maxDist Disntance to check if the wall point is desired.
791 * @param iterTol Tolerence for iteration.
792 * @return Inside (true) or not (false)
793 */
796 const Array<OneD, const NekDouble> &gloCoord,
797 const Array<OneD, const NekDouble> &projDir,
798 Array<OneD, NekDouble> &locCoord, NekDouble &projDist,
799 const NekDouble maxDist, const NekDouble iterTol)
800{
801 // Get variables
802 StdRegions::StdExpansionSharedPtr bndXmap = bndGeom->GetXmap();
803 const int npts = bndXmap->GetTotPoints();
804 const int nCoordDim = m_f->m_exp[0]->GetCoordim(0); // =2 for 2.5D cases
805
806 Array<OneD, Array<OneD, const NekDouble>> bndCoeffs(nCoordDim);
808 Array<OneD, Array<OneD, NekDouble>> projPts(nCoordDim);
809 Array<OneD, NekDouble> projGloCoord(3, 0.0);
810
811 for (int i = 0; i < nCoordDim; ++i)
812 {
813 pts[i] = Array<OneD, NekDouble>(npts);
814 projPts[i] = Array<OneD, NekDouble>(npts);
815 bndCoeffs[i] = bndGeom->GetCoeffs(i); // 0/1/2 for x/y/z
816 bndXmap->BwdTrans(bndCoeffs[i], pts[i]);
817 }
818
819 // Project the point and vertices of the element in the input direction
820 ProjectPoint(gloCoord, projDir, 0.0, projGloCoord);
821 ProjectVertices(pts, projDir, 0.0, projPts);
822
823 // Set the main direction(s) and the minor direction
824 // The gloCoord for minor direction will not be used for locCoord iteration
825 // dirUse[0] is the main dir for 2D/2.5D cases, dirUse[1]/[2] is the minor
826 // one dirUse[0] and dirUse[1] are the main dir for 3D cases, dirUse[2] is
827 // hte minor one
828 int dirMaxId = 0; // id to get the dir with largest projDir component
829 for (int i = 1; i < nCoordDim; ++i)
830 {
831 if (fabs(projDir[i]) > fabs(projDir[dirMaxId]))
832 {
833 dirMaxId = i;
834 }
835 }
836
837 Array<OneD, int> dirUse(3, 0);
838 if (nCoordDim == 2)
839 {
840 // 2D or 2.5D cases
841 if (dirMaxId == 0)
842 {
843 dirUse[0] = 1;
844 dirUse[1] = 0;
845 dirUse[2] = 2;
846 }
847 else
848 {
849 dirUse[0] = 0;
850 dirUse[1] = 1;
851 dirUse[2] = 2;
852 }
853 }
854 else
855 {
856 // 3D cases
857 if (dirMaxId == 0)
858 {
859 dirUse[0] = 1;
860 dirUse[1] = 2;
861 dirUse[2] = 0;
862 }
863 else if (dirMaxId == 1)
864 {
865 dirUse[0] = 2;
866 dirUse[1] = 0;
867 dirUse[2] = 1;
868 }
869 else
870 {
871 dirUse[0] = 0;
872 dirUse[1] = 1;
873 dirUse[2] = 2;
874 }
875 }
876
877 // Check if the projected point is in the projected elmt
878 // If yes, then compute the locCoord and check if desired point is found
879 if (nCoordDim == 2)
880 {
881 if (isInProjectedArea2D(projGloCoord, projPts, 1.0e-12))
882 {
883 bool isConverge, isDesired;
884
886 bndGeom, projGloCoord, projPts, dirUse, locCoord, iterTol);
887
888 Array<OneD, NekDouble> tmp(2, 0.0);
889 tmp[0] = bndXmap->PhysEvaluate(locCoord, pts[0]) - gloCoord[0];
890 tmp[1] = bndXmap->PhysEvaluate(locCoord, pts[1]) - gloCoord[1];
891 projDist = Vmath::Dot(2, tmp, 1, projDir, 1); // can be negative
892
893 isDesired = (projDist > 0.0) && (projDist < maxDist);
894
895 return isConverge && isDesired;
896 }
897 else
898 {
899 return false;
900 }
901 }
902 else
903 {
904 if (isInProjectedArea3D(projGloCoord, projPts, projDir, 1.0e-12,
905 1.0e-6))
906 {
907 NekDouble dist;
908 bool isConverge, isDesired;
909
910 isConverge =
911 NewtonIterForLocCoordOnBndElmt(bndGeom, projGloCoord, projPts,
912 dirUse, locCoord, dist, iterTol);
913
914 if (dist > iterTol)
915 {
916 std::ostringstream ss;
917 ss << "Collapsed locCoord out of range.\n"
918 << "Newton iteration gives the distance: " << dist;
919 WARNINGL1(false, ss.str());
920 }
921
922 Array<OneD, NekDouble> tmp(3, 0.0);
923 tmp[0] = bndXmap->PhysEvaluate(locCoord, pts[0]) - gloCoord[0];
924 tmp[1] = bndXmap->PhysEvaluate(locCoord, pts[1]) - gloCoord[1];
925 tmp[2] = bndXmap->PhysEvaluate(locCoord, pts[2]) - gloCoord[2];
926 projDist = Vmath::Dot(3, tmp, 1, projDir, 1); // can be negative
927
928 isDesired = (projDist > 0.0) && (projDist < maxDist);
929
930 return isConverge && isDesired;
931 }
932 else
933 {
934 return false;
935 }
936 }
937}
938
939/**
940 * @brief Get the normals for a given locCoord
941 * @param bndGeom Pointer to the geometry of the boundary element.
942 * @param locCoord Iteration results for local coordinates (if inside).
943 * @param normals Wall normal as the result
944 */
947 const Array<OneD, const NekDouble> &locCoord,
948 Array<OneD, NekDouble> &normals)
949{
950 const int nCoordDim = m_f->m_exp[0]->GetCoordim(0); // =2 for 2.5D cases
951 StdRegions::StdExpansionSharedPtr bndXmap = bndGeom->GetXmap();
954 int npts = bndXmap->GetTotPoints();
955
956 // Get pts in the element
957 for (int i = 0; i < nCoordDim; ++i)
958 {
959 pts[i] = Array<OneD, NekDouble>(npts);
960 bndCoeffs[i] = bndGeom->GetCoeffs(i); // 0/1/2 for x/y/z
961 bndXmap->BwdTrans(bndCoeffs[i], pts[i]);
962 }
963
964 // Get the outward-pointing normals according to the given locCoord
965 if (nCoordDim == 2)
966 {
967 Array<OneD, NekDouble> DxD1(pts[0].size());
968 Array<OneD, NekDouble> DyD1(pts[0].size());
969
970 bndXmap->PhysDeriv(pts[0], DxD1);
971 bndXmap->PhysDeriv(pts[1], DyD1);
972
973 NekDouble dxd1, dyd1, norm;
974 dxd1 = bndXmap->PhysEvaluate(locCoord, DxD1);
975 dyd1 = bndXmap->PhysEvaluate(locCoord, DyD1);
976 norm = sqrt(dxd1 * dxd1 + dyd1 * dyd1);
977
978 normals[0] = dyd1 / norm;
979 normals[1] = -dxd1 / norm;
980 normals[2] = 0.0;
981 }
982 else
983 {
984 Array<OneD, NekDouble> DxD1(pts[0].size());
985 Array<OneD, NekDouble> DxD2(pts[0].size());
986 Array<OneD, NekDouble> DyD1(pts[0].size());
987 Array<OneD, NekDouble> DyD2(pts[0].size());
988 Array<OneD, NekDouble> DzD1(pts[0].size());
989 Array<OneD, NekDouble> DzD2(pts[0].size());
990
991 bndXmap->PhysDeriv(pts[0], DxD1, DxD2);
992 bndXmap->PhysDeriv(pts[1], DyD1, DyD2);
993 bndXmap->PhysDeriv(pts[2], DzD1, DzD2);
994
995 NekDouble dxd1, dyd1, dzd1, dxd2, dyd2, dzd2;
996 dxd1 = bndXmap->PhysEvaluate(locCoord, DxD1);
997 dxd2 = bndXmap->PhysEvaluate(locCoord, DxD2);
998 dyd1 = bndXmap->PhysEvaluate(locCoord, DyD1);
999 dyd2 = bndXmap->PhysEvaluate(locCoord, DyD2);
1000 dzd1 = bndXmap->PhysEvaluate(locCoord, DzD1);
1001 dzd2 = bndXmap->PhysEvaluate(locCoord, DzD2);
1002
1003 NekDouble n1, n2, n3, norm;
1004 n1 = dyd1 * dzd2 - dyd2 * dzd1;
1005 n2 = dzd1 * dxd2 - dzd2 * dxd1;
1006 n3 = dxd1 * dyd2 - dxd2 * dyd1;
1007 norm = sqrt(n1 * n1 + n2 * n2 + n3 * n3);
1008
1009 normals[0] = n1 / norm;
1010 normals[1] = n2 / norm;
1011 normals[2] = n3 / norm;
1012 }
1013}
1014
1015} // namespace FieldUtils
1016} // namespace Nektar
#define WARNINGL1(condition, msg)
Definition: ErrorUtil.hpp:250
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
A class that contains algorithms for interpolation between pts fields, expansions and different meshe...
FIELD_UTILS_EXPORT void Interpolate(const T expInField, T &expOutField, NekDouble def_value=0.0)
Interpolate from an expansion to an expansion.
FieldSharedPtr m_f
Field object.
Definition: Module.h:234
std::map< std::string, ConfigOption > m_config
List of configuration values.
Definition: Module.h:263
This processing module sets up for the boundary field to be extracted.
virtual void v_Process(po::variables_map &vm) override
void ProjectPoint(const Array< OneD, const NekDouble > &gloCoord, const Array< OneD, const NekDouble > &projDir, const NekDouble distToOrig, Array< OneD, NekDouble > &projGloCoord)
Project a single point along the given direction to a plane.
bool BndElmtContainsPoint(SpatialDomains::GeometrySharedPtr bndGeom, const Array< OneD, const NekDouble > &gloCoord, const Array< OneD, const NekDouble > &projDir, Array< OneD, NekDouble > &locCoord, NekDouble &projDist, const NekDouble maxDist=1.0, const NekDouble iterTol=1.0e-8)
Check if a point can be projected onto an oundary element in a given direction. If yes,...
bool NewtonIterForLocCoordOnBndElmt(SpatialDomains::GeometrySharedPtr bndGeom, const Array< OneD, const NekDouble > &gloCoord, const Array< OneD, const Array< OneD, NekDouble > > &pts, const Array< OneD, const int > &dirUse, Array< OneD, NekDouble > &locCoord, NekDouble &dist, const NekDouble iterTol=1.0e-8, const int iterMax=51)
bool isInProjectedArea2D(const Array< OneD, const NekDouble > &projGloCoord, const Array< OneD, const Array< OneD, NekDouble > > &projPts, const NekDouble paralTol=1.0e-12)
Determine if the projected point is inside the projected element.
bool isInProjectedArea3D(const Array< OneD, const NekDouble > &projGloCoord, const Array< OneD, const Array< OneD, NekDouble > > &projPts, const Array< OneD, const NekDouble > &projDir, const NekDouble paralTol=1.0e-12, const NekDouble angleTol=1.0e-6)
bool BisectionForLocCoordOnBndElmt(SpatialDomains::GeometrySharedPtr bndGeom, const Array< OneD, const NekDouble > &gloCoord, const Array< OneD, const Array< OneD, NekDouble > > &pts, const Array< OneD, const int > &dirUse, Array< OneD, NekDouble > &locCoord, const NekDouble iterTol=1.0e-8, const int iterMax=51)
Use iteration to get the locCoord. This routine should be used after we have checked the projected po...
static std::shared_ptr< Module > create(FieldSharedPtr f)
Creates an instance of this class.
void GetNormals(SpatialDomains::GeometrySharedPtr bndGeom, const Array< OneD, const NekDouble > &locCoord, Array< OneD, NekDouble > &normals)
Get the normals for a given locCoord.
virtual void v_Process(po::variables_map &vm) override
Write mesh to output file.
void ProjectVertices(const Array< OneD, const Array< OneD, NekDouble > > &pts, const Array< OneD, const NekDouble > &projDir, const NekDouble distToOrig, Array< OneD, Array< OneD, NekDouble > > &projPts)
Project a single point along the given direction to a plane.
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:198
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
static bool GenerateVector(const std::string &str, std::vector< T > &out)
Takes a comma-separated string and converts it to entries in a vector.
Definition: ParseUtils.cpp:131
const BoundaryRegionCollection & GetBoundaryRegions(void) const
Definition: Conditions.h:234
std::shared_ptr< Field > FieldSharedPtr
Definition: Field.hpp:991
std::pair< ModuleType, std::string > ModuleKey
Definition: Module.h:317
ModuleFactory & GetModuleFactory()
Definition: Module.cpp:49
static std::map< PtsInfo, int > NullPtsInfoMap
Definition: PtsField.h:70
static const NekDouble kNekUnsetDouble
std::map< int, BoundaryRegionShPtr > BoundaryRegionCollection
Definition: Conditions.h:210
std::shared_ptr< Geometry > GeometrySharedPtr
Definition: Geometry.h:54
std::shared_ptr< StdExpansion > StdExpansionSharedPtr
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
double NekDouble
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:513
T Vsum(int n, const T *x, const int incx)
Subtract return sum(x)
Definition: Vmath.cpp:890
T Dot(int n, const T *w, const T *x)
dot (vector times vector): z = w*x
Definition: Vmath.cpp:1095
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.cpp:354
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition: Vmath.cpp:245
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294
Represents a command-line configuration option.
Definition: Module.h:131