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-2022 PCOpt/NTUA
9 Copyright (C) 2013-2022 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"
36#include "coupledFvPatch.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);
97
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 (
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
341)
342{
343 const vectorField& points = tPoints();
344 computeParametricCoordinates(points);
345}
346
347
349(
350 vector& vec,
351 scalar minValue,
352 scalar maxValue
353)
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{
394 if (Pstream::master())
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 }
429 }
430}
431
432
434{
435 const label nCPsU = basisU_.nCPs();
436 const label nCPsV = basisV_.nCPs();
437 const label nCPsW = basisW_.nCPs();
438
439 // Zero movement of the boundary control points. Active by default
440 if (confineBoundaryControlPoints_)
441 {
442 // Side patches
443 for (label iCPw = 0; iCPw < nCPsW; iCPw += nCPsW - 1)
444 {
445 for (label iCPv = 0; iCPv < nCPsV; iCPv++)
446 {
447 for (label iCPu = 0; iCPu < nCPsU; iCPu++)
448 {
449 confineControlPoint(getCPID(iCPu, iCPv, iCPw));
450 }
451 }
452 }
453 // Front-back patches
454 for (label iCPw = 0; iCPw < nCPsW; iCPw++)
455 {
456 for (label iCPv = 0; iCPv < nCPsV; iCPv++)
457 {
458 for (label iCPu = 0; iCPu < nCPsU; iCPu += nCPsU - 1)
459 {
460 confineControlPoint(getCPID(iCPu, iCPv, iCPw));
461 }
462 }
463 }
464 // Top-bottom patches
465 for (label iCPw = 0; iCPw < nCPsW; iCPw++)
466 {
467 for (label iCPv = 0; iCPv < nCPsV; iCPv += nCPsV - 1)
468 {
469 for (label iCPu = 0; iCPu < nCPsU; iCPu++)
470 {
471 confineControlPoint(getCPID(iCPu, iCPv, iCPw));
472 }
473 }
474 }
475 }
476}
477
478
480{
481 const label nCPsU = basisU_.nCPs();
482 const label nCPsV = basisV_.nCPs();
483 const label nCPsW = basisW_.nCPs();
484
485 // Zero movement to a number of x-constant slices of cps in order to
486 // preserve continuity at the boundary of the parameterized space
487 forAll(confineUMinCPs_, iCPu)
488 {
489 const boolVector& confineSlice = confineUMinCPs_[iCPu];
490 // Control points at the start of the parameterized space
491 for (label iCPw = 0; iCPw < nCPsW; iCPw++)
492 {
493 for (label iCPv = 0; iCPv < nCPsV; iCPv++)
494 {
495 confineControlPoint(getCPID(iCPu, iCPv, iCPw), confineSlice);
496 }
497 }
498 }
499
500 forAll(confineUMaxCPs_, sliceI)
501 {
502 const boolVector& confineSlice = confineUMaxCPs_[sliceI];
503 label iCPu = nCPsU - 1 - sliceI;
504 // Control points at the end of the parameterized space
505 for (label iCPw = 0; iCPw < nCPsW; iCPw++)
506 {
507 for (label iCPv = 0; iCPv < nCPsV; iCPv++)
508 {
509 confineControlPoint(getCPID(iCPu, iCPv, iCPw), confineSlice);
510 }
511 }
512 }
513
514 // Zero movement to a number of y-constant slices of cps in order to
515 // preserve continuity at the boundary of the parameterized space
516 forAll(confineVMinCPs_, iCPv)
517 {
518 const boolVector& confineSlice = confineVMinCPs_[iCPv];
519 // Control points at the start of the parameterized space
520 for (label iCPw = 0; iCPw < nCPsW; iCPw++)
521 {
522 for (label iCPu = 0; iCPu < nCPsU; iCPu++)
523 {
524 confineControlPoint(getCPID(iCPu, iCPv, iCPw), confineSlice);
525 }
526 }
527 }
528
529 forAll(confineVMaxCPs_, sliceI)
530 {
531 const boolVector& confineSlice = confineVMaxCPs_[sliceI];
532 label iCPv = nCPsV - 1 - sliceI;
533 // Control points at the end of the parameterized space
534 for (label iCPw = 0; iCPw < nCPsW; iCPw++)
535 {
536 for (label iCPu = 0; iCPu < nCPsU; iCPu++)
537 {
538 confineControlPoint(getCPID(iCPu, iCPv, iCPw), confineSlice);
539 }
540 }
541 }
542
543 // Zero movement to a number of w-constant slices of cps in order to
544 // preserve continuity at the boundary of the parameterized space
545 forAll(confineWMinCPs_, iCPw)
546 {
547 const boolVector& confineSlice = confineWMinCPs_[iCPw];
548 // Control points at the start of the parameterized space
549 for (label iCPv = 0; iCPv < nCPsV; iCPv++)
550 {
551 for (label iCPu = 0; iCPu < nCPsU; iCPu++)
552 {
553 confineControlPoint(getCPID(iCPu, iCPv, iCPw), confineSlice);
554 }
555 }
556 }
557
558 forAll(confineWMaxCPs_, sliceI)
559 {
560 const boolVector& confineSlice = confineWMaxCPs_[sliceI];
561 label iCPw = nCPsW - 1 - sliceI;
562 // Control points at the end of the parameterized space
563 for (label iCPv = 0; iCPv < nCPsV; iCPv++)
564 {
565 for (label iCPu = 0; iCPu < nCPsU; iCPu++)
566 {
567 confineControlPoint(getCPID(iCPu, iCPv, iCPw), confineSlice);
568 }
569 }
570 }
571}
572
573
575{
576 for (label cpI = 0; cpI < cps_.size(); ++cpI)
577 {
578 if (confineUMovement_) activeDesignVariables_[3*cpI] = false;
579 if (confineVMovement_) activeDesignVariables_[3*cpI + 1] = false;
580 if (confineWMovement_) activeDesignVariables_[3*cpI + 2] = false;
581 }
582}
583
584
586{
587 if (cpI < 0 || cpI > cps_.size() -1)
588 {
590 << "Attempted to confine control point movement for a control point "
591 << " ID which is out of bounds"
592 << exit(FatalError);
593 }
594 else
595 {
596 activeDesignVariables_[3*cpI] = false;
597 activeDesignVariables_[3*cpI + 1] = false;
598 activeDesignVariables_[3*cpI + 2] = false;
599 }
600}
601
602
604(
605 const label cpI,
606 const boolVector& confineDirections
607)
608{
609 if (cpI < 0 || cpI > cps_.size() -1)
610 {
612 << "Attempted to confine control point movement for a control point "
613 << " ID which is out of bounds"
614 << exit(FatalError);
615 }
616 else
617 {
618 if (confineDirections.x()) activeDesignVariables_[3*cpI] = false;
619 if (confineDirections.y()) activeDesignVariables_[3*cpI + 1] = false;
620 if (confineDirections.z()) activeDesignVariables_[3*cpI + 2] = false;
621 }
622}
623
624
625// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
626
628(
629 const dictionary& dict,
630 const fvMesh& mesh,
631 bool computeParamCoors
632)
633:
635 (
637 (
638 dict.dictName() + "cpsBsplines",
639 mesh.time().timeName(),
640 fileName("uniform")/fileName("volumetricBSplines"),
641 mesh,
642 IOobject::READ_IF_PRESENT,
643 IOobject::AUTO_WRITE
644 ),
645 word::null
646 ),
647 mesh_(mesh),
648 dict_(dict),
649 name_(dict.dictName()),
650 basisU_(dict.get<label>("nCPsU"), dict.get<label>("degreeU")),
651 basisV_(dict.get<label>("nCPsV"), dict.get<label>("degreeV")),
652 basisW_(dict.get<label>("nCPsW"), dict.get<label>("degreeW")),
653 maxIter_(dict.getOrDefault<label>("maxIterations", 10)),
654 tolerance_(dict.getOrDefault<scalar>("tolerance", 1.e-10)),
655 nMaxBound_(dict.getOrDefault<scalar>("nMaxBoundIterations", 4)),
656 cps_(),
657 mapPtr_(nullptr),
658 reverseMapPtr_(nullptr),
659 parametricCoordinatesPtr_(nullptr),
660 localSystemCoordinates_(mesh_.nPoints(), Zero),
661 confineUMovement_
662 (
663 dict.getOrDefaultCompat<bool>
664 (
665 "confineUMovement", {{"confineX1movement", 1912}}, false
666 )
667 ),
668 confineVMovement_
669 (
671 (
672 "confineVMovement", {{"confineX2movement", 1912}}, false
673 )
674 ),
675 confineWMovement_
676 (
678 (
679 "confineWMovement", {{"confineX3movement", 1912}}, false
680 )
681 ),
682 confineBoundaryControlPoints_
683 (
684 dict.getOrDefault<bool>("confineBoundaryControlPoints", true)
685 ),
686 confineUMinCPs_
687 (
688 dict.getOrDefaultCompat<boolVectorList>
689 (
690 "confineUMinCPs", {{"boundUMinCPs", 1912}}, boolVectorList()
691 )
692 ),
693 confineUMaxCPs_
694 (
695 dict.getOrDefaultCompat<boolVectorList>
696 (
697 "confineUMaxCPs", {{"boundUMaxCPs", 1912}}, boolVectorList()
698 )
699 ),
700 confineVMinCPs_
701 (
702 dict.getOrDefaultCompat<boolVectorList>
703 (
704 "confineVMinCPs", {{"boundVMinCPs", 1912}}, boolVectorList()
705 )
706 ),
707 confineVMaxCPs_
708 (
709 dict.getOrDefaultCompat<boolVectorList>
710 (
711 "confineVMaxCPs", {{"boundVMaxCPs", 1912}}, boolVectorList()
712 )
713 ),
714 confineWMinCPs_
715 (
716 dict.getOrDefaultCompat<boolVectorList>
717 (
718 "confineWMinCPs", {{"boundWMinCPs", 1912}}, boolVectorList()
719 )
720 ),
721 confineWMaxCPs_
722 (
723 dict.getOrDefaultCompat<boolVectorList>
724 (
725 "confineWMaxCPs", {{"boundWMaxCPs", 1912}}, boolVectorList()
726 )
727 ),
728 activeControlPoints_(0), //zero here, execute sanity checks first
729 activeDesignVariables_(0), //zero here, execute sanity checks first
730 cpsFolder_("controlPoints"),
731 readStoredData_(dict.getOrDefault<bool>("readStoredData", true))
732{
733 // Create folders
734 makeFolders();
735
736 // Sanity checks
737 if
738 (
739 (confineUMinCPs_.size() + confineUMaxCPs_.size() >= basisU_.nCPs())
740 || (confineVMinCPs_.size() + confineVMaxCPs_.size() >= basisV_.nCPs())
741 || (confineWMinCPs_.size() + confineWMaxCPs_.size() >= basisW_.nCPs())
742 )
743 {
745 << "Number of control point slices to be kept frozen at "
746 << "the boundaries is invalid \n"
747 << "Number of control points in u " << basisU_.nCPs() << "\n"
748 << "Number of control points in v " << basisV_.nCPs() << "\n"
749 << "Number of control points in w " << basisW_.nCPs() << "\n"
750 << exit(FatalError);
751 }
752
753 // Construct control points, if not already read from file
754 if (found("controlPoints"))
755 {
756 cps_ =
758 (
759 "controlPoints",
760 *this,
761 basisU_.nCPs()*basisV_.nCPs()*basisW_.nCPs()
762 );
763 }
764 else
765 {
767 }
768 determineActiveDesignVariablesAndPoints();
769}
770
771
772// * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
773
775(
776 const dictionary& dict,
777 const fvMesh& mesh,
778 bool computeParamCoors
779)
780{
781 const word modelType(dict.get<word>("type"));
782
783 Info<< "NURBS3DVolume type : " << modelType << endl;
784
785 auto* ctorPtr = dictionaryConstructorTable(modelType);
786
787 if (!ctorPtr)
788 {
790 (
791 dict,
792 "type",
793 modelType,
794 *dictionaryConstructorTablePtr_
795 ) << exit(FatalIOError);
796 }
797
798 return autoPtr<NURBS3DVolume>(ctorPtr(dict, mesh, computeParamCoors));
799}
800
801
802// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
803
805(
806 const scalar u,
807 const scalar v,
808 const scalar w
809) const
810{
811 const label degreeU = basisU_.degree();
812 const label degreeV = basisV_.degree();
813 const label degreeW = basisW_.degree();
814
815 const label nCPsU = basisU_.nCPs();
816 const label nCPsV = basisV_.nCPs();
817 const label nCPsW = basisW_.nCPs();
818
819 vector derivative(Zero);
820
821 for (label iCPw = 0; iCPw < nCPsW; ++iCPw)
822 {
823 const scalar basisW(basisW_.basisValue(iCPw, degreeW, w));
824 for (label iCPv = 0; iCPv < nCPsV; ++iCPv)
825 {
826 const scalar basisVW = basisW*basisV_.basisValue(iCPv, degreeV, v);
827 for (label iCPu = 0; iCPu < nCPsU; ++iCPu)
828 {
829 derivative +=
830 cps_[getCPID(iCPu, iCPv, iCPw)]
831 *basisU_.basisDerivativeU(iCPu, degreeU, u)
832 *basisVW;
833 }
834 }
835 }
836
837 return derivative;
838}
839
840
842(
843 const scalar u,
844 const scalar v,
845 const scalar w
846) const
847{
848 const label degreeU = basisU_.degree();
849 const label degreeV = basisV_.degree();
850 const label degreeW = basisW_.degree();
851
852 const label nCPsU = basisU_.nCPs();
853 const label nCPsV = basisV_.nCPs();
854 const label nCPsW = basisW_.nCPs();
855
856 vector derivative(Zero);
857
858 for (label iCPw = 0; iCPw < nCPsW; ++iCPw)
859 {
860 const scalar basisW(basisW_.basisValue(iCPw, degreeW, w));
861 for (label iCPv = 0; iCPv < nCPsV; ++iCPv)
862 {
863 const scalar basisWDeriV =
864 basisW*basisV_.basisDerivativeU(iCPv, degreeV, v);
865 for (label iCPu = 0; iCPu < nCPsU; ++iCPu)
866 {
867 derivative +=
868 cps_[getCPID(iCPu, iCPv, iCPw)]
869 *basisU_.basisValue(iCPu, degreeU, u)
870 *basisWDeriV;
871 }
872 }
873 }
874
875 return derivative;
876}
877
878
880(
881 const scalar u,
882 const scalar v,
883 const scalar w
884) const
885{
886 const label degreeU = basisU_.degree();
887 const label degreeV = basisV_.degree();
888 const label degreeW = basisW_.degree();
889
890 const label nCPsU = basisU_.nCPs();
891 const label nCPsV = basisV_.nCPs();
892 const label nCPsW = basisW_.nCPs();
893
894 vector derivative(Zero);
895
896 for (label iCPw = 0; iCPw < nCPsW; iCPw++)
897 {
898 const scalar derivW(basisW_.basisDerivativeU(iCPw, degreeW, w));
899 for (label iCPv = 0; iCPv < nCPsV; iCPv++)
900 {
901 const scalar derivWBasisV =
902 derivW*basisV_.basisValue(iCPv, degreeV, v);
903 for (label iCPu = 0; iCPu < nCPsU; iCPu++)
904 {
905 derivative +=
906 cps_[getCPID(iCPu, iCPv, iCPw)]
907 *basisU_.basisValue(iCPu, degreeU, u)
908 *derivWBasisV;
909 }
910 }
911 }
912
913 return derivative;
914}
915
916
918(
919 const vector& uVector
920) const
921{
922 const scalar u = uVector.x();
923 const scalar v = uVector.y();
924 const scalar w = uVector.z();
925
926 vector uDeriv = volumeDerivativeU(u, v, w);
927 vector vDeriv = volumeDerivativeV(u, v, w);
928 vector wDeriv = volumeDerivativeW(u, v, w);
929
930 tensor Jacobian(uDeriv, vDeriv, wDeriv, true);
931
932 return Jacobian;
933}
934
935
937(
938 const vector& uVector,
939 const label cpI
940) const
941{
942 const scalar u = uVector.x();
943 const scalar v = uVector.y();
944 const scalar w = uVector.z();
945
946 const label nCPsU = basisU_.nCPs();
947 const label nCPsV = basisV_.nCPs();
948
949 const label degreeU = basisU_.degree();
950 const label degreeV = basisV_.degree();
951 const label degreeW = basisW_.degree();
952
953 label iCPw = cpI/label(nCPsU*nCPsV);
954 label iCPv = (cpI - iCPw*nCPsU*nCPsV)/nCPsU;
955 label iCPu = (cpI - iCPw*nCPsU*nCPsV - iCPv*nCPsU);
956
957 // Normally, this should be a tensor, however the parameterization is
958 // isotropic. Hence the tensor degenerates to a diagonal tensor with all
959 // diagonal elements being equal. This returns the (unique) diag element
960 scalar derivative =
961 basisU_.basisValue(iCPu, degreeU, u)
962 *basisV_.basisValue(iCPv, degreeV, v)
963 *basisW_.basisValue(iCPw, degreeW, w);
964
965 return derivative;
966}
967
968
970(
971 const pointVectorField& pointSens,
972 const labelList& sensitivityPatchIDs
973)
974{
975 vectorField controlPointDerivs(cps_.size(), Zero);
976
977 // Get parametric coordinates
978 const vectorField& parametricCoordinates = getParametricCoordinates();
979
980 forAll(controlPointDerivs, cpI)
981 {
982 forAll(sensitivityPatchIDs, pI)
983 {
984 const label patchI = sensitivityPatchIDs[pI];
985 const polyPatch& patch = mesh_.boundaryMesh()[patchI];
986 const labelList& meshPoints = patch.meshPoints();
987
988 forAll(meshPoints, mpI)
989 {
990 const label globalIndex = meshPoints[mpI];
991 const label whichPointInBox = reverseMapPtr_()[globalIndex];
992
993 // If point resides within control points box,
994 // add contribution to cp derivative
995 if (whichPointInBox != -1)
996 {
997 controlPointDerivs[cpI] +=
998 (
999 pointSens[globalIndex]
1000 & transformationTensorDxDb(globalIndex)
1001 )
1002 *volumeDerivativeCP
1003 (
1004 parametricCoordinates[globalIndex],
1005 cpI
1006 );
1007 }
1008 }
1009 }
1010 }
1011
1012 // Sum contributions from all processors
1013 Pstream::listCombineAllGather(controlPointDerivs, plusEqOp<vector>());
1014
1015 return controlPointDerivs;
1016}
1017
1018
1020(
1021 const volVectorField& faceSens,
1022 const labelList& sensitivityPatchIDs
1023)
1024{
1025 return
1026 computeControlPointSensitivities
1027 (
1028 faceSens.boundaryField(),
1029 sensitivityPatchIDs
1030 );
1031}
1032
1033
1035(
1036 const boundaryVectorField& faceSens,
1037 const labelList& sensitivityPatchIDs
1038)
1039{
1040 // Return field
1041 vectorField controlPointDerivs(cps_.size(), Zero);
1042
1043 // Get parametric coordinates
1044 const vectorField& parametricCoordinates = getParametricCoordinates();
1045
1046 // Auxiliary quantities
1048 const labelList& reverseMap = reverseMapPtr_();
1049
1050 forAll(controlPointDerivs, cpI)
1051 {
1052 forAll(sensitivityPatchIDs, pI)
1053 {
1054 const label patchI = sensitivityPatchIDs[pI];
1055 const polyPatch& patch = mesh_.boundaryMesh()[patchI];
1056 const label patchStart = patch.start();
1057 const fvPatchVectorField& patchSens = faceSens[patchI];
1058
1059 // loop over patch faces
1060 forAll(patch, fI)
1061 {
1062 const face& fGlobal = mesh_.faces()[fI + patchStart];
1063 const pointField facePoints = fGlobal.points(mesh_.points());
1064 // loop over face points
1065 tensorField facePointDerivs(facePoints.size(), Zero);
1066 forAll(fGlobal, pI)
1067 {
1068 const label globalIndex = fGlobal[pI]; //global point index
1069 const label whichPointInBox = reverseMap[globalIndex];
1070 // if point resides within control points box,
1071 // add contribution to d( facePoints )/db
1072 if (whichPointInBox != -1)
1073 {
1074 // TENSOR-BASED
1075 //~~~~~~~~~~~~~
1076 facePointDerivs[pI] =
1077 transformationTensorDxDb(globalIndex)
1078 * volumeDerivativeCP
1079 (
1080 parametricCoordinates[globalIndex],
1081 cpI
1082 );
1083
1084 }
1085 }
1086
1087 tensor fCtrs_d =
1089 (
1090 facePoints,
1091 facePointDerivs
1092 )[0];
1093 controlPointDerivs[cpI] += patchSens[fI] & fCtrs_d;
1094 }
1095 }
1096 }
1097 // Sum contributions from all processors
1098 Pstream::listCombineAllGather(controlPointDerivs, plusEqOp<vector>());
1099
1100 return controlPointDerivs;
1101}
1102
1103
1105(
1106 const vectorField& faceSens,
1107 const label patchI,
1108 const label cpI
1109)
1110{
1111 // Return vector
1112 vector cpSens(Zero);
1113 // Get parametric coordinates
1114 const vectorField& parametricCoordinates = getParametricCoordinates();
1115
1116 // Auxiliary quantities
1118 const labelList& reverseMap = reverseMapPtr_();
1119
1120 const polyPatch& patch = mesh_.boundaryMesh()[patchI];
1121 const label patchStart = patch.start();
1122 // Loop over patch faces
1123 forAll(patch, fI)
1124 {
1125 const face& fGlobal = mesh_.faces()[fI + patchStart];
1126 const pointField facePoints = fGlobal.points(mesh_.points());
1127 // Loop over face points
1128 tensorField facePointDerivs(facePoints.size(), Zero);
1129 forAll(fGlobal, pI)
1130 {
1131 const label globalIndex = fGlobal[pI]; //global point index
1132 const label whichPointInBox = reverseMap[globalIndex];
1133 // If point resides within control points box,
1134 // add contribution to d( facePoints )/db
1135 if (whichPointInBox != -1)
1136 {
1137 // TENSOR-BASED
1138 //~~~~~~~~~~~~~
1139 facePointDerivs[pI] =
1140 transformationTensorDxDb(globalIndex)
1141 *volumeDerivativeCP
1142 (
1143 parametricCoordinates[globalIndex],
1144 cpI
1145 );
1146 }
1147 }
1148
1149 tensor fCtrs_d =
1151 (
1152 facePoints,
1153 facePointDerivs
1154 )[0];
1155 cpSens += faceSens[fI] & fCtrs_d;
1156 }
1157 // Sum contributions from all processors
1158 reduce(cpSens, sumOp<vector>());
1159
1160 return cpSens;
1161}
1162
1163
1165(
1166 const label patchI,
1167 const label cpI,
1168 bool DimensionedNormalSens
1169)
1170{
1171 const fvPatch& patch = mesh_.boundary()[patchI];
1172 const polyPatch& ppatch = patch.patch();
1173 // Return field
1174 tmp<tensorField> tdndbSens(new tensorField(patch.size(), Zero));
1175 tensorField& dndbSens = tdndbSens.ref();
1176 // Auxiliary quantities
1178 const label patchStart = ppatch.start();
1179 const labelList& reverseMap = reverseMapPtr_();
1180
1181 // Get parametric coordinates
1182 const vectorField& parametricCoordinates = getParametricCoordinates();
1183
1184 // Loop over patch faces
1185 forAll(patch, fI)
1186 {
1187 const face& fGlobal = mesh_.faces()[fI + patchStart];
1188 const pointField facePoints = fGlobal.points(mesh_.points());
1189 // Loop over face points
1190 tensorField facePointDerivs(facePoints.size(), Zero);
1191 forAll(fGlobal, pI)
1192 {
1193 const label globalIndex = fGlobal[pI]; //global point index
1194 const label whichPointInBox = reverseMap[globalIndex];
1195 // If point resides within control points box,
1196 // add contribution to d( facePoints )/db
1197 if (whichPointInBox != -1)
1198 {
1199 // TENSOR-BASED
1200 //~~~~~~~~~~~~~
1201 facePointDerivs[pI] =
1202 transformationTensorDxDb(globalIndex)
1203 *volumeDerivativeCP
1204 (
1205 parametricCoordinates[globalIndex],
1206 cpI
1207 );
1208 }
1209 }
1210
1211 // Determine whether to return variance of dimensioned or unit normal
1212 tensorField dNdbSens =
1214 (
1215 facePoints,
1216 facePointDerivs
1217 );
1218
1219 if (DimensionedNormalSens)
1220 {
1221 dndbSens[fI] = dNdbSens[1];
1222 }
1223 else
1224 {
1225 dndbSens[fI] = dNdbSens[2];
1226 }
1227 }
1228
1229 return tdndbSens;
1230}
1231
1232
1234(
1235 const label patchI,
1236 const label cpI
1237)
1238{
1239 // Get parametric coordinates
1240 const vectorField& parametricCoordinates = getParametricCoordinates();
1241
1242 // Patch data
1243 const polyPatch& patch = mesh_.boundaryMesh()[patchI];
1244 const labelList& meshPoints = patch.meshPoints();
1245
1246 // Return field
1247 auto tdxdb = tmp<tensorField>::New(patch.nPoints(), Zero);
1248 auto& dxdb = tdxdb.ref();
1249
1250 forAll(meshPoints, pI)
1251 {
1252 const label globalIndex = meshPoints[pI]; //global point index
1253 const label whichPointInBox = reverseMapPtr_()[globalIndex];
1254
1255 // If point resides within control points box, find dxdb
1256 if (whichPointInBox != -1)
1257 {
1258 dxdb[pI] =
1259 transformationTensorDxDb(globalIndex)
1260 *volumeDerivativeCP
1261 (
1262 parametricCoordinates[globalIndex],
1263 cpI
1264 );
1265 }
1266 }
1267
1268 return tdxdb;
1269}
1270
1271
1273(
1274 const label patchI,
1275 const label cpI
1276)
1277{
1278 // get parametric coordinates
1279 const vectorField& parametricCoordinates = getParametricCoordinates();
1280
1281 // Patch data
1282 const polyPatch& patch = mesh_.boundaryMesh()[patchI];
1283 const label patchStart = patch.start();
1284
1285 // Return field
1286 auto tdxdb = tmp<tensorField>::New(patch.size(), Zero);
1287 auto& dxdb = tdxdb.ref();
1288
1289 // Mesh differentiation engine
1290 deltaBoundary deltaBound(mesh_);
1291
1292 forAll(patch, fI)
1293 {
1294 const face& fGlobal = mesh_.faces()[fI + patchStart];
1295 const pointField facePoints = fGlobal.points(mesh_.points());
1296 // Loop over face points
1297 tensorField facePointDerivs(facePoints.size(), Zero);
1298 forAll(fGlobal, pI)
1299 {
1300 const label globalIndex = fGlobal[pI]; //global point index
1301 const label whichPointInBox = reverseMapPtr_()[globalIndex];
1302 // If point resides within control points box,
1303 // add contribution to d( facePoints )/db
1304 if (whichPointInBox != -1)
1305 {
1306 // TENSOR-BASED
1307 //~~~~~~~~~~~~~
1308 facePointDerivs[pI] =
1309 transformationTensorDxDb(globalIndex)
1310 *volumeDerivativeCP
1311 (
1312 parametricCoordinates[globalIndex],
1313 cpI
1314 );
1315
1316 }
1317 }
1318 dxdb[fI] =
1319 deltaBound.makeFaceCentresAndAreas_d
1320 (
1321 facePoints,
1322 facePointDerivs
1323 )[0];
1324 }
1325
1326 return tdxdb;
1327}
1328
1329
1331(
1332 const vectorField& uVector
1333) const
1334{
1335 const label nPoints = mapPtr_().size();
1336 auto tpoints = tmp<vectorField>::New(nPoints, Zero);
1337 auto& points = tpoints.ref();
1338
1339 forAll(points, pI)
1340 {
1341 const label globalPI = mapPtr_()[pI];
1342 points[pI] = coordinates(uVector[globalPI]);
1343 }
1344
1345 return tpoints;
1346}
1347
1348
1350(
1351 const vector& uVector
1352) const
1353{
1354 const label degreeU = basisU_.degree();
1355 const label degreeV = basisV_.degree();
1356 const label degreeW = basisW_.degree();
1357
1358 const label nCPsU = basisU_.nCPs();
1359 const label nCPsV = basisV_.nCPs();
1360 const label nCPsW = basisW_.nCPs();
1361
1362 const scalar u = uVector.x();
1363 const scalar v = uVector.y();
1364 const scalar w = uVector.z();
1365
1366 vector point(Zero);
1367 for (label iCPw = 0; iCPw < nCPsW; iCPw++)
1368 {
1369 const scalar basisW(basisW_.basisValue(iCPw, degreeW, w));
1370 for (label iCPv = 0; iCPv < nCPsV; iCPv++)
1371 {
1372 const scalar basisVW =
1373 basisW*basisV_.basisValue(iCPv, degreeV, v);
1374 for (label iCPu = 0; iCPu < nCPsU; iCPu++)
1375 {
1376 point +=
1377 cps_[getCPID(iCPu, iCPv, iCPw)]
1378 *basisU_.basisValue(iCPu, degreeU, u)
1379 *basisVW;
1380 }
1381 }
1382 }
1383
1384 return point;
1385}
1386
1387
1389(
1390 const vectorField& controlPointsMovement
1391)
1392{
1393 // Get parametric coordinates and map
1394 const vectorField& paramCoors = getParametricCoordinates();
1395 const labelList& map = mapPtr_();
1396
1397 // Update control points position
1398 cps_ += controlPointsMovement;
1399 writeCps("cpsBsplines"+mesh_.time().timeName());
1400
1401 // Compute new mesh points based on updated control points
1402 tmp<vectorField> tparameterizedPoints = coordinates(paramCoors);
1403 const vectorField& parameterizedPoints = tparameterizedPoints();
1404
1405 // Return field. Initialized with current mesh points
1406 tmp<vectorField> tnewPoints(new vectorField(mesh_.points()));
1407 vectorField& newPoints = tnewPoints.ref();
1408
1409 // Update position of parameterized points
1410 forAll(parameterizedPoints, pI)
1411 {
1412 newPoints[map[pI]] = transformPointToCartesian(parameterizedPoints[pI]);
1413 }
1414
1415 // Update coordinates in the local system based on the cartesian points
1416 updateLocalCoordinateSystem(newPoints);
1417 DebugInfo
1418 << "Max mesh movement equal to "
1419 << gMax(mag(newPoints - mesh_.points())) << endl;
1420
1421 return tnewPoints;
1422}
1423
1424
1426(
1427 const vectorField& controlPointsMovement,
1428 const labelList& patchesToBeMoved,
1429 const bool updateCPs
1430)
1431{
1432 // Get parametric coordinates
1433 const vectorField& paramCoors = getParametricCoordinates();
1434
1435 // Update control points position
1436 cps_ += controlPointsMovement;
1437
1438 if (updateCPs)
1439 {
1440 writeCps("cpsBsplines"+mesh_.time().timeName());
1441 }
1442
1443 // Return field. Initialized with current mesh points
1444 tmp<vectorField> tnewPoints(new vectorField(mesh_.points()));
1445 vectorField& newPoints = tnewPoints.ref();
1446
1447 // Update position of parameterized boundary points
1448 for (const label patchI : patchesToBeMoved)
1449 {
1450 const polyPatch& patch = mesh_.boundaryMesh()[patchI];
1451 const labelList& meshPoints = patch.meshPoints();
1452
1453 for (const label globalIndex : meshPoints)
1454 {
1455 const label whichPointInBox = reverseMapPtr_()[globalIndex];
1456 // If point resides within control points box,
1457 // compute new cartesian coordinates
1458 if (whichPointInBox != -1)
1459 {
1460 newPoints[globalIndex] =
1461 transformPointToCartesian
1462 (
1464 (
1465 paramCoors[globalIndex]
1466 )
1467 );
1468 }
1469 }
1470 }
1471
1472 if (updateCPs)
1473 {
1474 // Update coordinates in the local system based on the cartesian points
1475 updateLocalCoordinateSystem(newPoints);
1476 }
1477 else
1478 {
1479 // Move control points to their initial position
1480 cps_ -= controlPointsMovement;
1481 }
1482
1483 DebugInfo
1484 << "Max mesh movement equal to "
1485 << gMax(mag(newPoints - mesh_.points())) << endl;
1486
1487 return tnewPoints;
1488}
1489
1490
1492(
1493 const label i,
1494 const label j,
1495 const label k
1496) const
1497{
1498 const label nCPsU = basisU_.nCPs();
1499 const label nCPsV = basisV_.nCPs();
1500
1501 return k*nCPsU*nCPsV + j*nCPsU + i;
1502}
1503
1504
1506{
1507 if (cps_.size() != newCps.size())
1508 {
1510 << "Attempting to replace control points with a set of "
1511 << "different size"
1512 << exit(FatalError);
1513 }
1514 cps_ = newCps;
1515}
1516
1517
1519(
1520 vectorField& controlPointsMovement
1521) const
1522{
1523 forAll(controlPointsMovement, cpI)
1524 {
1525 if (!activeDesignVariables_[3*cpI])
1526 {
1527 controlPointsMovement[cpI].x() = Zero;
1528 }
1529 if (!activeDesignVariables_[3*cpI + 1])
1530 {
1531 controlPointsMovement[cpI].y() = Zero;
1532 }
1533 if (!activeDesignVariables_[3*cpI + 2])
1534 {
1535 controlPointsMovement[cpI].z() = Zero;
1536 }
1537 }
1538}
1539
1540
1542(
1543 const vectorField& controlPointsMovement,
1544 const labelList& patchesToBeMoved
1545)
1546{
1547 // Backup old cps
1548 vectorField oldCPs = cps_;
1549 // Get parametric coordinates
1550 const vectorField& paramCoors = getParametricCoordinates();
1551 // Update control points position
1552 cps_ += controlPointsMovement;
1553 // Update position of parameterized boundary points
1554 scalar maxDisplacement(Zero);
1555 for (const label patchI : patchesToBeMoved)
1556 {
1557 const polyPatch& patch = mesh_.boundaryMesh()[patchI];
1558 const labelList& meshPoints = patch.meshPoints();
1559
1560 for (const label globalIndex : meshPoints)
1561 {
1562 const label whichPointInBox = reverseMapPtr_()[globalIndex];
1563 // If point resides within control points box,
1564 // compute new cartesian coordinates
1565 if (whichPointInBox != -1)
1566 {
1567 vector newPoint =
1568 transformPointToCartesian
1569 (
1571 (
1572 paramCoors[globalIndex]
1573 )
1574 );
1575 maxDisplacement =
1576 max
1577 (
1578 maxDisplacement,
1579 mag(newPoint - mesh_.points()[globalIndex])
1580 );
1581 }
1582 }
1583 }
1584 reduce(maxDisplacement, maxOp<scalar>());
1585 cps_ = oldCPs;
1586
1587 return maxDisplacement;
1588}
1589
1590
1592{
1593 if (!mapPtr_)
1594 {
1595 findPointsInBox(localSystemCoordinates_);
1596 }
1597 tmp<vectorField> pointsInBox
1598 (
1599 new vectorField(localSystemCoordinates_, mapPtr_())
1600 );
1601
1602 return pointsInBox;
1603}
1604
1605
1607{
1608 if (!mapPtr_)
1609 {
1610 findPointsInBox(localSystemCoordinates_);
1611 }
1612
1613 return mapPtr_();
1614}
1615
1616
1618{
1619 if (!reverseMapPtr_)
1620 {
1621 findPointsInBox(localSystemCoordinates_);
1622 }
1623
1624 return reverseMapPtr_();
1625}
1626
1627
1629{
1630 // If not computed yet, compute parametric coordinates
1631 if (!parametricCoordinatesPtr_)
1632 {
1633 // Find mesh points inside control points box
1634 // if they have been identified yet
1635 if (!mapPtr_)
1636 {
1637 findPointsInBox(localSystemCoordinates_);
1638 }
1639 computeParametricCoordinates(getPointsInBox()());
1640 }
1641
1642 return parametricCoordinatesPtr_();
1643}
1644
1645
1647{
1648 // Get parametric coordinates
1649 const vectorField& parametricCoordinates = getParametricCoordinates();
1650
1651 // Set return field to zero
1653 (
1655 (
1656 IOobject
1657 (
1658 "DxDb",
1659 mesh_.time().timeName(),
1660 mesh_,
1663 ),
1664 pointMesh::New(mesh_),
1666 )
1667 );
1668
1669 pointTensorField& DxDb = tDxDb.ref();
1670
1671 // All points outside the control box remain unmoved.
1672 // Loop over only points within the control box
1673 const labelList& map = mapPtr_();
1674 for (const label globalIndex : map)
1675 {
1676 DxDb[globalIndex] =
1677 transformationTensorDxDb(globalIndex)
1678 *volumeDerivativeCP
1679 (
1680 parametricCoordinates[globalIndex],
1681 cpI
1682 );
1683 }
1684
1685 return tDxDb;
1686}
1687
1688
1690(
1691 const label cpI
1692)
1693{
1694 // Get parametric coordinates
1695 const vectorField& parametricCoordinates = getParametricCoordinates();
1696
1697 // Set return field to zero
1699 (
1700 new volTensorField
1701 (
1702 IOobject
1703 (
1704 "DxDb",
1705 mesh_.time().timeName(),
1706 mesh_,
1709 ),
1710 mesh_,
1712 )
1713 );
1714
1715 volTensorField& DxDb = tDxDb.ref();
1716 deltaBoundary deltaBound(mesh_);
1717 const labelListList& pointCells = mesh_.pointCells();
1718
1719 // All points outside the control box remain unmoved.
1720 // Loop over only points within the control box
1721 const labelList& map = mapPtr_();
1722 for (const label globalIndex : map)
1723 {
1724 tensor pointDxDb =
1725 transformationTensorDxDb(globalIndex)
1726 *volumeDerivativeCP
1727 (
1728 parametricCoordinates[globalIndex],
1729 cpI
1730 );
1731 const labelList& pointCellsI = pointCells[globalIndex];
1732 tmp<tensorField> tC_d = deltaBound.cellCenters_d(globalIndex);
1733 const tensorField& C_d = tC_d();
1734 forAll(pointCellsI, cI)
1735 {
1736 const label cellI = pointCellsI[cI];
1737 DxDb[cellI] += C_d[cI] & pointDxDb;
1738 }
1739 }
1740
1741 // Assign boundary values since the grad of this field is often needed
1742 forAll(mesh_.boundary(), pI)
1743 {
1744 const fvPatch& patch = mesh_.boundary()[pI];
1745 if (!isA<coupledFvPatch>(patch))
1746 {
1747 DxDb.boundaryFieldRef()[pI] = patchDxDbFace(pI, cpI);
1748 }
1749 }
1750
1751 // Correct coupled boundaries
1753
1754 return tDxDb;
1755}
1756
1757
1759{
1760 label nU(basisU_.nCPs());
1761 return label(nU % 2 == 0 ? 0.5*nU : 0.5*(nU - 1) + 1);
1762}
1763
1764
1766{
1767 label nV(basisV_.nCPs());
1768 return label(nV % 2 == 0 ? 0.5*nV : 0.5*(nV - 1) + 1);
1769}
1770
1771
1773{
1774 label nW(basisW_.nCPs());
1775 return label(nW % 2 == 0 ? 0.5*nW : 0.5*(nW - 1) + 1);
1776}
1777
1778
1780{
1781 return Vector<label>(nUSymmetry(), nVSymmetry(), nWSymmetry());
1782}
1783
1784
1786(
1787 const fileName& baseName,
1788 const bool transform
1789) const
1790{
1791 const label nCPsU = basisU_.nCPs();
1792 const label nCPsV = basisV_.nCPs();
1793
1794 vectorField cpsInCartesian(cps_);
1795 if (transform)
1796 {
1797 forAll(cpsInCartesian, cpI)
1798 {
1799 cpsInCartesian[cpI] = transformPointToCartesian(cps_[cpI]);
1800 }
1801 }
1802
1803 Info<< "Writing control point positions to file" << endl;
1804
1805 if (Pstream::master())
1806 {
1807 OFstream cpsFile("optimisation"/cpsFolder_/name_ + baseName + ".csv");
1808 // Write header
1809 cpsFile
1810 << "\"Points : 0\", \"Points : 1\", \"Points : 2\","
1811 << "\"i\", \"j\", \"k\","
1812 << "\"active : 0\", \"active : 1\", \"active : 2\"" << endl;
1813
1814 forAll(cpsInCartesian, cpI)
1815 {
1816 const label iCPw = cpI/label(nCPsU*nCPsV);
1817 const label iCPv = (cpI - iCPw*nCPsU*nCPsV)/nCPsU;
1818 const label iCPu = (cpI - iCPw*nCPsU*nCPsV - iCPv*nCPsU);
1819
1820 cpsFile
1821 << cpsInCartesian[cpI].x() << ", "
1822 << cpsInCartesian[cpI].y() << ", "
1823 << cpsInCartesian[cpI].z() << ", "
1824 << iCPu << ", "
1825 << iCPv << ", "
1826 << iCPw << ", "
1827 << activeDesignVariables_[3*cpI] << ", "
1828 << activeDesignVariables_[3*cpI + 1] << ", "
1829 << activeDesignVariables_[3*cpI + 2] << endl;
1830 }
1831 }
1832}
1833
1834
1836{
1837 parametricCoordinatesPtr_().write();
1838}
1839
1840
1842{
1843 cps_.writeEntry("controlPoints", os);
1844 return true;
1845}
1846
1847
1848// ************************************************************************* //
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:545
Generic GeometricBoundaryField class.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
void correctBoundaryConditions()
Correct boundary field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
void setSize(const label n)
Alias for resize()
Definition: List.H:218
NURBS3DVolume morpher. Includes support functions for gradient computations Base class providing supp...
Definition: NURBS3DVolume.H:76
tmp< tensorField > patchDxDbFace(const label patchI, const label cpI)
Get patch dx/db.
void computeParametricCoordinates(const vectorField &points)
label nUSymmetry() const
Get number of variables if CPs are moved symmetrically in U.
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.
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.
void findPointsInBox(const vectorField &meshPoints)
Find points within control points box.
Definition: NURBS3DVolume.C:50
tmp< vectorField > computeNewBoundaryPoints(const vectorField &controlPointsMovement, const labelList &patchesToBeMoved, const bool moveCPs=true)
Boundary mesh movement based on given control point movement.
vectorField computeControlPointSensitivities(const pointVectorField &pointSens, const labelList &sensitivityPatchIDs)
const pointVectorField & getParametricCoordinates()
Get parametric coordinates.
tmp< volTensorField > getDxCellsDb(const label cpI)
Get dxCartesiandb for a certain control point on cells.
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.
void makeFolders()
Create folders to store cps and derivatives.
void writeCps(const fileName &="cpsFile", const bool transform=true) const
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.
tmp< tensorField > dndbBasedSensitivities(const label patchI, const label cpI, bool DimensionedNormalSens=true)
void confineControlPoint(const label cpI)
Confine all three movements for a prescribed control point.
Vector< label > nSymmetry() const
vectorField cps_
The volumetric B-Splines control points.
void boundControlPointMovement(vectorField &controlPointsMovement) const
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.
scalar volumeDerivativeCP(const vector &u, const label cpI) const
Volume point derivative wrt to control point cpI at point u,v,w.
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 > 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.
const labelList & getReverseMap()
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, using an OSstream.
Definition: OFstream.H:57
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
static void listCombineAllGather(const List< commsStruct > &comms, List< T > &values, const CombineOp &cop, const int tag, const label comm)
After completion all processors have the same data.
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
const Cmpt & component(const direction) const
Definition: VectorSpaceI.H:87
const Cmpt & z() const
Access to the vector z component.
Definition: VectorI.H:85
const Cmpt & y() const
Access to the vector y component.
Definition: VectorI.H:79
const Cmpt & x() const
Access to the vector x component.
Definition: VectorI.H:73
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
Specialized bundling of boolean values as a vector of 3 components, element access using x(),...
Definition: boolVector.H:59
bool x() const
The x component.
Definition: boolVectorI.H:110
bool y() const
The y component.
Definition: boolVectorI.H:111
bool z() const
The z component.
Definition: boolVectorI.H:112
Differentiation of the mesh data structure.
Definition: deltaBoundary.H:59
tmp< tensorField > cellCenters_d(const label pointI)
vectorField makeFaceCentresAndAreas_d(const pointField &p, const pointField &p_d)
Definition: deltaBoundary.C:72
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
T getOrDefaultCompat(const word &keyword, std::initializer_list< std::pair< const char *, int > > compat, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
static const char *const coordinates
The keyword "coordinates".
Definition: ensightFile.H:88
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
pointField points(const UList< point > &pts) const
Return the points corresponding to this face.
Definition: faceI.H:87
A class for handling file names.
Definition: fileName.H:76
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:68
localIOdictionary is derived from IOdictionary but excludes parallel master reading.
Smooth ATC in cells having a point to a set of patches supplied by type.
Definition: pointCells.H:59
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:75
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:364
splitCell * master() const
Definition: splitCell.H:113
A class for managing temporary objects.
Definition: tmp.H:65
T & ref() const
Definition: tmpI.H:227
const Switch & bound() const
Return the bound flag.
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
bool
Definition: EEqn.H:20
dynamicFvMesh & mesh
PtrList< coordinateSystem > coordinates(solidRegions.size())
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition: error.H:478
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
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.
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:47
const dimensionSet dimless
Dimensionless.
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:536
List< label > labelList
A List of labels.
Definition: List.H:66
bool mkDir(const fileName &pathName, mode_t mode=0777)
Make a directory and return an error if it could not be created.
Definition: MSwindows.C:515
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:35
messageStream Info
Information stream (stdout output on master, null elsewhere)
GeometricField< vector, pointPatchField, pointMesh > pointVectorField
vector point
Point is a vector.
Definition: point.H:43
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
scalar diff(const triad &A, const triad &B)
Return a quantity of the difference between two triads.
Definition: triad.C:378
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
Field< vector > vectorField
Specialisation of Field<T> for vector.
IOerror FatalIOError
List< bool > boolList
A List of bools.
Definition: List.H:64
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Field< tensor > tensorField
Specialisation of Field<T> for tensor.
error FatalError
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
void cmptMag(FieldField< Field, Type > &cf, const FieldField< Field, Type > &f)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Type gMax(const FieldField< Field, Type > &f)
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
A non-counting (dummy) refCount.
Definition: refCount.H:59