Loading...
Searching...
No Matches
NURBS3DVolume.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2007-2023 PCOpt/NTUA
9 Copyright (C) 2013-2023 FOSS GP
10 Copyright (C) 2019-2021 OpenCFD Ltd.
11-------------------------------------------------------------------------------
12License
13 This file is part of OpenFOAM.
14
15 OpenFOAM is free software: you can redistribute it and/or modify it
16 under the terms of the GNU General Public License as published by
17 the Free Software Foundation, either version 3 of the License, or
18 (at your option) any later version.
19
20 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
21 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
22 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
23 for more details.
24
25 You should have received a copy of the GNU General Public License
26 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
27
28\*---------------------------------------------------------------------------*/
29
31#include "NURBS3DVolume.H"
32
33#include "OFstream.H"
34#include "Time.H"
35#include "deltaBoundary.H"
38
39// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40
41namespace Foam
42{
45}
46
47
48// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
49
51{
52 // It is considered an error to recompute points in the control boxes
54 {
56 << "Attempting to recompute points residing within control boxes"
57 << exit(FatalError);
58 }
59
60 mapPtr_.reset(new labelList(meshPoints.size(), -1));
61 reverseMapPtr_.reset(new labelList(meshPoints.size(), -1));
62 labelList& map = mapPtr_();
63 labelList& reverseMap = reverseMapPtr_();
64
65 // Identify points inside morphing boxes
66 scalar lowerX = min(cps_.component(0));
67 scalar upperX = max(cps_.component(0));
68 scalar lowerY = min(cps_.component(1));
69 scalar upperY = max(cps_.component(1));
70 scalar lowerZ = min(cps_.component(2));
71 scalar upperZ = max(cps_.component(2));
72
73 Info<< "Control Points bounds \n"
74 << "\tX1 : (" << lowerX << " " << upperX << ")\n"
75 << "\tX2 : (" << lowerY << " " << upperY << ")\n"
76 << "\tX3 : (" << lowerZ << " " << upperZ << ")\n" << endl;
77
78 label count(0);
79 forAll(meshPoints, pI)
80 {
81 const vector& pointI = meshPoints[pI];
82 if
83 (
84 pointI.x() >= lowerX && pointI.x() <= upperX
85 && pointI.y() >= lowerY && pointI.y() <= upperY
86 && pointI.z() >= lowerZ && pointI.z() <= upperZ
87 )
88 {
89 map[count] = pI;
90 reverseMap[pI] = count;
91 ++count;
92 }
93 }
94
95 // Resize lists
96 map.setSize(count);
98 reduce(count, sumOp<label>());
99 Info<< "Initially found " << count << " points inside control boxes"
100 << endl;
101}
102
103
105(
106 const vectorField& points
107)
108{
109 scalar timeBef = mesh_.time().elapsedCpuTime();
110
111 if (parametricCoordinatesPtr_)
112 {
114 << "Attempting to recompute parametric coordinates"
115 << exit(FatalError);
116 }
117
118 labelList& map = mapPtr_();
119 labelList& reverseMap = reverseMapPtr_();
120
121 parametricCoordinatesPtr_.reset
122 (
124 (
125 IOobject
126 (
127 "parametricCoordinates" + name_,
128 mesh_.time().timeName(),
129 mesh_,
132 ),
133 pointMesh::New(mesh_),
135 )
136 );
137 vectorField& paramCoors = parametricCoordinatesPtr_().primitiveFieldRef();
138
139 // If already present, read from file
140 if
141 (
142 parametricCoordinatesPtr_().typeHeaderOk<pointVectorField>(true)
143 && readStoredData_
144 )
145 {
146 // These are the parametric coordinates that have been computed
147 // correctly (hopefully!) in a previous run.
148 // findPointsInBox has possibly found more points and lists need
149 // resizing
150 Info<< "Reading parametric coordinates from file" << endl;
151 IOobject& header = parametricCoordinatesPtr_().ref();
152 parametricCoordinatesPtr_() =
154 (
155 header,
156 pointMesh::New(mesh_)
157 );
158
159 // Initialize intermediate fields with sizes from findPointInBox
160 labelList actualMap(map.size());
161
162 // Read and store
163 label curIndex(0);
164 forAll(map, pI)
165 {
166 const label globalPointIndex = map[pI];
167 if (paramCoors[globalPointIndex] != vector::zero)
168 {
169 actualMap[curIndex] = map[pI];
170 reverseMap[globalPointIndex] = curIndex;
171 ++curIndex;
172 }
173 else
174 {
175 reverseMap[globalPointIndex] = -1;
176 }
177 }
178
179 // Resize intermediate
180 actualMap.setSize(curIndex);
181
182 reduce(curIndex, sumOp<label>());
183 Info<< "Read non-zero parametric coordinates for " << curIndex
184 << " points" << endl;
185
186 // Update lists with the appropriate entries
187 map = actualMap;
188 }
189 // Else, compute parametric coordinates iteratively
190 else
191 {
192 // Initialize parametric coordinates based on min/max of control points
193 scalar minX1 = min(cps_.component(0));
194 scalar maxX1 = max(cps_.component(0));
195 scalar minX2 = min(cps_.component(1));
196 scalar maxX2 = max(cps_.component(1));
197 scalar minX3 = min(cps_.component(2));
198 scalar maxX3 = max(cps_.component(2));
199
200 scalar oneOverDenomX(1./(maxX1 - minX1));
201 scalar oneOverDenomY(1./(maxX2 - minX2));
202 scalar oneOverDenomZ(1./(maxX3 - minX3));
203
204 forAll(map, pI)
205 {
206 const label globalPI = map[pI];
207 paramCoors[globalPI].x() = (points[pI].x() - minX1)*oneOverDenomX;
208 paramCoors[globalPI].y() = (points[pI].y() - minX2)*oneOverDenomY;
209 paramCoors[globalPI].z() = (points[pI].z() - minX3)*oneOverDenomZ;
210 }
211
212 // Indices of points that failed to converge
213 // (i.e. are bounded for nMaxBound iters)
214 boolList dropOffPoints(map.size(), false);
215 label nDropedPoints(0);
216
217 // Initial cartesian coordinates
218 tmp<vectorField> tsplinesBasedCoors(coordinates(paramCoors));
219 vectorField& splinesBasedCoors = tsplinesBasedCoors.ref();
220
221 // Newton-Raphson loop to compute parametric coordinates
222 // based on cartesian coordinates and the known control points
223 Info<< "Mapping of mesh points to parametric space for box " << name_
224 << " ..." << endl;
225 // Do loop on a point-basis and check residual of each point equation
226 label maxIterNeeded(0);
227 forAll(points, pI)
228 {
229 label iter(0);
230 label nBoundIters(0);
231 vector res(GREAT, GREAT, GREAT);
232 do
233 {
234 const label globalPI = map[pI];
235 vector& uVec = paramCoors[globalPI];
236 vector& coorPointI = splinesBasedCoors[pI];
237 uVec += ((inv(JacobianUVW(uVec))) & (points[pI] - coorPointI));
238 // Bounding might be needed for the first iterations
239 // If multiple bounds happen, point is outside of the control
240 // boxes and should be discarded
241 if (bound(uVec))
242 {
243 ++nBoundIters;
244 }
245 if (nBoundIters > nMaxBound_)
246 {
247 dropOffPoints[pI] = true;
248 ++nDropedPoints;
249 break;
250 }
251 // Update current cartesian coordinates based on parametric ones
252 coorPointI = coordinates(uVec);
253 // Compute residual
254 res = cmptMag(points[pI] - coorPointI);
255 }
256 while
257 (
258 (iter++ < maxIter_)
259 && (
260 res.component(0) > tolerance_
261 || res.component(1) > tolerance_
262 || res.component(2) > tolerance_
263 )
264 );
265 if (iter > maxIter_)
266 {
268 << "Mapping to parametric space for point " << pI
269 << " failed." << endl
270 << "Residual after " << maxIter_ + 1 << " iterations : "
271 << res << endl
272 << "parametric coordinates " << paramCoors[map[pI]]
273 << endl
274 << "Local system coordinates " << points[pI] << endl
275 << "Threshold residual per direction : " << tolerance_
276 << endl;
277 }
278 maxIterNeeded = max(maxIterNeeded, iter);
279 }
280 reduce(maxIterNeeded, maxOp<label>());
281
282 label nParameterizedPoints = map.size() - nDropedPoints;
283
284 // Resize mapping lists and parametric coordinates after dropping some
285 // points
286 labelList mapOld(map);
287
288 map.setSize(nParameterizedPoints);
289
290 label curIndex(0);
291 forAll(dropOffPoints, pI)
292 {
293 if (!dropOffPoints[pI])
294 {
295 map[curIndex] = mapOld[pI];
296 reverseMap[mapOld[pI]] = curIndex;
297 ++curIndex;
298 }
299 else
300 {
301 paramCoors[mapOld[pI]] = vector::zero;
302 reverseMap[mapOld[pI]] = -1;
303 }
304 }
305
306 reduce(nDropedPoints, sumOp<label>());
307 reduce(nParameterizedPoints, sumOp<label>());
308 Info<< "Found " << nDropedPoints
309 << " to discard from morphing boxes" << endl;
310 Info<< "Keeping " << nParameterizedPoints
311 << " parameterized points in boxes" << endl;
312
313 splinesBasedCoors = coordinates(paramCoors)();
314 scalar maxDiff(-GREAT);
315 forAll(splinesBasedCoors, pI)
316 {
317 scalar diff =
318 mag(splinesBasedCoors[pI] - localSystemCoordinates_[map[pI]]);
319 if (diff > maxDiff)
320 {
321 maxDiff = diff;
322 }
323 }
324 reduce(maxDiff, maxOp<scalar>());
325 scalar timeAft = mesh_.time().elapsedCpuTime();
326 Info<< "\tMapping completed in " << timeAft - timeBef << " seconds"
327 << endl;
328 Info<< "\tMax iterations per point needed to compute parametric "
329 << "coordinates : "
330 << maxIterNeeded << endl;
331 Info<< "\tMax difference between original mesh points and "
332 << "parameterized ones "
333 << maxDiff << endl;
334 }
335}
336
337
339(
340 tmp<vectorField> tPoints
342{
343 const vectorField& points = tPoints();
345}
346
347
349(
350 vector& vec,
351 scalar minValue,
352 scalar maxValue
353) const
354{
355 bool boundPoint(false);
356 // Lower value bounding
357 if (vec.x() < scalar(0))
358 {
359 vec.x() = minValue;
360 boundPoint = true;
361 }
362 if (vec.y() < scalar(0))
363 {
364 vec.y() = minValue;
365 boundPoint = true;
366 }
367 if (vec.z() < scalar(0))
368 {
369 vec.z() = minValue;
370 boundPoint = true;
371 }
372 // Upper value bounding
373 if (vec.x() > 1)
374 {
375 vec.x() = maxValue;
376 boundPoint = true;
377 }
378 if (vec.y() > 1)
379 {
380 vec.y() = maxValue;
381 boundPoint = true;
382 }
383 if (vec.z() > 1)
384 {
385 vec.z() = maxValue;
386 boundPoint = true;
387 }
388 return boundPoint;
389}
390
391
393{
395 {
396 mkDir(mesh_.time().globalPath()/"optimisation"/cpsFolder_);
397 }
398}
399
400
402{
403 label nCPs = cps_.size();
404 activeControlPoints_ = boolList(nCPs, true);
405 activeDesignVariables_ = boolList(3*nCPs, true);
406
407 // Check whether all boundary control points should be confined
408 confineBoundaryControlPoints();
409
410 // Apply confinement to maintain continuity
411 continuityRealatedConfinement();
412
413 // Confine user-specified directions
414 confineControlPointsDirections();
415
416 // Determine active control points. A control point is considered active
417 // if at least one of its components is free to move
418 forAll(activeControlPoints_, cpI)
419 {
420 if
421 (
422 !activeDesignVariables_[3*cpI]
423 && !activeDesignVariables_[3*cpI + 1]
424 && !activeDesignVariables_[3*cpI + 2]
425 )
426 {
427 activeControlPoints_[cpI] = false;
428 }
430
431 // Deactivate design varibles if they affect no mesh point
432 confineInertControlPoints();
433}
434
435
437{
438 // Loop over all active control points and check whether they parameterize
439 // at least one mesh point of the ones laying within the morphing box
440 const labelList& map = getMap();
441 const pointVectorField& parametricCoors = getParametricCoordinates();
442
443 const scalarField& knotsU = basisU_.knots();
444 const scalarField& knotsV = basisV_.knots();
445 const scalarField& knotsW = basisW_.knots();
446
447 const label pU = basisU_.degree();
448 const label pV = basisV_.degree();
449 const label pW = basisW_.degree();
450 forAll(activeControlPoints_, cpI)
451 {
452 if (activeControlPoints_[cpI])
453 {
454 // Get i, j, k corresponding to this control point
455 label i(-1), j(-1), k(-1);
456 getIJK(i, j, k, cpI);
457
458 const scalar uMin(knotsU[i]);
459 const scalar uMax(knotsU[i + pU + 1]);
460 const scalar vMin(knotsV[j]);
461 const scalar vMax(knotsV[j + pV + 1]);
462 const scalar wMin(knotsW[k]);
463 const scalar wMax(knotsW[k + pW + 1]);
464 bool foundParamPt(false);
465 for (const label paramI : map)
466 {
467 const vector& paramCoors = parametricCoors()[paramI];
468 if
469 (
470 paramCoors.x() >= uMin && paramCoors.x() < uMax
471 && paramCoors.y() >= vMin && paramCoors.y() < vMax
472 && paramCoors.z() >= wMin && paramCoors.z() < wMax
473 )
474 {
475 foundParamPt = true;
476 break;
477 }
478 }
479 UPstream::reduceOr(foundParamPt);
480 if (!foundParamPt)
481 {
482 activeControlPoints_[cpI] = false;
483 activeDesignVariables_[3*cpI] = false;
484 activeDesignVariables_[3*cpI + 1] = false;
485 activeDesignVariables_[3*cpI + 2] = false;
486 Info<< "Disabling control " << cpI
487 << " and variables "
488 << "(" << 3*cpI << ", " << 3*cpI+1 << ", " << 3*cpI+2 << ")"
489 << " since they does not parameterize any mesh point"
491 }
492 }
493 }
494}
495
496
498{
499 const label nCPsU = basisU_.nCPs();
500 const label nCPsV = basisV_.nCPs();
501 const label nCPsW = basisW_.nCPs();
502
503 // Zero movement of the boundary control points. Active by default
504 if (confineBoundaryControlPoints_)
505 {
506 // Side patches
507 for (label iCPw = 0; iCPw < nCPsW; iCPw += nCPsW - 1)
508 {
509 for (label iCPv = 0; iCPv < nCPsV; iCPv++)
510 {
511 for (label iCPu = 0; iCPu < nCPsU; iCPu++)
512 {
513 confineControlPoint(getCPID(iCPu, iCPv, iCPw));
514 }
515 }
516 }
517 // Front-back patches
518 for (label iCPw = 0; iCPw < nCPsW; iCPw++)
519 {
520 for (label iCPv = 0; iCPv < nCPsV; iCPv++)
521 {
522 for (label iCPu = 0; iCPu < nCPsU; iCPu += nCPsU - 1)
523 {
524 confineControlPoint(getCPID(iCPu, iCPv, iCPw));
525 }
526 }
527 }
528 // Top-bottom patches
529 for (label iCPw = 0; iCPw < nCPsW; iCPw++)
530 {
531 for (label iCPv = 0; iCPv < nCPsV; iCPv += nCPsV - 1)
532 {
533 for (label iCPu = 0; iCPu < nCPsU; iCPu++)
534 {
535 confineControlPoint(getCPID(iCPu, iCPv, iCPw));
537 }
538 }
539 }
540}
541
542
544{
545 const label nCPsU = basisU_.nCPs();
546 const label nCPsV = basisV_.nCPs();
547 const label nCPsW = basisW_.nCPs();
548
549 // Zero movement to a number of x-constant slices of cps in order to
550 // preserve continuity at the boundary of the parameterized space
551 forAll(confineUMinCPs_, iCPu)
552 {
553 const boolVector& confineSlice = confineUMinCPs_[iCPu];
554 // Control points at the start of the parameterized space
555 for (label iCPw = 0; iCPw < nCPsW; iCPw++)
556 {
557 for (label iCPv = 0; iCPv < nCPsV; iCPv++)
558 {
559 confineControlPoint(getCPID(iCPu, iCPv, iCPw), confineSlice);
560 }
561 }
562 }
563
564 forAll(confineUMaxCPs_, sliceI)
565 {
566 const boolVector& confineSlice = confineUMaxCPs_[sliceI];
567 label iCPu = nCPsU - 1 - sliceI;
568 // Control points at the end of the parameterized space
569 for (label iCPw = 0; iCPw < nCPsW; iCPw++)
570 {
571 for (label iCPv = 0; iCPv < nCPsV; iCPv++)
572 {
573 confineControlPoint(getCPID(iCPu, iCPv, iCPw), confineSlice);
574 }
575 }
576 }
577
578 // Zero movement to a number of y-constant slices of cps in order to
579 // preserve continuity at the boundary of the parameterized space
580 forAll(confineVMinCPs_, iCPv)
581 {
582 const boolVector& confineSlice = confineVMinCPs_[iCPv];
583 // Control points at the start of the parameterized space
584 for (label iCPw = 0; iCPw < nCPsW; iCPw++)
585 {
586 for (label iCPu = 0; iCPu < nCPsU; iCPu++)
587 {
588 confineControlPoint(getCPID(iCPu, iCPv, iCPw), confineSlice);
589 }
590 }
591 }
592
593 forAll(confineVMaxCPs_, sliceI)
594 {
595 const boolVector& confineSlice = confineVMaxCPs_[sliceI];
596 label iCPv = nCPsV - 1 - sliceI;
597 // Control points at the end of the parameterized space
598 for (label iCPw = 0; iCPw < nCPsW; iCPw++)
599 {
600 for (label iCPu = 0; iCPu < nCPsU; iCPu++)
601 {
602 confineControlPoint(getCPID(iCPu, iCPv, iCPw), confineSlice);
603 }
604 }
605 }
606
607 // Zero movement to a number of w-constant slices of cps in order to
608 // preserve continuity at the boundary of the parameterized space
609 forAll(confineWMinCPs_, iCPw)
610 {
611 const boolVector& confineSlice = confineWMinCPs_[iCPw];
612 // Control points at the start of the parameterized space
613 for (label iCPv = 0; iCPv < nCPsV; iCPv++)
614 {
615 for (label iCPu = 0; iCPu < nCPsU; iCPu++)
616 {
617 confineControlPoint(getCPID(iCPu, iCPv, iCPw), confineSlice);
618 }
619 }
620 }
621
622 forAll(confineWMaxCPs_, sliceI)
623 {
624 const boolVector& confineSlice = confineWMaxCPs_[sliceI];
625 label iCPw = nCPsW - 1 - sliceI;
626 // Control points at the end of the parameterized space
627 for (label iCPv = 0; iCPv < nCPsV; iCPv++)
628 {
629 for (label iCPu = 0; iCPu < nCPsU; iCPu++)
630 {
631 confineControlPoint(getCPID(iCPu, iCPv, iCPw), confineSlice);
632 }
633 }
634 }
635}
636
637
639{
640 for (label cpI = 0; cpI < cps_.size(); ++cpI)
641 {
643 if (confineVMovement_) activeDesignVariables_[3*cpI + 1] = false;
644 if (confineWMovement_) activeDesignVariables_[3*cpI + 2] = false;
645 }
646}
647
648
650{
651 if (cpI < 0 || cpI > cps_.size() -1)
652 {
654 << "Attempted to confine control point movement for a control point "
655 << " ID which is out of bounds"
656 << exit(FatalError);
657 }
658 else
659 {
661 activeDesignVariables_[3*cpI + 1] = false;
662 activeDesignVariables_[3*cpI + 2] = false;
663 }
664}
665
666
668(
669 const label cpI,
670 const boolVector& confineDirections
671)
672{
673 if (cpI < 0 || cpI > cps_.size() -1)
674 {
676 << "Attempted to confine control point movement for a control point "
677 << " ID which is out of bounds"
678 << exit(FatalError);
679 }
680 else
681 {
682 if (confineDirections.x()) activeDesignVariables_[3*cpI] = false;
683 if (confineDirections.y()) activeDesignVariables_[3*cpI + 1] = false;
684 if (confineDirections.z()) activeDesignVariables_[3*cpI + 2] = false;
685 }
686}
687
688
689// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
690
692(
693 const dictionary& dict,
694 const fvMesh& mesh,
695 bool computeParamCoors
696)
697:
698 localIOdictionary
699 (
700 IOobject
701 (
702 dict.dictName() + "cpsBsplines",
703 mesh.time().timeName(),
704 fileName("uniform")/fileName("volumetricBSplines"),
705 mesh,
706 IOobject::READ_IF_PRESENT,
707 IOobject::AUTO_WRITE
708 ),
709 word::null
710 ),
711 mesh_(mesh),
712 dict_(dict),
713 name_(dict.dictName()),
714 basisU_(dict.get<label>("nCPsU"), dict.get<label>("degreeU")),
715 basisV_(dict.get<label>("nCPsV"), dict.get<label>("degreeV")),
716 basisW_(dict.get<label>("nCPsW"), dict.get<label>("degreeW")),
717 maxIter_(dict.getOrDefault<label>("maxIterations", 10)),
718 tolerance_(dict.getOrDefault<scalar>("tolerance", 1.e-10)),
719 nMaxBound_(dict.getOrDefault<scalar>("nMaxBoundIterations", 4)),
720 cps_(),
721 mapPtr_(nullptr),
722 reverseMapPtr_(nullptr),
723 parametricCoordinatesPtr_(nullptr),
724 localSystemCoordinates_(mesh_.nPoints(), Zero),
725 confineUMovement_
726 (
727 dict.getOrDefaultCompat<bool>
728 (
729 "confineUMovement", {{"confineX1movement", 1912}}, false
730 )
731 ),
732 confineVMovement_
733 (
735 (
736 "confineVMovement", {{"confineX2movement", 1912}}, false
737 )
738 ),
739 confineWMovement_
740 (
742 (
743 "confineWMovement", {{"confineX3movement", 1912}}, false
744 )
745 ),
746 confineBoundaryControlPoints_
747 (
748 dict.getOrDefault<bool>("confineBoundaryControlPoints", true)
749 ),
750 confineUMinCPs_
751 (
752 dict.getOrDefaultCompat<boolVectorList>
753 (
754 "confineUMinCPs", {{"boundUMinCPs", 1912}}, boolVectorList()
755 )
756 ),
757 confineUMaxCPs_
758 (
759 dict.getOrDefaultCompat<boolVectorList>
760 (
761 "confineUMaxCPs", {{"boundUMaxCPs", 1912}}, boolVectorList()
762 )
763 ),
764 confineVMinCPs_
765 (
766 dict.getOrDefaultCompat<boolVectorList>
767 (
768 "confineVMinCPs", {{"boundVMinCPs", 1912}}, boolVectorList()
769 )
770 ),
771 confineVMaxCPs_
772 (
773 dict.getOrDefaultCompat<boolVectorList>
774 (
775 "confineVMaxCPs", {{"boundVMaxCPs", 1912}}, boolVectorList()
776 )
777 ),
778 confineWMinCPs_
779 (
780 dict.getOrDefaultCompat<boolVectorList>
781 (
782 "confineWMinCPs", {{"boundWMinCPs", 1912}}, boolVectorList()
783 )
784 ),
785 confineWMaxCPs_
786 (
787 dict.getOrDefaultCompat<boolVectorList>
788 (
789 "confineWMaxCPs", {{"boundWMaxCPs", 1912}}, boolVectorList()
790 )
791 ),
792 activeControlPoints_(0), //zero here, execute sanity checks first
793 activeDesignVariables_(0), //zero here, execute sanity checks first
794 cpsFolder_("controlPoints"),
795 readStoredData_(dict.getOrDefault<bool>("readStoredData", true))
796{
797 // Create folders
798 makeFolders();
799
800 // Sanity checks
801 if
802 (
803 (confineUMinCPs_.size() + confineUMaxCPs_.size() >= basisU_.nCPs())
804 || (confineVMinCPs_.size() + confineVMaxCPs_.size() >= basisV_.nCPs())
805 || (confineWMinCPs_.size() + confineWMaxCPs_.size() >= basisW_.nCPs())
806 )
807 {
809 << "Number of control point slices to be kept frozen at "
810 << "the boundaries is invalid \n"
811 << "Number of control points in u " << basisU_.nCPs() << "\n"
812 << "Number of control points in v " << basisV_.nCPs() << "\n"
813 << "Number of control points in w " << basisW_.nCPs() << "\n"
814 << exit(FatalError);
815 }
816
817 // Construct control points, if not already read from file
818 if (found("controlPoints"))
819 {
820 cps_ =
822 (
823 "controlPoints",
824 *this,
825 basisU_.nCPs()*basisV_.nCPs()*basisW_.nCPs()
826 );
827 }
828 else
829 {
831 }
832}
833
834
835// * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
836
838(
839 const dictionary& dict,
840 const fvMesh& mesh,
841 bool computeParamCoors
842)
843{
844 const word modelType(dict.get<word>("type"));
845
846 Info<< "NURBS3DVolume type : " << modelType << endl;
847
848 auto* ctorPtr = dictionaryConstructorTable(modelType);
849
850 if (!ctorPtr)
851 {
853 (
854 dict,
855 "type",
856 modelType,
857 *dictionaryConstructorTablePtr_
858 ) << exit(FatalIOError);
859 }
860
861 return autoPtr<NURBS3DVolume>(ctorPtr(dict, mesh, computeParamCoors));
862}
863
864
865// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
866
867
869(
870 const vectorField& points
871) const
872{
874 vectorField& paramCoors = tparamCoors.ref();
875 // Initialize parametric coordinates based on min/max of control points
876 scalar minX1 = min(cps_.component(0));
877 scalar maxX1 = max(cps_.component(0));
878 scalar minX2 = min(cps_.component(1));
879 scalar maxX2 = max(cps_.component(1));
880 scalar minX3 = min(cps_.component(2));
881 scalar maxX3 = max(cps_.component(2));
882
883 scalar oneOverDenomX(1./(maxX1 - minX1));
884 scalar oneOverDenomY(1./(maxX2 - minX2));
885 scalar oneOverDenomZ(1./(maxX3 - minX3));
886
887 forAll(points, pI)
888 {
889 paramCoors[pI].x() = (points[pI].x() - minX1)*oneOverDenomX;
890 paramCoors[pI].y() = (points[pI].y() - minX2)*oneOverDenomY;
891 paramCoors[pI].z() = (points[pI].z() - minX3)*oneOverDenomZ;
892 }
893
894 // Indices of points that failed to converge
895 // (i.e. are bounded for nMaxBound iters)
896 boolList dropOffPoints(points.size(), false);
897 label nDropedPoints(0);
898
899 // Initial cartesian coordinates
900 tmp<vectorField> tsplinesBasedCoors(coordinates(paramCoors));
901 vectorField& splinesBasedCoors = tsplinesBasedCoors.ref();
902
903 // Newton-Raphson loop to compute parametric coordinates
904 // based on cartesian coordinates and the known control points
905 Info<< "Mapping of mesh points to parametric space for box " << name_
906 << " ..." << endl;
907 // Do loop on a point-basis and check residual of each point equation
908 label maxIterNeeded(0);
909 forAll(points, pI)
910 {
911 label iter(0);
912 label nBoundIters(0);
913 vector res(GREAT, GREAT, GREAT);
914 do
915 {
916 vector& uVec = paramCoors[pI];
917 vector& coorPointI = splinesBasedCoors[pI];
918 uVec += ((inv(JacobianUVW(uVec))) & (points[pI] - coorPointI));
919 // Bounding might be needed for the first iterations
920 // If multiple bounds happen, point is outside of the control
921 // boxes and should be discarded
922 if (bound(uVec))
923 {
924 ++nBoundIters;
925 }
926 if (nBoundIters > nMaxBound_)
927 {
928 dropOffPoints[pI] = true;
929 ++nDropedPoints;
930 break;
931 }
932 // Update current cartesian coordinates based on parametric ones
933 coorPointI = coordinates(uVec);
934 // Compute residual
935 res = cmptMag(points[pI] - coorPointI);
936 }
937 while
938 (
939 (iter++ < maxIter_)
940 && (
941 res.component(0) > tolerance_
942 || res.component(1) > tolerance_
943 || res.component(2) > tolerance_
944 )
945 );
946 if (iter > maxIter_)
947 {
949 << "Mapping to parametric space for point " << pI
950 << " failed." << endl
951 << "Residual after " << maxIter_ + 1 << " iterations : "
952 << res << endl
953 << "parametric coordinates " << paramCoors[pI]
954 << endl
955 << "Local system coordinates " << points[pI] << endl
956 << "Threshold residual per direction : " << tolerance_
957 << endl;
958 }
959 maxIterNeeded = max(maxIterNeeded, iter);
960 }
961 reduce(maxIterNeeded, maxOp<label>());
962
963 label nParameterizedPoints = points.size() - nDropedPoints;
964
965 reduce(nDropedPoints, sumOp<label>());
966 reduce(nParameterizedPoints, sumOp<label>());
967 Info<< "Found " << nDropedPoints
968 << " to discard from morphing boxes" << endl;
969 Info<< "Keeping " << nParameterizedPoints
970 << " parameterized points in boxes" << endl;
971 return tparamCoors;
972}
973
974
976(
977 const scalar u,
978 const scalar v,
979 const scalar w
980) const
981{
982 const label degreeU = basisU_.degree();
983 const label degreeV = basisV_.degree();
984 const label degreeW = basisW_.degree();
985
986 const label nCPsU = basisU_.nCPs();
987 const label nCPsV = basisV_.nCPs();
988 const label nCPsW = basisW_.nCPs();
989
990 vector derivative(Zero);
991
992 for (label iCPw = 0; iCPw < nCPsW; ++iCPw)
993 {
994 const scalar basisW(basisW_.basisValue(iCPw, degreeW, w));
995 for (label iCPv = 0; iCPv < nCPsV; ++iCPv)
996 {
997 const scalar basisVW = basisW*basisV_.basisValue(iCPv, degreeV, v);
998 for (label iCPu = 0; iCPu < nCPsU; ++iCPu)
999 {
1000 derivative +=
1001 cps_[getCPID(iCPu, iCPv, iCPw)]
1002 *basisU_.basisDerivativeU(iCPu, degreeU, u)
1003 *basisVW;
1004 }
1006 }
1007
1008 return derivative;
1009}
1010
1011
1013(
1014 const scalar u,
1015 const scalar v,
1016 const scalar w
1017) const
1018{
1019 const label degreeU = basisU_.degree();
1020 const label degreeV = basisV_.degree();
1021 const label degreeW = basisW_.degree();
1022
1023 const label nCPsU = basisU_.nCPs();
1024 const label nCPsV = basisV_.nCPs();
1025 const label nCPsW = basisW_.nCPs();
1026
1027 vector derivative(Zero);
1028
1029 for (label iCPw = 0; iCPw < nCPsW; ++iCPw)
1030 {
1031 const scalar basisW(basisW_.basisValue(iCPw, degreeW, w));
1032 for (label iCPv = 0; iCPv < nCPsV; ++iCPv)
1033 {
1034 const scalar basisWDeriV =
1035 basisW*basisV_.basisDerivativeU(iCPv, degreeV, v);
1036 for (label iCPu = 0; iCPu < nCPsU; ++iCPu)
1037 {
1038 derivative +=
1039 cps_[getCPID(iCPu, iCPv, iCPw)]
1040 *basisU_.basisValue(iCPu, degreeU, u)
1041 *basisWDeriV;
1042 }
1044 }
1045
1046 return derivative;
1047}
1048
1049
1051(
1052 const scalar u,
1053 const scalar v,
1054 const scalar w
1055) const
1056{
1057 const label degreeU = basisU_.degree();
1058 const label degreeV = basisV_.degree();
1059 const label degreeW = basisW_.degree();
1060
1061 const label nCPsU = basisU_.nCPs();
1062 const label nCPsV = basisV_.nCPs();
1063 const label nCPsW = basisW_.nCPs();
1064
1065 vector derivative(Zero);
1066
1067 for (label iCPw = 0; iCPw < nCPsW; iCPw++)
1068 {
1069 const scalar derivW(basisW_.basisDerivativeU(iCPw, degreeW, w));
1070 for (label iCPv = 0; iCPv < nCPsV; iCPv++)
1071 {
1072 const scalar derivWBasisV =
1073 derivW*basisV_.basisValue(iCPv, degreeV, v);
1074 for (label iCPu = 0; iCPu < nCPsU; iCPu++)
1075 {
1076 derivative +=
1077 cps_[getCPID(iCPu, iCPv, iCPw)]
1078 *basisU_.basisValue(iCPu, degreeU, u)
1079 *derivWBasisV;
1080 }
1082 }
1083
1084 return derivative;
1085}
1086
1087
1089(
1090 const vector& uVector
1091) const
1092{
1093 const scalar u = uVector.x();
1094 const scalar v = uVector.y();
1095 const scalar w = uVector.z();
1096
1097 vector uDeriv = volumeDerivativeU(u, v, w);
1098 vector vDeriv = volumeDerivativeV(u, v, w);
1099 vector wDeriv = volumeDerivativeW(u, v, w);
1101 tensor Jacobian(uDeriv, vDeriv, wDeriv, true);
1102
1103 return Jacobian;
1104}
1105
1106
1108(
1109 const vector& uVector,
1110 const label cpI
1111) const
1112{
1113 const scalar u = uVector.x();
1114 const scalar v = uVector.y();
1115 const scalar w = uVector.z();
1116
1117 const label nCPsU = basisU_.nCPs();
1118 const label nCPsV = basisV_.nCPs();
1119
1120 const label degreeU = basisU_.degree();
1121 const label degreeV = basisV_.degree();
1122 const label degreeW = basisW_.degree();
1123
1124 label iCPw = cpI/label(nCPsU*nCPsV);
1125 label iCPv = (cpI - iCPw*nCPsU*nCPsV)/nCPsU;
1126 label iCPu = (cpI - iCPw*nCPsU*nCPsV - iCPv*nCPsU);
1127
1128 // Normally, this should be a tensor, however the parameterization is
1129 // isotropic. Hence the tensor degenerates to a diagonal tensor with all
1130 // diagonal elements being equal. This returns the (unique) diag element
1131 scalar derivative =
1132 basisU_.basisValue(iCPu, degreeU, u)
1133 *basisV_.basisValue(iCPv, degreeV, v)
1134 *basisW_.basisValue(iCPw, degreeW, w);
1135
1136 return derivative;
1137}
1138
1139
1141(
1142 const pointVectorField& pointSens,
1143 const labelList& sensitivityPatchIDs
1144)
1145{
1146 vectorField controlPointDerivs(cps_.size(), Zero);
1147
1148 // Get parametric coordinates
1149 const vectorField& parametricCoordinates = getParametricCoordinates();
1150
1151 forAll(controlPointDerivs, cpI)
1152 {
1153 forAll(sensitivityPatchIDs, pI)
1154 {
1155 const label patchI = sensitivityPatchIDs[pI];
1156 const polyPatch& patch = mesh_.boundaryMesh()[patchI];
1157 const labelList& meshPoints = patch.meshPoints();
1158
1159 forAll(meshPoints, mpI)
1160 {
1161 const label globalIndex = meshPoints[mpI];
1162 const label whichPointInBox = reverseMapPtr_()[globalIndex];
1163
1164 // If point resides within control points box,
1165 // add contribution to cp derivative
1166 if (whichPointInBox != -1)
1167 {
1168 controlPointDerivs[cpI] +=
1169 (
1170 pointSens[globalIndex]
1171 & transformationTensorDxDb(globalIndex)
1172 )
1173 *volumeDerivativeCP
1174 (
1175 parametricCoordinates[globalIndex],
1176 cpI
1177 );
1178 }
1179 }
1180 }
1181 }
1182
1183 // Sum contributions from all processors
1184 Pstream::listReduce(controlPointDerivs, sumOp<vector>());
1185
1186 return controlPointDerivs;
1187}
1188
1189
1191(
1192 const volVectorField& faceSens,
1193 const labelList& sensitivityPatchIDs
1194)
1195{
1196 return
1197 computeControlPointSensitivities
1199 faceSens.boundaryField(),
1200 sensitivityPatchIDs
1201 );
1202}
1203
1204
1206(
1207 const boundaryVectorField& faceSens,
1208 const labelList& sensitivityPatchIDs
1209)
1210{
1211 // Return field
1212 vectorField controlPointDerivs(cps_.size(), Zero);
1213
1214 // Get parametric coordinates
1215 const vectorField& parametricCoordinates = getParametricCoordinates();
1216
1217 // Auxiliary quantities
1219 const labelList& reverseMap = reverseMapPtr_();
1220
1221 forAll(controlPointDerivs, cpI)
1222 {
1223 forAll(sensitivityPatchIDs, pI)
1224 {
1225 const label patchI = sensitivityPatchIDs[pI];
1226 const polyPatch& patch = mesh_.boundaryMesh()[patchI];
1227 const label patchStart = patch.start();
1228 const fvPatchVectorField& patchSens = faceSens[patchI];
1229
1230 // loop over patch faces
1231 forAll(patch, fI)
1232 {
1233 const face& fGlobal = mesh_.faces()[fI + patchStart];
1234 const pointField facePoints = fGlobal.points(mesh_.points());
1235 // loop over face points
1236 tensorField facePointDerivs(facePoints.size(), Zero);
1237 forAll(fGlobal, pI)
1238 {
1239 const label globalIndex = fGlobal[pI]; //global point index
1240 const label whichPointInBox = reverseMap[globalIndex];
1241 // if point resides within control points box,
1242 // add contribution to d( facePoints )/db
1243 if (whichPointInBox != -1)
1244 {
1245 // TENSOR-BASED
1246 //~~~~~~~~~~~~~
1247 facePointDerivs[pI] =
1248 transformationTensorDxDb(globalIndex)
1249 * volumeDerivativeCP
1250 (
1251 parametricCoordinates[globalIndex],
1252 cpI
1253 );
1254
1255 }
1256 }
1257
1258 tensor fCtrs_d =
1260 (
1261 facePoints,
1262 facePointDerivs
1263 )[0];
1264 controlPointDerivs[cpI] += patchSens[fI] & fCtrs_d;
1265 }
1266 }
1267 }
1268 // Sum contributions from all processors
1269 Pstream::listReduce(controlPointDerivs, sumOp<vector>());
1270
1271 return controlPointDerivs;
1272}
1273
1274
1276(
1277 const vectorField& faceSens,
1278 const label patchI,
1279 const label cpI
1280)
1281{
1282 // Return vector
1283 vector cpSens(Zero);
1284 // Get parametric coordinates
1285 const vectorField& parametricCoordinates = getParametricCoordinates();
1286
1287 // Auxiliary quantities
1289 const labelList& reverseMap = reverseMapPtr_();
1290
1291 const polyPatch& patch = mesh_.boundaryMesh()[patchI];
1292 const label patchStart = patch.start();
1293 // Loop over patch faces
1294 forAll(patch, fI)
1295 {
1296 const face& fGlobal = mesh_.faces()[fI + patchStart];
1297 const pointField facePoints = fGlobal.points(mesh_.points());
1298 // Loop over face points
1299 tensorField facePointDerivs(facePoints.size(), Zero);
1300 forAll(fGlobal, pI)
1301 {
1302 const label globalIndex = fGlobal[pI]; //global point index
1303 const label whichPointInBox = reverseMap[globalIndex];
1304 // If point resides within control points box,
1305 // add contribution to d( facePoints )/db
1306 if (whichPointInBox != -1)
1307 {
1308 // TENSOR-BASED
1309 //~~~~~~~~~~~~~
1310 facePointDerivs[pI] =
1311 transformationTensorDxDb(globalIndex)
1312 *volumeDerivativeCP
1313 (
1314 parametricCoordinates[globalIndex],
1315 cpI
1316 );
1317 }
1318 }
1319
1320 tensor fCtrs_d =
1322 (
1323 facePoints,
1324 facePointDerivs
1325 )[0];
1326 cpSens += faceSens[fI] & fCtrs_d;
1327 }
1328 // Sum contributions from all processors
1329 reduce(cpSens, sumOp<vector>());
1330
1331 return cpSens;
1332}
1333
1334
1336(
1337 const label patchI,
1338 const label cpI,
1339 bool DimensionedNormalSens
1340)
1341{
1342 const fvPatch& patch = mesh_.boundary()[patchI];
1343 const polyPatch& ppatch = patch.patch();
1344 // Return field
1345 auto tdndbSens = tmp<tensorField>::New(patch.size(), Zero);
1346 auto& dndbSens = tdndbSens.ref();
1347
1348 // Auxiliary quantities
1350 const label patchStart = ppatch.start();
1351 const labelList& reverseMap = reverseMapPtr_();
1352
1353 // Get parametric coordinates
1354 const vectorField& parametricCoordinates = getParametricCoordinates();
1355
1356 // Loop over patch faces
1357 forAll(patch, fI)
1358 {
1359 const face& fGlobal = mesh_.faces()[fI + patchStart];
1360 const pointField facePoints = fGlobal.points(mesh_.points());
1361 // Loop over face points
1362 tensorField facePointDerivs(facePoints.size(), Zero);
1363 forAll(fGlobal, pI)
1364 {
1365 const label globalIndex = fGlobal[pI]; //global point index
1366 const label whichPointInBox = reverseMap[globalIndex];
1367 // If point resides within control points box,
1368 // add contribution to d( facePoints )/db
1369 if (whichPointInBox != -1)
1370 {
1371 // TENSOR-BASED
1372 //~~~~~~~~~~~~~
1373 facePointDerivs[pI] =
1374 transformationTensorDxDb(globalIndex)
1375 *volumeDerivativeCP
1376 (
1377 parametricCoordinates[globalIndex],
1378 cpI
1379 );
1380 }
1381 }
1382
1383 // Determine whether to return variance of dimensioned or unit normal
1384 tensorField dNdbSens =
1386 (
1387 facePoints,
1388 facePointDerivs
1389 );
1390
1391 if (DimensionedNormalSens)
1392 {
1393 dndbSens[fI] = dNdbSens[1];
1394 }
1395 else
1396 {
1397 dndbSens[fI] = dNdbSens[2];
1399 }
1400
1401 return tdndbSens;
1402}
1403
1404
1405Foam::tmp<Foam::tensorField> Foam::NURBS3DVolume::patchDxDb
1406(
1407 const label patchI,
1408 const label cpI
1409)
1410{
1411 // Get parametric coordinates
1412 const vectorField& parametricCoordinates = getParametricCoordinates();
1413
1414 // Patch data
1415 const polyPatch& patch = mesh_.boundaryMesh()[patchI];
1416 const labelList& meshPoints = patch.meshPoints();
1417
1418 // Return field
1419 auto tdxdb = tmp<tensorField>::New(patch.nPoints(), Zero);
1420 auto& dxdb = tdxdb.ref();
1421
1422 forAll(meshPoints, pI)
1423 {
1424 const label globalIndex = meshPoints[pI]; //global point index
1425 const label whichPointInBox = reverseMapPtr_()[globalIndex];
1426
1427 // If point resides within control points box, find dxdb
1428 if (whichPointInBox != -1)
1429 {
1430 dxdb[pI] =
1431 transformationTensorDxDb(globalIndex)
1432 *volumeDerivativeCP
1433 (
1434 parametricCoordinates[globalIndex],
1435 cpI
1436 );
1438 }
1439
1440 return tdxdb;
1441}
1442
1443
1444Foam::tmp<Foam::tensorField> Foam::NURBS3DVolume::patchDxDbFace
1445(
1446 const label patchI,
1447 const label cpI
1448)
1449{
1450 // get parametric coordinates
1451 const vectorField& parametricCoordinates = getParametricCoordinates();
1452
1453 // Patch data
1454 const polyPatch& patch = mesh_.boundaryMesh()[patchI];
1455 const label patchStart = patch.start();
1456
1457 // Return field
1458 auto tdxdb = tmp<tensorField>::New(patch.size(), Zero);
1459 auto& dxdb = tdxdb.ref();
1460
1461 // Mesh differentiation engine
1462 deltaBoundary deltaBound(mesh_);
1463
1464 forAll(patch, fI)
1465 {
1466 const face& fGlobal = mesh_.faces()[fI + patchStart];
1467 const pointField facePoints = fGlobal.points(mesh_.points());
1468 // Loop over face points
1469 tensorField facePointDerivs(facePoints.size(), Zero);
1470 forAll(fGlobal, pI)
1471 {
1472 const label globalIndex = fGlobal[pI]; //global point index
1473 const label whichPointInBox = reverseMapPtr_()[globalIndex];
1474 // If point resides within control points box,
1475 // add contribution to d( facePoints )/db
1476 if (whichPointInBox != -1)
1477 {
1478 // TENSOR-BASED
1479 //~~~~~~~~~~~~~
1480 facePointDerivs[pI] =
1481 transformationTensorDxDb(globalIndex)
1482 *volumeDerivativeCP
1483 (
1484 parametricCoordinates[globalIndex],
1485 cpI
1486 );
1487
1488 }
1489 }
1490 dxdb[fI] =
1491 deltaBound.makeFaceCentresAndAreas_d
1492 (
1493 facePoints,
1494 facePointDerivs
1495 )[0];
1496 }
1497
1498 return tdxdb;
1499}
1500
1501
1503(
1504 const vectorField& uVector
1505) const
1506{
1507 const label nPoints = mapPtr_().size();
1508 auto tpoints = tmp<vectorField>::New(nPoints, Zero);
1509 auto& points = tpoints.ref();
1510
1511 forAll(points, pI)
1512 {
1513 const label globalPI = mapPtr_()[pI];
1514 points[pI] = coordinates(uVector[globalPI]);
1515 }
1516
1517 return tpoints;
1518}
1519
1520
1522(
1523 const vector& uVector
1524) const
1525{
1526 const label degreeU = basisU_.degree();
1527 const label degreeV = basisV_.degree();
1528 const label degreeW = basisW_.degree();
1529
1530 const label nCPsU = basisU_.nCPs();
1531 const label nCPsV = basisV_.nCPs();
1532 const label nCPsW = basisW_.nCPs();
1533
1534 const scalar u = uVector.x();
1535 const scalar v = uVector.y();
1536 const scalar w = uVector.z();
1537
1538 vector point(Zero);
1539 for (label iCPw = 0; iCPw < nCPsW; iCPw++)
1540 {
1541 const scalar basisW(basisW_.basisValue(iCPw, degreeW, w));
1542 for (label iCPv = 0; iCPv < nCPsV; iCPv++)
1543 {
1544 const scalar basisVW =
1545 basisW*basisV_.basisValue(iCPv, degreeV, v);
1546 for (label iCPu = 0; iCPu < nCPsU; iCPu++)
1547 {
1548 point +=
1549 cps_[getCPID(iCPu, iCPv, iCPw)]
1550 *basisU_.basisValue(iCPu, degreeU, u)
1551 *basisVW;
1552 }
1554 }
1555
1556 return point;
1557}
1558
1559
1560Foam::tmp<Foam::vectorField> Foam::NURBS3DVolume::computeNewPoints
1561(
1562 const vectorField& controlPointsMovement
1563)
1564{
1565 // Get parametric coordinates and map
1566 const vectorField& paramCoors = getParametricCoordinates();
1567 const labelList& map = mapPtr_();
1568
1569 // Update control points position
1570 cps_ += controlPointsMovement;
1571 writeCps("cpsBsplines"+mesh_.time().timeName());
1572
1573 // Compute new mesh points based on updated control points
1574 tmp<vectorField> tparameterizedPoints = coordinates(paramCoors);
1575 const vectorField& parameterizedPoints = tparameterizedPoints();
1576
1577 // Return field. Initialized with current mesh points
1578 auto tnewPoints = tmp<vectorField>::New(mesh_.points());
1579 auto& newPoints = tnewPoints.ref();
1580
1581 // Update position of parameterized points
1582 forAll(parameterizedPoints, pI)
1583 {
1584 newPoints[map[pI]] = transformPointToCartesian(parameterizedPoints[pI]);
1585 }
1586
1587 // Update coordinates in the local system based on the cartesian points
1588 updateLocalCoordinateSystem(newPoints);
1589 DebugInfo
1590 << "Max mesh movement equal to "
1591 << gMax(mag(newPoints - mesh_.points())) << endl;
1592
1593 return tnewPoints;
1594}
1595
1596
1598(
1599 const vectorField& controlPointsMovement,
1600 const labelList& patchesToBeMoved,
1601 const bool updateCPs
1602)
1603{
1604 // Get parametric coordinates
1605 const vectorField& paramCoors = getParametricCoordinates();
1606
1607 // Update control points position
1608 cps_ += controlPointsMovement;
1609
1610 if (updateCPs)
1611 {
1612 writeCps("cpsBsplines"+mesh_.time().timeName());
1613 }
1614
1615 // Return field. Initialized with current mesh points
1616 auto tnewPoints = tmp<vectorField>::New(mesh_.points());
1617 auto& newPoints = tnewPoints.ref();
1618
1619 // Update position of parameterized boundary points
1620 for (const label patchI : patchesToBeMoved)
1621 {
1622 const polyPatch& patch = mesh_.boundaryMesh()[patchI];
1623 const labelList& meshPoints = patch.meshPoints();
1624
1625 for (const label globalIndex : meshPoints)
1626 {
1627 const label whichPointInBox = reverseMapPtr_()[globalIndex];
1628 // If point resides within control points box,
1629 // compute new cartesian coordinates
1630 if (whichPointInBox != -1)
1631 {
1632 newPoints[globalIndex] =
1633 transformPointToCartesian
1634 (
1636 (
1637 paramCoors[globalIndex]
1638 )
1639 );
1640 }
1641 }
1642 }
1643
1644 if (updateCPs)
1645 {
1646 // Update coordinates in the local system based on the cartesian points
1647 updateLocalCoordinateSystem(newPoints);
1648 }
1649 else
1650 {
1651 // Move control points to their initial position
1652 cps_ -= controlPointsMovement;
1653 }
1654
1655 DebugInfo
1656 << "Max mesh movement equal to "
1657 << gMax(mag(newPoints - mesh_.points())) << endl;
1658
1659 return tnewPoints;
1660}
1661
1662
1664(
1665 const label i,
1666 const label j,
1667 const label k
1668) const
1669{
1670 const label nCPsU = basisU_.nCPs();
1671 const label nCPsV = basisV_.nCPs();
1672
1673 return k*nCPsU*nCPsV + j*nCPsU + i;
1674}
1675
1676
1678(
1679 label& i,
1680 label& j,
1681 label& k,
1682 const label cpID
1683) const
1684{
1685 const label nCPsU = basisU_.nCPs();
1686 const label nCPsV = basisV_.nCPs();
1687 k = cpID/(nCPsU*nCPsV);
1688 const label inKplaneID = (cpID - k*(nCPsU*nCPsV));
1689 j = inKplaneID/nCPsU;
1690 i = inKplaneID - j*nCPsU;
1691}
1692
1693
1695{
1696 if (cps_.size() != newCps.size())
1697 {
1699 << "Attempting to replace control points with a set of "
1700 << "different size"
1701 << exit(FatalError);
1702 }
1703 cps_ = newCps;
1704}
1705
1706
1708(
1709 vectorField& controlPointsMovement
1710) const
1711{
1712 forAll(controlPointsMovement, cpI)
1713 {
1714 if (!activeDesignVariables_[3*cpI])
1715 {
1716 controlPointsMovement[cpI].x() = Zero;
1717 }
1718 if (!activeDesignVariables_[3*cpI + 1])
1719 {
1720 controlPointsMovement[cpI].y() = Zero;
1721 }
1722 if (!activeDesignVariables_[3*cpI + 2])
1724 controlPointsMovement[cpI].z() = Zero;
1725 }
1726 }
1727}
1728
1729
1731(
1732 const vectorField& controlPointsMovement,
1733 const labelList& patchesToBeMoved
1734)
1735{
1736 // Backup old cps
1737 vectorField oldCPs = cps_;
1738 // Get parametric coordinates
1739 const vectorField& paramCoors = getParametricCoordinates();
1740 // Update control points position
1741 cps_ += controlPointsMovement;
1742 // Update position of parameterized boundary points
1743 scalar maxDisplacement(Zero);
1744 for (const label patchI : patchesToBeMoved)
1745 {
1746 const polyPatch& patch = mesh_.boundaryMesh()[patchI];
1747 const labelList& meshPoints = patch.meshPoints();
1748
1749 for (const label globalIndex : meshPoints)
1750 {
1751 const label whichPointInBox = reverseMapPtr_()[globalIndex];
1752 // If point resides within control points box,
1753 // compute new cartesian coordinates
1754 if (whichPointInBox != -1)
1755 {
1756 vector newPoint =
1757 transformPointToCartesian
1758 (
1760 (
1761 paramCoors[globalIndex]
1762 )
1763 );
1764 maxDisplacement =
1765 max
1766 (
1767 maxDisplacement,
1768 mag(newPoint - mesh_.points()[globalIndex])
1769 );
1770 }
1771 }
1772 }
1773 reduce(maxDisplacement, maxOp<scalar>());
1774 cps_ = oldCPs;
1775
1776 return maxDisplacement;
1777}
1778
1779
1781{
1782 if (!mapPtr_)
1783 {
1784 findPointsInBox(localSystemCoordinates_);
1785 }
1786 tmp<vectorField> pointsInBox
1787 (
1789 );
1790
1791 return pointsInBox;
1792}
1793
1794
1796{
1797 if (!mapPtr_)
1798 {
1800 }
1801
1802 return mapPtr_();
1803}
1804
1805
1807{
1808 if (!reverseMapPtr_)
1809 {
1811 }
1812
1813 return reverseMapPtr_();
1814}
1815
1816
1818{
1819 // If not computed yet, compute parametric coordinates
1820 if (!parametricCoordinatesPtr_)
1821 {
1822 // Find mesh points inside control points box
1823 // if they have been identified yet
1824 if (!mapPtr_)
1825 {
1826 findPointsInBox(localSystemCoordinates_);
1827 }
1829 }
1830
1831 return parametricCoordinatesPtr_();
1832}
1833
1834
1836{
1837 // Get parametric coordinates
1838 const vectorField& parametricCoordinates = getParametricCoordinates();
1839
1840 // Set return field to zero
1841 tmp<pointTensorField> tDxDb
1842 (
1844 (
1845 IOobject
1846 (
1847 "DxDb",
1848 mesh_.time().timeName(),
1849 mesh_,
1852 ),
1853 pointMesh::New(mesh_),
1855 )
1856 );
1857
1858 pointTensorField& DxDb = tDxDb.ref();
1859
1860 // All points outside the control box remain unmoved.
1861 // Loop over only points within the control box
1862 const labelList& map = mapPtr_();
1863 for (const label globalIndex : map)
1864 {
1865 DxDb[globalIndex] =
1866 transformationTensorDxDb(globalIndex)
1867 *volumeDerivativeCP
1868 (
1869 parametricCoordinates[globalIndex],
1870 cpI
1872 }
1873
1874 return tDxDb;
1875}
1876
1877
1879(
1880 const label cpI
1881)
1882{
1883 // Get parametric coordinates
1884 const vectorField& parametricCoordinates = getParametricCoordinates();
1885
1886 // Set return field to zero
1887 auto tDxDb = volTensorField::New
1888 (
1889 "DxDb",
1891 mesh_,
1893 );
1894
1895 volTensorField& DxDb = tDxDb.ref();
1896 deltaBoundary deltaBound(mesh_);
1897 const labelListList& pointCells = mesh_.pointCells();
1898
1899 // All points outside the control box remain unmoved.
1900 // Loop over only points within the control box
1901 const labelList& map = mapPtr_();
1902 for (const label globalIndex : map)
1903 {
1904 tensor pointDxDb =
1905 transformationTensorDxDb(globalIndex)
1906 *volumeDerivativeCP
1907 (
1908 parametricCoordinates[globalIndex],
1909 cpI
1910 );
1911 const labelList& pointCellsI = pointCells[globalIndex];
1912 tmp<tensorField> tC_d = deltaBound.cellCenters_d(globalIndex);
1913 const tensorField& C_d = tC_d();
1914 forAll(pointCellsI, cI)
1915 {
1916 const label cellI = pointCellsI[cI];
1917 DxDb[cellI] += C_d[cI] & pointDxDb;
1918 }
1919 }
1920
1921 // Assign boundary values since the grad of this field is often needed
1922 forAll(mesh_.boundary(), pI)
1923 {
1924 const fvPatch& patch = mesh_.boundary()[pI];
1925 if (!isA<coupledFvPatch>(patch))
1926 {
1927 DxDb.boundaryFieldRef()[pI] = patchDxDbFace(pI, cpI);
1928 }
1929 }
1930
1931 // Correct coupled boundaries
1932 DxDb.correctBoundaryConditions();
1933
1934 return tDxDb;
1935}
1936
1937
1939{
1940 label nU(basisU_.nCPs());
1941 return label(nU % 2 == 0 ? 0.5*nU : 0.5*(nU - 1) + 1);
1942}
1943
1944
1946{
1947 label nV(basisV_.nCPs());
1948 return label(nV % 2 == 0 ? 0.5*nV : 0.5*(nV - 1) + 1);
1949}
1950
1951
1953{
1954 label nW(basisW_.nCPs());
1955 return label(nW % 2 == 0 ? 0.5*nW : 0.5*(nW - 1) + 1);
1956}
1957
1960{
1962}
1963
1964
1966(
1967 const fileName& baseName,
1968 const bool transform
1969) const
1970{
1971 const label nCPsU = basisU_.nCPs();
1972 const label nCPsV = basisV_.nCPs();
1973
1974 vectorField cpsInCartesian(cps_);
1975 if (transform)
1976 {
1977 forAll(cpsInCartesian, cpI)
1978 {
1979 cpsInCartesian[cpI] = transformPointToCartesian(cps_[cpI]);
1980 }
1981 }
1982
1983 Info<< "Writing control point positions to file" << endl;
1984
1985 if (Pstream::master())
1986 {
1987 OFstream cpsFile("optimisation"/cpsFolder_/name_ + baseName + ".csv");
1988 // Write header
1989 cpsFile
1990 << "\"Points : 0\", \"Points : 1\", \"Points : 2\","
1991 << "\"i\", \"j\", \"k\","
1992 << "\"active : 0\", \"active : 1\", \"active : 2\"" << endl;
1993
1994 forAll(cpsInCartesian, cpI)
1995 {
1996 const label iCPw = cpI/label(nCPsU*nCPsV);
1997 const label iCPv = (cpI - iCPw*nCPsU*nCPsV)/nCPsU;
1998 const label iCPu = (cpI - iCPw*nCPsU*nCPsV - iCPv*nCPsU);
1999
2000 cpsFile
2001 << cpsInCartesian[cpI].x() << ", "
2002 << cpsInCartesian[cpI].y() << ", "
2003 << cpsInCartesian[cpI].z() << ", "
2004 << iCPu << ", "
2005 << iCPv << ", "
2006 << iCPw << ", "
2007 << activeDesignVariables_[3*cpI] << ", "
2008 << activeDesignVariables_[3*cpI + 1] << ", "
2009 << activeDesignVariables_[3*cpI + 2] << endl;
2010 }
2011 }
2012}
2013
2016{
2017 parametricCoordinatesPtr_().write();
2018}
2019
2020
2022{
2023 cps_.writeEntry("controlPoints", os);
2024 return true;
2025}
2026
2027
2028// ************************************************************************* //
scalar maxValue
scalar minValue
label k
bool found
Macros for easy insertion into run-time selection tables.
tmp< Field< cmptType > > component(const direction) const
Return a component field of the field.
Definition Field.C:607
static tmp< GeometricField< tensor, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< tensor >::calculatedType())
Internal & ref(const bool updateAccessTime=true)
Same as internalFieldRef().
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
void correctBoundaryConditions()
Correct boundary field.
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
@ NO_REGISTER
Do not request registration (bool: false).
@ NO_READ
Nothing to be read.
@ READ_IF_PRESENT
Reading is optional [identical to LAZY_READ].
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const Time & time() const noexcept
Return Time associated with the objectRegistry.
Definition IOobject.C:456
bool typeHeaderOk(const bool checkType=true, const bool search=true, const bool verbose=true)
Read header (respects is_globalIOobject trait) and check its info. A void type suppresses trait and t...
void setSize(label n)
Alias for resize().
Definition List.H:536
static FOAM_NO_DANGLING_REFERENCE const pointMesh & New(const polyMesh &mesh, Args &&... args)
NURBS3DVolume morpher. Includes support functions for gradient computations Base class providing supp...
bool readStoredData_
Read parametric coordinates from file if present in the folder.
void getIJK(label &i, label &j, label &k, const label cpID) const
Get I-J-K coordinates from control point ID.
label maxIter_
Max iterations for Newton-Raphson.
virtual void updateLocalCoordinateSystem(const vectorField &cartesianPoints)=0
Update coordinates in the local system based on the cartesian points.
boolVectorList confineWMinCPs_
static autoPtr< NURBS3DVolume > New(const dictionary &dict, const fvMesh &mesh, bool computeParamCoors=true)
Return a reference to the selected NURBS model.
boolVectorList confineVMinCPs_
tmp< tensorField > patchDxDbFace(const label patchI, const label cpI)
Get patch dx/db.
void computeParametricCoordinates(const vectorField &points)
Compute parametric coordinates in order to match a given set of coordinates, based on the cps of the ...
label nUSymmetry() const
Get number of variables if CPs are moved symmetrically in U.
void confineInertControlPoints()
Deactivate control points if they affect no mesh point.
const fvMesh & mesh_
autoPtr< labelList > mapPtr_
Map of points-in-box to mesh points.
vector volumeDerivativeU(const scalar u, const scalar v, const scalar w) const
Volume point derivative wrt u at point u,v,w.
void confineBoundaryControlPoints()
Confine movement in boundary control points if necessary.
const labelList & getMap()
Get map of points in box to mesh points.
scalar tolerance_
Tolerance for Newton-Raphson.
virtual bool writeData(Ostream &os) const
Write the control points to support restart.
scalar computeMaxBoundaryDisplacement(const vectorField &controlPointsMovement, const labelList &patchesToBeMoved)
Compute max. displacement at the boundary.
boolVectorList confineVMaxCPs_
void findPointsInBox(const vectorField &meshPoints)
Find points within control points box.
tmp< vectorField > computeNewBoundaryPoints(const vectorField &controlPointsMovement, const labelList &patchesToBeMoved, const bool moveCPs=true)
Boundary mesh movement based on given control point movement.
boolVectorList confineWMaxCPs_
vectorField localSystemCoordinates_
Coordinates in the local system for which CPs are defined.
boolVectorList confineUMaxCPs_
vectorField computeControlPointSensitivities(const pointVectorField &pointSens, const labelList &sensitivityPatchIDs)
Control point sensitivities computed using point-based surface sensitivities.
boolVectorList confineUMinCPs_
Which movement components to freeze in each plane.
const pointVectorField & getParametricCoordinates()
Get parametric coordinates.
tmp< volTensorField > getDxCellsDb(const label cpI)
Get dxCartesiandb for a certain control point on cells.
virtual tensor transformationTensorDxDb(label globalPointIndex)=0
Transformation tensor for dxdb, from local coordinate system to cartesian.
autoPtr< labelList > reverseMapPtr_
Map of mesh points to points-in-box.
label nVSymmetry() const
Get number of variables if CPs are moved symmetrically in V.
tmp< pointTensorField > getDxDb(const label cpI)
Get dxCartesiandb for a certain control point.
bool bound(vector &vec, scalar minValue=1e-7, scalar maxValue=0.999999) const
Bound components to certain limits.
void makeFolders()
Create folders to store cps and derivatives.
void writeCps(const fileName &="cpsFile", const bool transform=true) const
Write control points on a cartesian coordinates system for visualization.
NURBS3DVolume(const dictionary &dict, const fvMesh &mesh, bool computeParamCoors=true)
Construct from dictionary.
boolList activeDesignVariables_
Which design variables are changed in an optimisation.
tmp< tensorField > patchDxDb(const label patchI, const label cpI)
Get patch dx/db.
tmp< vectorField > getPointsInBox()
Get mesh points that reside within the control points box.
vector volumeDerivativeV(const scalar u, const scalar v, const scalar w) const
Volume point derivative wrt v at point u,v,w.
virtual vector transformPointToCartesian(const vector &localCoordinates) const =0
Transform a point from its coordinate system to a cartesian system.
tmp< tensorField > dndbBasedSensitivities(const label patchI, const label cpI, bool DimensionedNormalSens=true)
Part of control point sensitivities related to the face normal variations.
string cpsFolder_
Folder to store control points.
NURBSbasis basisU_
NURBS basis functions.
label confineUMovement_
Confine movement in certain directions and control points. Refers to the local system.
void confineControlPoint(const label cpI)
Confine all three movements for a prescribed control point.
const dictionary & dict() const
Get dictionary.
Vector< label > nSymmetry() const
Get number of variables per direction, if CPs are moved symmetrically.
vectorField cps_
The volumetric B-Splines control points.
void boundControlPointMovement(vectorField &controlPointsMovement) const
Bound control points movement in the boundary control points and in certain directions if needed.
void continuityRealatedConfinement()
Confine control point movement to maintain user-defined continuity.
label getCPID(const label i, const label j, const label k) const
Get control point ID from its I-J-K coordinates.
autoPtr< pointVectorField > parametricCoordinatesPtr_
Parametric coordinates of pointsInBox.
scalar volumeDerivativeCP(const vector &u, const label cpI) const
Volume point derivative wrt to control point cpI at point u,v,w.
label confineBoundaryControlPoints_
label nMaxBound_
How many times to bound parametric coordinates until deciding it is outside the box.
const NURBSbasis & basisW() const
const fvMesh & mesh() const
Get mesh.
label nWSymmetry() const
Get number of variables if CPs are moved symmetrically in W.
void writeParamCoordinates() const
Write parametric coordinates.
void determineActiveDesignVariablesAndPoints()
Create lists with active design variables and control points.
tmp< vectorField > coordinates(const vectorField &uVector) const
Compute cartesian coordinates based on control points and parametric coordinates.
tmp< vectorField > computeNewPoints(const vectorField &controlPointsMovement)
Mesh movement based on given control point movement.
tensor JacobianUVW(const vector &u) const
Jacobian matrix wrt to the volume parametric coordinates.
void setControlPoints(const vectorField &newCps)
Set new control points.
boolList activeControlPoints_
Which of the cps are moved in an optimisation.
const labelList & getReverseMap()
Get map of mesh points to points in box. Return -1 if point is outside the box.
vector volumeDerivativeW(const scalar u, const scalar v, const scalar w) const
Volume point derivative wrt w at point u,v,w.
void confineControlPointsDirections()
Confine movement in all control points for user-defined directions.
Output to file stream as an OSstream, normally using std::ofstream for the actual output.
Definition OFstream.H:75
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
static void listReduce(UList< T > &values, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce list elements (list must be equal size on all ranks), applying bop to each list element.
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition UPstream.H:1714
static void reduceOr(bool &value, const int communicator=worldComm)
Logical (or) reduction (MPI_AllReduce).
static const Form zero
const Cmpt & component(const direction) const
Templated 3D Vector derived from VectorSpace adding construction from 3 components,...
Definition Vector.H:61
const Cmpt & x() const noexcept
Access to the vector x component.
Definition Vector.H:135
const Cmpt & z() const noexcept
Access to the vector z component.
Definition Vector.H:145
const Cmpt & y() const noexcept
Access to the vector y component.
Definition Vector.H:140
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
Specialized bundling of boolean values as a vector of 3 components, element access using x(),...
Definition boolVector.H:55
bool z() const noexcept
The z component.
Definition boolVector.H:194
bool x() const noexcept
The x component.
Definition boolVector.H:184
bool y() const noexcept
The y component.
Definition boolVector.H:189
static autoPtr< controlPointsDefinition > New(NURBS3DVolume &box)
Return a reference to the selected controlPointsDefinition model.
Differentiation of the mesh data structure.
tmp< tensorField > cellCenters_d(const label pointI)
Compute the change of the cell centers of the pointCells of pointI, for a unitary movement of pointI ...
vectorField makeFaceCentresAndAreas_d(const pointField &p, const pointField &p_d)
Given a face and the points to be moved in the normal direction, find faceArea, faceCentre and unitVe...
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
T getOrDefaultCompat(const word &keyword, std::initializer_list< std::pair< const char *, int > > compat, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value using any compatibility names if needed.
static const dictionary null
An empty dictionary, which is also the parent for all dictionaries.
Definition dictionary.H:487
word dictName() const
The local dictionary name (final part of scoped name).
Definition dictionaryI.H:53
A face is a list of labels corresponding to mesh vertices.
Definition face.H:71
pointField points(const UList< point > &pts) const
Return the points corresponding to this face.
Definition faceI.H:80
A class for handling file names.
Definition fileName.H:75
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
Calculates a non-overlapping list of offsets based on an input size (eg, number of cells) from differ...
Definition globalIndex.H:77
localIOdictionary(const IOobject &io, const dictionary *fallback=nullptr)
Construct given an IOobject and optional fallback dictionary content.
Smooth ATC in cells having a point to a set of patches supplied by type.
Definition pointCells.H:55
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
label start() const noexcept
Return start label of this patch in the polyMesh face list.
Definition polyPatch.H:446
Tensor of scalars, i.e. Tensor<scalar>.
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition tmpI.H:235
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
dynamicFvMesh & mesh
PtrList< coordinateSystem > coordinates(solidRegions.size())
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition error.H:637
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
const word dictName("faMeshDefinition")
const pointField & points
label nPoints
word timeName
Definition getTimeIndex.H:3
#define DebugInfo
Report an information message using Foam::Info.
#define WarningInFunction
Report a warning using Foam::Warning.
const std::string patch
OpenFOAM patch number as a std::string.
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimless
Dimensionless.
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
GeometricField< tensor, pointPatchField, pointMesh > pointTensorField
List< label > labelList
A List of labels.
Definition List.H:62
bool mkDir(const fileName &pathName, mode_t mode=0777)
Make a directory and return an error if it could not be created.
Definition POSIX.C:616
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition bound.C:29
refinementData transform(const tensor &, const refinementData val)
No-op rotational transform for base types.
messageStream Info
Information stream (stdout output on master, null elsewhere).
GeometricField< vector, pointPatchField, pointMesh > pointVectorField
Tensor< scalar > tensor
Definition symmTensor.H:57
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
GeometricField< tensor, fvPatchField, volMesh > volTensorField
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
scalar diff(const triad &A, const triad &B)
Return a quantity of the difference between two triads.
Definition triad.C:373
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
Field< vector > vectorField
Specialisation of Field<T> for vector.
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
vector point
Point is a vector.
Definition point.H:37
List< bool > boolList
A List of bools.
Definition List.H:60
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
Field< tensor > tensorField
Specialisation of Field<T> for tensor.
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
void cmptMag(FieldField< Field, Type > &cf, const FieldField< Field, Type > &f)
vectorField pointField
pointField is a vectorField.
Vector< scalar > vector
Definition vector.H:57
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
volVectorField::Boundary boundaryVectorField
Type gMax(const FieldField< Field, Type > &f)
fvPatchField< vector > fvPatchVectorField
dimensioned< tensor > dimensionedTensor
Dimensioned tensor obtained from generic dimensioned type.
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299