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-2020 PCOpt/NTUA
9  Copyright (C) 2013-2020 FOSS GP
10  Copyright (C) 2019-2021 OpenCFD Ltd.
11 -------------------------------------------------------------------------------
12 License
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 
41 namespace Foam
42 {
43  defineTypeNameAndDebug(NURBS3DVolume, 0);
44  defineRunTimeSelectionTable(NURBS3DVolume, dictionary);
45 }
46 
47 
48 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
49 
51 {
52  // It is considered an error to recompute points in the control boxes
53  if (mapPtr_ || reverseMapPtr_)
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 
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  (
123  new pointVectorField
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
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 :
634  mesh_(mesh),
635  dict_(dict),
636  name_(dict.dictName()),
637  basisU_(dict.get<label>("nCPsU"), dict.get<label>("degreeU")),
638  basisV_(dict.get<label>("nCPsV"), dict.get<label>("degreeV")),
639  basisW_(dict.get<label>("nCPsW"), dict.get<label>("degreeW")),
640  maxIter_(dict.getOrDefault<label>("maxIterations", 10)),
641  tolerance_(dict.getOrDefault<scalar>("tolerance", 1.e-10)),
642  nMaxBound_(dict.getOrDefault<scalar>("nMaxBoundIterations", 4)),
643  cps_(0),
644  mapPtr_(nullptr),
645  reverseMapPtr_(nullptr),
646  parametricCoordinatesPtr_(nullptr),
647  localSystemCoordinates_(mesh_.nPoints(), Zero),
648  confineUMovement_
649  (
651  (
652  "confineUMovement", {{"confineX1movement", 1912}}, false
653  )
654  ),
655  confineVMovement_
656  (
658  (
659  "confineVMovement", {{"confineX2movement", 1912}}, false
660  )
661  ),
662  confineWMovement_
663  (
665  (
666  "confineWMovement", {{"confineX3movement", 1912}}, false
667  )
668  ),
669  confineBoundaryControlPoints_
670  (
671  dict.getOrDefault<bool>("confineBoundaryControlPoints", true)
672  ),
673  confineUMinCPs_
674  (
675  dict.getOrDefaultCompat<boolVectorList>
676  (
677  "confineUMinCPs", {{"boundUMinCPs", 1912}}, boolVectorList()
678  )
679  ),
680  confineUMaxCPs_
681  (
682  dict.getOrDefaultCompat<boolVectorList>
683  (
684  "confineUMaxCPs", {{"boundUMaxCPs", 1912}}, boolVectorList()
685  )
686  ),
687  confineVMinCPs_
688  (
689  dict.getOrDefaultCompat<boolVectorList>
690  (
691  "confineVMinCPs", {{"boundVMinCPs", 1912}}, boolVectorList()
692  )
693  ),
694  confineVMaxCPs_
695  (
696  dict.getOrDefaultCompat<boolVectorList>
697  (
698  "confineVMaxCPs", {{"boundVMaxCPs", 1912}}, boolVectorList()
699  )
700  ),
701  confineWMinCPs_
702  (
703  dict.getOrDefaultCompat<boolVectorList>
704  (
705  "confineWMinCPs", {{"boundWMinCPs", 1912}}, boolVectorList()
706  )
707  ),
708  confineWMaxCPs_
709  (
710  dict.getOrDefaultCompat<boolVectorList>
711  (
712  "confineWMaxCPs", {{"boundWMaxCPs", 1912}}, boolVectorList()
713  )
714  ),
715  activeControlPoints_(0), //zero here, execute sanity checks first
716  activeDesignVariables_(0), //zero here, execute sanity checks first
717  cpsFolder_("controlPoints"),
718  readStoredData_(dict.getOrDefault<bool>("readStoredData", true))
719 {
720  // Create folders
721  makeFolders();
722 
723  // Sanity checks
724  if
725  (
726  (confineUMinCPs_.size() + confineUMaxCPs_.size() >= basisU_.nCPs())
727  || (confineVMinCPs_.size() + confineVMaxCPs_.size() >= basisV_.nCPs())
728  || (confineWMinCPs_.size() + confineWMaxCPs_.size() >= basisW_.nCPs())
729  )
730  {
732  << "Number of control point slices to be kept frozen at "
733  << "the boundaries is invalid \n"
734  << "Number of control points in u " << basisU_.nCPs() << "\n"
735  << "Number of control points in v " << basisV_.nCPs() << "\n"
736  << "Number of control points in w " << basisW_.nCPs() << "\n"
737  << exit(FatalError);
738  }
739 
740  // Define control points
742  determineActiveDesignVariablesAndPoints();
743  writeCpsInDict();
744 }
745 
746 
747 // * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
748 
750 (
751  const dictionary& dict,
752  const fvMesh& mesh,
753  bool computeParamCoors
754 )
755 {
756  const word modelType(dict.get<word>("type"));
757 
758  Info<< "NURBS3DVolume type : " << modelType << endl;
759 
760  auto* ctorPtr = dictionaryConstructorTable(modelType);
761 
762  if (!ctorPtr)
763  {
765  (
766  dict,
767  "type",
768  modelType,
769  *dictionaryConstructorTablePtr_
770  ) << exit(FatalIOError);
771  }
772 
773  return autoPtr<NURBS3DVolume>(ctorPtr(dict, mesh, computeParamCoors));
774 }
775 
776 
777 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
778 
780 (
781  const scalar u,
782  const scalar v,
783  const scalar w
784 ) const
785 {
786  const label degreeU = basisU_.degree();
787  const label degreeV = basisV_.degree();
788  const label degreeW = basisW_.degree();
789 
790  const label nCPsU = basisU_.nCPs();
791  const label nCPsV = basisV_.nCPs();
792  const label nCPsW = basisW_.nCPs();
793 
794  vector derivative(Zero);
795 
796  for (label iCPw = 0; iCPw < nCPsW; ++iCPw)
797  {
798  const scalar basisW(basisW_.basisValue(iCPw, degreeW, w));
799  for (label iCPv = 0; iCPv < nCPsV; ++iCPv)
800  {
801  const scalar basisVW = basisW*basisV_.basisValue(iCPv, degreeV, v);
802  for (label iCPu = 0; iCPu < nCPsU; ++iCPu)
803  {
804  derivative +=
805  cps_[getCPID(iCPu, iCPv, iCPw)]
806  *basisU_.basisDerivativeU(iCPu, degreeU, u)
807  *basisVW;
808  }
809  }
810  }
811 
812  return derivative;
813 }
814 
815 
817 (
818  const scalar u,
819  const scalar v,
820  const scalar w
821 ) const
822 {
823  const label degreeU = basisU_.degree();
824  const label degreeV = basisV_.degree();
825  const label degreeW = basisW_.degree();
826 
827  const label nCPsU = basisU_.nCPs();
828  const label nCPsV = basisV_.nCPs();
829  const label nCPsW = basisW_.nCPs();
830 
831  vector derivative(Zero);
832 
833  for (label iCPw = 0; iCPw < nCPsW; ++iCPw)
834  {
835  const scalar basisW(basisW_.basisValue(iCPw, degreeW, w));
836  for (label iCPv = 0; iCPv < nCPsV; ++iCPv)
837  {
838  const scalar basisWDeriV =
839  basisW*basisV_.basisDerivativeU(iCPv, degreeV, v);
840  for (label iCPu = 0; iCPu < nCPsU; ++iCPu)
841  {
842  derivative +=
843  cps_[getCPID(iCPu, iCPv, iCPw)]
844  *basisU_.basisValue(iCPu, degreeU, u)
845  *basisWDeriV;
846  }
847  }
848  }
849 
850  return derivative;
851 }
852 
853 
855 (
856  const scalar u,
857  const scalar v,
858  const scalar w
859 ) const
860 {
861  const label degreeU = basisU_.degree();
862  const label degreeV = basisV_.degree();
863  const label degreeW = basisW_.degree();
864 
865  const label nCPsU = basisU_.nCPs();
866  const label nCPsV = basisV_.nCPs();
867  const label nCPsW = basisW_.nCPs();
868 
869  vector derivative(Zero);
870 
871  for (label iCPw = 0; iCPw < nCPsW; iCPw++)
872  {
873  const scalar derivW(basisW_.basisDerivativeU(iCPw, degreeW, w));
874  for (label iCPv = 0; iCPv < nCPsV; iCPv++)
875  {
876  const scalar derivWBasisV =
877  derivW*basisV_.basisValue(iCPv, degreeV, v);
878  for (label iCPu = 0; iCPu < nCPsU; iCPu++)
879  {
880  derivative +=
881  cps_[getCPID(iCPu, iCPv, iCPw)]
882  *basisU_.basisValue(iCPu, degreeU, u)
883  *derivWBasisV;
884  }
885  }
886  }
887 
888  return derivative;
889 }
890 
891 
893 (
894  const vector& uVector
895 ) const
896 {
897  const scalar u = uVector.x();
898  const scalar v = uVector.y();
899  const scalar w = uVector.z();
900 
901  vector uDeriv = volumeDerivativeU(u, v, w);
902  vector vDeriv = volumeDerivativeV(u, v, w);
903  vector wDeriv = volumeDerivativeW(u, v, w);
904 
905  tensor Jacobian(Zero);
906 
907  Jacobian[0] = uDeriv.component(0);
908  Jacobian[1] = vDeriv.component(0);
909  Jacobian[2] = wDeriv.component(0);
910  Jacobian[3] = uDeriv.component(1);
911  Jacobian[4] = vDeriv.component(1);
912  Jacobian[5] = wDeriv.component(1);
913  Jacobian[6] = uDeriv.component(2);
914  Jacobian[7] = vDeriv.component(2);
915  Jacobian[8] = wDeriv.component(2);
916 
917  return Jacobian;
918 }
919 
920 
922 (
923  const vector& uVector,
924  const label cpI
925 ) const
926 {
927  const scalar u = uVector.x();
928  const scalar v = uVector.y();
929  const scalar w = uVector.z();
930 
931  const label nCPsU = basisU_.nCPs();
932  const label nCPsV = basisV_.nCPs();
933 
934  const label degreeU = basisU_.degree();
935  const label degreeV = basisV_.degree();
936  const label degreeW = basisW_.degree();
937 
938  label iCPw = cpI/label(nCPsU*nCPsV);
939  label iCPv = (cpI - iCPw*nCPsU*nCPsV)/nCPsU;
940  label iCPu = (cpI - iCPw*nCPsU*nCPsV - iCPv*nCPsU);
941 
942  // Normally, this should be a tensor, however the parameterization is
943  // isotropic. Hence the tensor degenerates to a diagonal tensor with all
944  // diagonal elements being equal. This returns the (unique) diag element
945  scalar derivative =
946  basisU_.basisValue(iCPu, degreeU, u)
947  *basisV_.basisValue(iCPv, degreeV, v)
948  *basisW_.basisValue(iCPw, degreeW, w);
949 
950  return derivative;
951 }
952 
953 
955 (
956  const pointVectorField& pointSens,
957  const labelList& sensitivityPatchIDs
958 )
959 {
960  vectorField controlPointDerivs(cps_.size(), Zero);
961 
962  // Get parametric coordinates
963  const vectorField& parametricCoordinates = getParametricCoordinates();
964 
965  forAll(controlPointDerivs, cpI)
966  {
967  forAll(sensitivityPatchIDs, pI)
968  {
969  const label patchI = sensitivityPatchIDs[pI];
970  const polyPatch& patch = mesh_.boundaryMesh()[patchI];
971  const labelList& meshPoints = patch.meshPoints();
972 
973  forAll(meshPoints, mpI)
974  {
975  const label globalIndex = meshPoints[mpI];
976  const label whichPointInBox = reverseMapPtr_()[globalIndex];
977 
978  // If point resides within control points box,
979  // add contribution to cp derivative
980  if (whichPointInBox != -1)
981  {
982  controlPointDerivs[cpI] +=
983  (
984  pointSens[globalIndex]
985  & transformationTensorDxDb(globalIndex)
986  )
987  *volumeDerivativeCP
988  (
989  parametricCoordinates[globalIndex],
990  cpI
991  );
992  }
993  }
994  }
995  }
996 
997  // Sum contributions from all processors
998  Pstream::listCombineGather(controlPointDerivs, plusEqOp<vector>());
999  Pstream::listCombineScatter(controlPointDerivs);
1000 
1001  return controlPointDerivs;
1002 }
1003 
1004 
1007  const volVectorField& faceSens,
1008  const labelList& sensitivityPatchIDs
1009 )
1010 {
1011  return
1012  computeControlPointSensitivities
1013  (
1014  faceSens.boundaryField(),
1015  sensitivityPatchIDs
1016  );
1017 }
1018 
1019 
1022  const boundaryVectorField& faceSens,
1023  const labelList& sensitivityPatchIDs
1024 )
1025 {
1026  // Return field
1027  vectorField controlPointDerivs(cps_.size(), Zero);
1028 
1029  // Get parametric coordinates
1030  const vectorField& parametricCoordinates = getParametricCoordinates();
1031 
1032  // Auxiliary quantities
1034  const labelList& reverseMap = reverseMapPtr_();
1035 
1036  forAll(controlPointDerivs, cpI)
1037  {
1038  forAll(sensitivityPatchIDs, pI)
1039  {
1040  const label patchI = sensitivityPatchIDs[pI];
1041  const polyPatch& patch = mesh_.boundaryMesh()[patchI];
1042  const label patchStart = patch.start();
1043  const fvPatchVectorField& patchSens = faceSens[patchI];
1044 
1045  // loop over patch faces
1046  forAll(patch, fI)
1047  {
1048  const face& fGlobal = mesh_.faces()[fI + patchStart];
1049  const pointField facePoints = fGlobal.points(mesh_.points());
1050  // loop over face points
1051  tensorField facePointDerivs(facePoints.size(), Zero);
1052  forAll(fGlobal, pI)
1053  {
1054  const label globalIndex = fGlobal[pI]; //global point index
1055  const label whichPointInBox = reverseMap[globalIndex];
1056  // if point resides within control points box,
1057  // add contribution to d( facePoints )/db
1058  if (whichPointInBox != -1)
1059  {
1060  // TENSOR-BASED
1061  //~~~~~~~~~~~~~
1062  facePointDerivs[pI] =
1063  transformationTensorDxDb(globalIndex)
1064  * volumeDerivativeCP
1065  (
1066  parametricCoordinates[globalIndex],
1067  cpI
1068  );
1069 
1070  }
1071  }
1072 
1073  tensor fCtrs_d =
1075  (
1076  facePoints,
1077  facePointDerivs
1078  )[0];
1079  controlPointDerivs[cpI] += patchSens[fI] & fCtrs_d;
1080  }
1081  }
1082  }
1083  // Sum contributions from all processors
1084  Pstream::listCombineGather(controlPointDerivs, plusEqOp<vector>());
1085  Pstream::listCombineScatter(controlPointDerivs);
1086 
1087  return controlPointDerivs;
1088 }
1089 
1090 
1093  const vectorField& faceSens,
1094  const label patchI,
1095  const label cpI
1096 )
1097 {
1098  // Return vector
1099  vector cpSens(Zero);
1100  // Get parametric coordinates
1101  const vectorField& parametricCoordinates = getParametricCoordinates();
1102 
1103  // Auxiliary quantities
1105  const labelList& reverseMap = reverseMapPtr_();
1106 
1107  const polyPatch& patch = mesh_.boundaryMesh()[patchI];
1108  const label patchStart = patch.start();
1109  // Loop over patch faces
1110  forAll(patch, fI)
1111  {
1112  const face& fGlobal = mesh_.faces()[fI + patchStart];
1113  const pointField facePoints = fGlobal.points(mesh_.points());
1114  // Loop over face points
1115  tensorField facePointDerivs(facePoints.size(), Zero);
1116  forAll(fGlobal, pI)
1117  {
1118  const label globalIndex = fGlobal[pI]; //global point index
1119  const label whichPointInBox = reverseMap[globalIndex];
1120  // If point resides within control points box,
1121  // add contribution to d( facePoints )/db
1122  if (whichPointInBox != -1)
1123  {
1124  // TENSOR-BASED
1125  //~~~~~~~~~~~~~
1126  facePointDerivs[pI] =
1127  transformationTensorDxDb(globalIndex)
1128  *volumeDerivativeCP
1129  (
1130  parametricCoordinates[globalIndex],
1131  cpI
1132  );
1133  }
1134  }
1135 
1136  tensor fCtrs_d =
1138  (
1139  facePoints,
1140  facePointDerivs
1141  )[0];
1142  cpSens += faceSens[fI] & fCtrs_d;
1143  }
1144  // Sum contributions from all processors
1145  reduce(cpSens, sumOp<vector>());
1146 
1147  return cpSens;
1148 }
1149 
1150 
1153  const label patchI,
1154  const label cpI,
1155  bool DimensionedNormalSens
1156 )
1157 {
1158  const fvPatch& patch = mesh_.boundary()[patchI];
1159  const polyPatch& ppatch = patch.patch();
1160  // Return field
1161  tmp<tensorField> tdndbSens(new tensorField(patch.size(), Zero));
1162  tensorField& dndbSens = tdndbSens.ref();
1163  // Auxiliary quantities
1165  const label patchStart = ppatch.start();
1166  const labelList& reverseMap = reverseMapPtr_();
1167 
1168  // Get parametric coordinates
1169  const vectorField& parametricCoordinates = getParametricCoordinates();
1170 
1171  // Loop over patch faces
1172  forAll(patch, fI)
1173  {
1174  const face& fGlobal = mesh_.faces()[fI + patchStart];
1175  const pointField facePoints = fGlobal.points(mesh_.points());
1176  // Loop over face points
1177  tensorField facePointDerivs(facePoints.size(), Zero);
1178  forAll(fGlobal, pI)
1179  {
1180  const label globalIndex = fGlobal[pI]; //global point index
1181  const label whichPointInBox = reverseMap[globalIndex];
1182  // If point resides within control points box,
1183  // add contribution to d( facePoints )/db
1184  if (whichPointInBox != -1)
1185  {
1186  // TENSOR-BASED
1187  //~~~~~~~~~~~~~
1188  facePointDerivs[pI] =
1189  transformationTensorDxDb(globalIndex)
1190  *volumeDerivativeCP
1191  (
1192  parametricCoordinates[globalIndex],
1193  cpI
1194  );
1195  }
1196  }
1197 
1198  // Determine whether to return variance of dimensioned or unit normal
1199  tensorField dNdbSens =
1201  (
1202  facePoints,
1203  facePointDerivs
1204  );
1205 
1206  if (DimensionedNormalSens)
1207  {
1208  dndbSens[fI] = dNdbSens[1];
1209  }
1210  else
1211  {
1212  dndbSens[fI] = dNdbSens[2];
1213  }
1214  }
1215 
1216  return tdndbSens;
1217 }
1218 
1219 
1222  const label patchI,
1223  const label cpI
1224 )
1225 {
1226  // Get parametric coordinates
1227  const vectorField& parametricCoordinates = getParametricCoordinates();
1228 
1229  // Patch data
1230  const polyPatch& patch = mesh_.boundaryMesh()[patchI];
1231  const labelList& meshPoints = patch.meshPoints();
1232 
1233  // Return field
1234  auto tdxdb = tmp<tensorField>::New(patch.nPoints(), Zero);
1235  auto& dxdb = tdxdb.ref();
1236 
1237  forAll(meshPoints, pI)
1238  {
1239  const label globalIndex = meshPoints[pI]; //global point index
1240  const label whichPointInBox = reverseMapPtr_()[globalIndex];
1241 
1242  // If point resides within control points box, find dxdb
1243  if (whichPointInBox != -1)
1244  {
1245  dxdb[pI] =
1246  transformationTensorDxDb(globalIndex)
1247  *volumeDerivativeCP
1248  (
1249  parametricCoordinates[globalIndex],
1250  cpI
1251  );
1252  }
1253  }
1254 
1255  return tdxdb;
1256 }
1257 
1258 
1261  const label patchI,
1262  const label cpI
1263 )
1264 {
1265  // get parametric coordinates
1266  const vectorField& parametricCoordinates = getParametricCoordinates();
1267 
1268  // Patch data
1269  const polyPatch& patch = mesh_.boundaryMesh()[patchI];
1270  const label patchStart = patch.start();
1271 
1272  // Return field
1273  auto tdxdb = tmp<tensorField>::New(patch.size(), Zero);
1274  auto& dxdb = tdxdb.ref();
1275 
1276  // Mesh differentiation engine
1277  deltaBoundary deltaBound(mesh_);
1278 
1279  forAll(patch, fI)
1280  {
1281  const face& fGlobal = mesh_.faces()[fI + patchStart];
1282  const pointField facePoints = fGlobal.points(mesh_.points());
1283  // Loop over face points
1284  tensorField facePointDerivs(facePoints.size(), Zero);
1285  forAll(fGlobal, pI)
1286  {
1287  const label globalIndex = fGlobal[pI]; //global point index
1288  const label whichPointInBox = reverseMapPtr_()[globalIndex];
1289  // If point resides within control points box,
1290  // add contribution to d( facePoints )/db
1291  if (whichPointInBox != -1)
1292  {
1293  // TENSOR-BASED
1294  //~~~~~~~~~~~~~
1295  facePointDerivs[pI] =
1296  transformationTensorDxDb(globalIndex)
1297  *volumeDerivativeCP
1298  (
1299  parametricCoordinates[globalIndex],
1300  cpI
1301  );
1302 
1303  }
1304  }
1305  dxdb[fI] =
1306  deltaBound.makeFaceCentresAndAreas_d
1307  (
1308  facePoints,
1309  facePointDerivs
1310  )[0];
1311  }
1312 
1313  return tdxdb;
1314 }
1315 
1316 
1319  const vectorField& uVector
1320 ) const
1321 {
1322  const label nPoints = mapPtr_().size();
1323  auto tpoints = tmp<vectorField>::New(nPoints, Zero);
1324  auto& points = tpoints.ref();
1325 
1326  forAll(points, pI)
1327  {
1328  const label globalPI = mapPtr_()[pI];
1329  points[pI] = coordinates(uVector[globalPI]);
1330  }
1331 
1332  return tpoints;
1333 }
1334 
1335 
1338  const vector& uVector
1339 ) const
1340 {
1341  const label degreeU = basisU_.degree();
1342  const label degreeV = basisV_.degree();
1343  const label degreeW = basisW_.degree();
1344 
1345  const label nCPsU = basisU_.nCPs();
1346  const label nCPsV = basisV_.nCPs();
1347  const label nCPsW = basisW_.nCPs();
1348 
1349  const scalar u = uVector.x();
1350  const scalar v = uVector.y();
1351  const scalar w = uVector.z();
1352 
1353  vector point(Zero);
1354  for (label iCPw = 0; iCPw < nCPsW; iCPw++)
1355  {
1356  const scalar basisW(basisW_.basisValue(iCPw, degreeW, w));
1357  for (label iCPv = 0; iCPv < nCPsV; iCPv++)
1358  {
1359  const scalar basisVW =
1360  basisW*basisV_.basisValue(iCPv, degreeV, v);
1361  for (label iCPu = 0; iCPu < nCPsU; iCPu++)
1362  {
1363  point +=
1364  cps_[getCPID(iCPu, iCPv, iCPw)]
1365  *basisU_.basisValue(iCPu, degreeU, u)
1366  *basisVW;
1367  }
1368  }
1369  }
1370 
1371  return point;
1372 }
1373 
1374 
1377  const vectorField& controlPointsMovement
1378 )
1379 {
1380  // Get parametric coordinates and map
1381  const vectorField& paramCoors = getParametricCoordinates();
1382  const labelList& map = mapPtr_();
1383 
1384  // Update control points position
1385  cps_ += controlPointsMovement;
1386  writeCps("cpsBsplines"+mesh_.time().timeName());
1387  writeCpsInDict();
1388 
1389  // Compute new mesh points based on updated control points
1390  tmp<vectorField> tparameterizedPoints = coordinates(paramCoors);
1391  const vectorField& parameterizedPoints = tparameterizedPoints();
1392 
1393  // Return field. Initialized with current mesh points
1394  tmp<vectorField> tnewPoints(new vectorField(mesh_.points()));
1395  vectorField& newPoints = tnewPoints.ref();
1396 
1397  // Update position of parameterized points
1398  forAll(parameterizedPoints, pI)
1399  {
1400  newPoints[map[pI]] = transformPointToCartesian(parameterizedPoints[pI]);
1401  }
1402 
1403  // Update coordinates in the local system based on the cartesian points
1404  updateLocalCoordinateSystem(newPoints);
1405  DebugInfo
1406  << "Max mesh movement equal to "
1407  << gMax(mag(newPoints - mesh_.points())) << endl;
1408 
1409  return tnewPoints;
1410 }
1411 
1412 
1415  const vectorField& controlPointsMovement,
1416  const labelList& patchesToBeMoved
1417 )
1418 {
1419  // Get parametric coordinates
1420  const vectorField& paramCoors = getParametricCoordinates();
1421 
1422  // Update control points position
1423  cps_ += controlPointsMovement;
1424 
1425  writeCps("cpsBsplines"+mesh_.time().timeName());
1426  writeCpsInDict();
1427 
1428  // Return field. Initialized with current mesh points
1429  tmp<vectorField> tnewPoints(new vectorField(mesh_.points()));
1430  vectorField& newPoints = tnewPoints.ref();
1431 
1432  // Update position of parameterized boundary points
1433  for (const label patchI : patchesToBeMoved)
1434  {
1435  const polyPatch& patch = mesh_.boundaryMesh()[patchI];
1436  const labelList& meshPoints = patch.meshPoints();
1437 
1438  for (const label globalIndex : meshPoints)
1439  {
1440  const label whichPointInBox = reverseMapPtr_()[globalIndex];
1441  // If point resides within control points box,
1442  // compute new cartesian coordinates
1443  if (whichPointInBox != -1)
1444  {
1445  newPoints[globalIndex] =
1446  transformPointToCartesian
1447  (
1448  coordinates
1449  (
1450  paramCoors[globalIndex]
1451  )
1452  );
1453  }
1454  }
1455  }
1456 
1457  // Update coordinates in the local system based on the cartesian points
1458  updateLocalCoordinateSystem(newPoints);
1459  DebugInfo
1460  << "Max mesh movement equal to "
1461  << gMax(mag(newPoints - mesh_.points())) << endl;
1462 
1463  return tnewPoints;
1464 }
1465 
1466 
1467 Foam::label Foam::NURBS3DVolume::getCPID
1469  const label i,
1470  const label j,
1471  const label k
1472 ) const
1473 {
1474  const label nCPsU = basisU_.nCPs();
1475  const label nCPsV = basisV_.nCPs();
1476 
1477  return k*nCPsU*nCPsV + j*nCPsU + i;
1478 }
1479 
1480 
1482 {
1483  if (cps_.size() != newCps.size())
1484  {
1486  << "Attempting to replace control points with a set of "
1487  << "different size"
1488  << exit(FatalError);
1489  }
1490  cps_ = newCps;
1491 }
1492 
1493 
1496  vectorField& controlPointsMovement
1497 )
1498 {
1499  forAll(controlPointsMovement, cpI)
1500  {
1501  if (!activeDesignVariables_[3*cpI])
1502  {
1503  controlPointsMovement[cpI].x() = Zero;
1504  }
1505  if (!activeDesignVariables_[3*cpI + 1])
1506  {
1507  controlPointsMovement[cpI].y() = Zero;
1508  }
1509  if (!activeDesignVariables_[3*cpI + 2])
1510  {
1511  controlPointsMovement[cpI].z() = Zero;
1512  }
1513  }
1514 }
1515 
1516 
1519  const vectorField& controlPointsMovement,
1520  const labelList& patchesToBeMoved
1521 )
1522 {
1523  // Backup old cps
1524  vectorField oldCPs = cps_;
1525  // Get parametric coordinates
1526  const vectorField& paramCoors = getParametricCoordinates();
1527  // Update control points position
1528  cps_ += controlPointsMovement;
1529  // Update position of parameterized boundary points
1530  scalar maxDisplacement(Zero);
1531  for (const label patchI : patchesToBeMoved)
1532  {
1533  const polyPatch& patch = mesh_.boundaryMesh()[patchI];
1534  const labelList& meshPoints = patch.meshPoints();
1535 
1536  for (const label globalIndex : meshPoints)
1537  {
1538  const label whichPointInBox = reverseMapPtr_()[globalIndex];
1539  // If point resides within control points box,
1540  // compute new cartesian coordinates
1541  if (whichPointInBox != -1)
1542  {
1543  vector newPoint =
1544  transformPointToCartesian
1545  (
1546  coordinates
1547  (
1548  paramCoors[globalIndex]
1549  )
1550  );
1551  maxDisplacement =
1552  max
1553  (
1554  maxDisplacement,
1555  mag(newPoint - mesh_.points()[globalIndex])
1556  );
1557  }
1558  }
1559  }
1560  reduce(maxDisplacement, maxOp<scalar>());
1561  cps_ = oldCPs;
1562 
1563  return maxDisplacement;
1564 }
1565 
1566 
1568 {
1569  if (!mapPtr_)
1570  {
1571  findPointsInBox(localSystemCoordinates_);
1572  }
1573  tmp<vectorField> pointsInBox
1574  (
1575  new vectorField(localSystemCoordinates_, mapPtr_())
1576  );
1577 
1578  return pointsInBox;
1579 }
1580 
1581 
1583 {
1584  if (!mapPtr_)
1585  {
1586  findPointsInBox(localSystemCoordinates_);
1587  }
1588 
1589  return mapPtr_();
1590 }
1591 
1592 
1594 {
1595  if (!reverseMapPtr_)
1596  {
1597  findPointsInBox(localSystemCoordinates_);
1598  }
1599 
1600  return reverseMapPtr_();
1601 }
1602 
1603 
1605 {
1606  // If not computed yet, compute parametric coordinates
1607  if (!parametricCoordinatesPtr_)
1608  {
1609  // Find mesh points inside control points box
1610  // if they have been identified yet
1611  if (!mapPtr_)
1612  {
1613  findPointsInBox(localSystemCoordinates_);
1614  }
1615  computeParametricCoordinates(getPointsInBox()());
1616  }
1617 
1618  return parametricCoordinatesPtr_();
1619 }
1620 
1621 
1623 {
1624  // Get parametric coordinates
1625  const vectorField& parametricCoordinates = getParametricCoordinates();
1626 
1627  // Set return field to zero
1628  tmp<pointTensorField> tDxDb
1629  (
1630  new pointTensorField
1631  (
1632  IOobject
1633  (
1634  "DxDb",
1635  mesh_.time().timeName(),
1636  mesh_,
1639  ),
1640  pointMesh::New(mesh_),
1642  )
1643  );
1644 
1645  pointTensorField& DxDb = tDxDb.ref();
1646 
1647  // All points outside the control box remain unmoved.
1648  // Loop over only points within the control box
1649  const labelList& map = mapPtr_();
1650  for (const label globalIndex : map)
1651  {
1652  DxDb[globalIndex] =
1653  transformationTensorDxDb(globalIndex)
1654  *volumeDerivativeCP
1655  (
1656  parametricCoordinates[globalIndex],
1657  cpI
1658  );
1659  }
1660 
1661  return tDxDb;
1662 }
1663 
1664 
1667  const label cpI
1668 )
1669 {
1670  // Get parametric coordinates
1671  const vectorField& parametricCoordinates = getParametricCoordinates();
1672 
1673  // Set return field to zero
1674  tmp<volTensorField> tDxDb
1675  (
1676  new volTensorField
1677  (
1678  IOobject
1679  (
1680  "DxDb",
1681  mesh_.time().timeName(),
1682  mesh_,
1685  ),
1686  mesh_,
1688  )
1689  );
1690 
1691  volTensorField& DxDb = tDxDb.ref();
1692  deltaBoundary deltaBound(mesh_);
1693  const labelListList& pointCells = mesh_.pointCells();
1694 
1695  // All points outside the control box remain unmoved.
1696  // Loop over only points within the control box
1697  const labelList& map = mapPtr_();
1698  for (const label globalIndex : map)
1699  {
1700  tensor pointDxDb =
1701  transformationTensorDxDb(globalIndex)
1702  *volumeDerivativeCP
1703  (
1704  parametricCoordinates[globalIndex],
1705  cpI
1706  );
1707  const labelList& pointCellsI = pointCells[globalIndex];
1708  tmp<tensorField> tC_d = deltaBound.cellCenters_d(globalIndex);
1709  const tensorField& C_d = tC_d();
1710  forAll(pointCellsI, cI)
1711  {
1712  const label cellI = pointCellsI[cI];
1713  DxDb[cellI] += C_d[cI] & pointDxDb;
1714  }
1715  }
1716 
1717  // Assign boundary values since the grad of this field is often needed
1718  forAll(mesh_.boundary(), pI)
1719  {
1720  const fvPatch& patch = mesh_.boundary()[pI];
1721  if (!isA<coupledFvPatch>(patch))
1722  {
1723  DxDb.boundaryFieldRef()[pI] = patchDxDbFace(pI, cpI);
1724  }
1725  }
1726 
1727  // Correct coupled boundaries
1729 
1730  return tDxDb;
1731 }
1732 
1733 
1735 {
1736  label nU(basisU_.nCPs());
1737  if (nU % 2 == 0)
1738  {
1739  nU /=2;
1740  }
1741  else
1742  {
1743  nU = (nU - 1)/2 + 1;
1744  }
1745  return nU;
1746 }
1747 
1748 
1750 {
1751  label nV(basisV_.nCPs());
1752  if (nV % 2 == 0)
1753  {
1754  nV /=2;
1755  }
1756  else
1757  {
1758  nV = (nV - 1)/2 + 1;
1759  }
1760  return nV;
1761 }
1762 
1763 
1765 {
1766  label nW(basisW_.nCPs());
1767  if (nW % 2 == 0)
1768  {
1769  nW /=2;
1770  }
1771  else
1772  {
1773  nW = (nW - 1)/2 + 1;
1774  }
1775  return nW;
1776 }
1777 
1778 
1781  const fileName& baseName,
1782  const bool transform
1783 ) const
1784 {
1785  const label nCPsU = basisU_.nCPs();
1786  const label nCPsV = basisV_.nCPs();
1787 
1788  vectorField cpsInCartesian(cps_);
1789  if (transform)
1790  {
1791  forAll(cpsInCartesian, cpI)
1792  {
1793  cpsInCartesian[cpI] = transformPointToCartesian(cps_[cpI]);
1794  }
1795  }
1796 
1797  Info<< "Writing control point positions to file" << endl;
1798 
1799  if (Pstream::master())
1800  {
1801  OFstream cpsFile("optimisation"/cpsFolder_/name_ + baseName + ".csv");
1802  // Write header
1803  cpsFile
1804  << "\"Points : 0\", \"Points : 1\", \"Points : 2\","
1805  << "\"i\", \"j\", \"k\","
1806  << "\"active : 0\", \"active : 1\", \"active : 2\"" << endl;
1807 
1808  forAll(cpsInCartesian, cpI)
1809  {
1810  const label iCPw = cpI/label(nCPsU*nCPsV);
1811  const label iCPv = (cpI - iCPw*nCPsU*nCPsV)/nCPsU;
1812  const label iCPu = (cpI - iCPw*nCPsU*nCPsV - iCPv*nCPsU);
1813 
1814  cpsFile
1815  << cpsInCartesian[cpI].x() << ", "
1816  << cpsInCartesian[cpI].y() << ", "
1817  << cpsInCartesian[cpI].z() << ", "
1818  << iCPu << ", "
1819  << iCPv << ", "
1820  << iCPw << ", "
1821  << activeDesignVariables_[3*cpI] << ", "
1822  << activeDesignVariables_[3*cpI + 1] << ", "
1823  << activeDesignVariables_[3*cpI + 2] << endl;
1824  }
1825  }
1826 }
1827 
1828 
1830 {
1831  IOdictionary cpsDict
1832  (
1833  IOobject
1834  (
1835  name_ + "cpsBsplines" + mesh_.time().timeName(),
1836  mesh_.time().caseConstant(),
1837  cpsFolder_,
1838  mesh_,
1841  false
1842  )
1843  );
1844 
1845  cpsDict.add("controlPoints", cps_);
1846 
1847  // Always write in ASCII, but allow compression
1848  cpsDict.regIOobject::writeObject
1849  (
1850  IOstreamOption(IOstream::ASCII, mesh_.time().writeCompression()),
1851  true
1852  );
1853 }
1854 
1855 
1857 {
1858  parametricCoordinatesPtr_().write();
1859 }
1860 
1861 
1862 // ************************************************************************* //
Foam::fvPatchField< vector >
Foam::NURBS3DVolume::computeParametricCoordinates
void computeParametricCoordinates(const vectorField &points)
Definition: NURBS3DVolume.C:105
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
Foam::pointVectorField
GeometricField< vector, pointPatchField, pointMesh > pointVectorField
Definition: pointFieldsFwd.H:61
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:54
Foam::Tensor< scalar >
Foam::maxOp
Definition: ops.H:223
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::Vector::x
const Cmpt & x() const
Access to the vector x component.
Definition: VectorI.H:73
Foam::IOobject::AUTO_WRITE
Definition: IOobject.H:194
Foam::deltaBoundary
Differentiation of the mesh data structure.
Definition: deltaBoundary.H:58
Foam::NURBS3DVolume::nWSymmetry
label nWSymmetry() const
Get number of variables if CPs are moved symmetrically in W.
Definition: NURBS3DVolume.C:1764
Foam::VectorSpace::component
const Cmpt & component(const direction) const
Definition: VectorSpaceI.H:87
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::fileName
A class for handling file names.
Definition: fileName.H:73
Foam::NURBS3DVolume::nUSymmetry
label nUSymmetry() const
Get number of variables if CPs are moved symmetrically in U.
Definition: NURBS3DVolume.C:1734
Foam::NURBS3DVolume::getMap
const labelList & getMap()
Get map of points in box to mesh points.
Definition: NURBS3DVolume.C:1582
Foam::NURBS3DVolume::volumeDerivativeU
vector volumeDerivativeU(const scalar u, const scalar v, const scalar w) const
Volume point derivative wrt u at point u,v,w.
Definition: NURBS3DVolume.C:780
Foam::NURBS3DVolume::write
void write() const
Write parametric coordinates.
Definition: NURBS3DVolume.C:1856
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::defineRunTimeSelectionTable
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
Foam::bound
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:35
Foam::NURBS3DVolume::volumeDerivativeW
vector volumeDerivativeW(const scalar u, const scalar v, const scalar w) const
Volume point derivative wrt w at point u,v,w.
Definition: NURBS3DVolume.C:855
Foam::dictionary::getOrDefaultCompat
T getOrDefaultCompat(const word &keyword, std::initializer_list< std::pair< const char *, int >> compat, const T &deflt, enum keyType::option=keyType::REGEX) const
Definition: dictionaryTemplates.C:432
Foam::MeshObject< polyMesh, UpdateableMeshObject, pointMesh >::New
static const pointMesh & New(const polyMesh &mesh, Args &&... args)
Get existing or create a new MeshObject.
Definition: MeshObject.C:48
Foam::tensorField
Field< tensor > tensorField
Specialisation of Field<T> for tensor.
Definition: primitiveFieldsFwd.H:57
Foam::UPstream::master
static bool master(const label communicator=worldComm)
Am I the master process.
Definition: UPstream.H:457
Foam::NURBS3DVolume::computeNewBoundaryPoints
tmp< vectorField > computeNewBoundaryPoints(const vectorField &controlPointsMovement, const labelList &patchesToBeMoved)
Boundary mesh movement based on given control point movement.
Definition: NURBS3DVolume.C:1414
Foam::boolList
List< bool > boolList
A List of bools.
Definition: List.H:65
Foam::FatalIOError
IOerror FatalIOError
Foam::NURBS3DVolume::makeFolders
void makeFolders()
Create folders to store cps and derivatives.
Definition: NURBS3DVolume.C:392
Foam::NURBS3DVolume::bound
bool bound(vector &vec, scalar minValue=1e-7, scalar maxValue=0.999999)
Bound components to certain limits.
Definition: NURBS3DVolume.C:349
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::NURBS3DVolume::confineControlPoint
void confineControlPoint(const label cpI)
Confine all three movements for a prescribed control point.
Definition: NURBS3DVolume.C:585
Foam::NURBS3DVolume::cps_
vectorField cps_
The volumetric B-Splines control points.
Definition: NURBS3DVolume.H:99
Foam::dictionary::get
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:107
Foam::transform
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:521
Foam::Vector::z
const Cmpt & z() const
Access to the vector z component.
Definition: VectorI.H:85
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
Foam::NURBS3DVolume::confineControlPointsDirections
void confineControlPointsDirections()
Confine movement in all control points for user-defined directions.
Definition: NURBS3DVolume.C:574
Foam::boolVector
Specialized bundling of boolean values as a vector of 3 components, element access using x(),...
Definition: boolVector.H:56
Foam::NURBS3DVolume::continuityRealatedConfinement
void continuityRealatedConfinement()
Confine control point movement to maintain user-defined continuity.
Definition: NURBS3DVolume.C:479
Foam::NURBS3DVolume::computeControlPointSensitivities
vectorField computeControlPointSensitivities(const pointVectorField &pointSens, const labelList &sensitivityPatchIDs)
Definition: NURBS3DVolume.C:955
Foam::sumOp
Definition: ops.H:213
Foam::NURBS3DVolume::reverseMapPtr_
autoPtr< labelList > reverseMapPtr_
Map of mesh points to points-in-box.
Definition: NURBS3DVolume.H:106
Foam::NURBS3DVolume::getCPID
label getCPID(const label i, const label j, const label k) const
Get control point ID from its I-J-K coordinates.
Definition: NURBS3DVolume.C:1468
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
coupledFvPatch.H
OFstream.H
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:54
Foam::dimensionedVector
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Definition: dimensionedVector.H:50
nPoints
label nPoints
Definition: gmvOutputHeader.H:2
controlPointsDefinition.H
Foam::NURBS3DVolume::getParametricCoordinates
const pointVectorField & getParametricCoordinates()
Get parametric coordinates.
Definition: NURBS3DVolume.C:1604
Foam::diff
scalar diff(const triad &A, const triad &B)
Return a quantity of the difference between two triads.
Definition: triad.C:378
Foam::NURBS3DVolume::confineBoundaryControlPoints
void confineBoundaryControlPoints()
Confine movement in boundary control points if necessary.
Definition: NURBS3DVolume.C:433
Foam::cmptMag
void cmptMag(FieldField< Field, Type > &cf, const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:400
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:51
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:227
FatalIOErrorInLookup
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition: error.H:478
Foam::Field< vector >
Foam::NURBS3DVolume::findPointsInBox
void findPointsInBox(const vectorField &meshPoints)
Find points within control points box.
Definition: NURBS3DVolume.C:50
Foam::inv
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
Definition: dimensionedSphericalTensor.C:73
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:68
Foam::deltaBoundary::cellCenters_d
tmp< tensorField > cellCenters_d(const label pointI)
Definition: deltaBoundary.C:285
Foam::List::setSize
void setSize(const label n)
Alias for resize()
Definition: List.H:222
coordinates
PtrList< coordinateSystem > coordinates(solidRegions.size())
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
Foam::IOstreamOption
The IOstreamOption is a simple container for options an IOstream can normally have.
Definition: IOstreamOption.H:63
Foam::NURBS3DVolume::determineActiveDesignVariablesAndPoints
void determineActiveDesignVariablesAndPoints()
Create lists with active design variables and control points.
Definition: NURBS3DVolume.C:401
Foam::NURBS3DVolume::mapPtr_
autoPtr< labelList > mapPtr_
Map of points-in-box to mesh points.
Definition: NURBS3DVolume.H:102
Foam::NURBS3DVolume::volumeDerivativeCP
scalar volumeDerivativeCP(const vector &u, const label cpI) const
Volume point derivative wrt to control point cpI at point u,v,w.
Definition: NURBS3DVolume.C:922
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::NURBS3DVolume::nVSymmetry
label nVSymmetry() const
Get number of variables if CPs are moved symmetrically in V.
Definition: NURBS3DVolume.C:1749
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::deltaBoundary::makeFaceCentresAndAreas_d
vectorField makeFaceCentresAndAreas_d(const pointField &p, const pointField &p_d)
Definition: deltaBoundary.C:72
Foam::controlPointsDefinition::New
static autoPtr< controlPointsDefinition > New(NURBS3DVolume &box)
Return a reference to the selected controlPointsDefinition model.
Definition: controlPointsDefinition.C:115
Foam::NURBS3DVolume::writeCpsInDict
void writeCpsInDict() const
Write control points on the local coordinate system.
Definition: NURBS3DVolume.C:1829
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::boolVector::z
bool z() const
The z component.
Definition: boolVectorI.H:112
Foam::globalIndex
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:68
Foam::Field::component
tmp< Field< cmptType > > component(const direction) const
Return a component field of the field.
Definition: Field.C:539
Foam::polyPatch::start
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:361
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::OFstream
Output to file stream, using an OSstream.
Definition: OFstream.H:53
Foam::NURBS3DVolume::getPointsInBox
tmp< vectorField > getPointsInBox()
Get mesh points that reside within the control points box.
Definition: NURBS3DVolume.C:1567
Foam::GeometricField::correctBoundaryConditions
void correctBoundaryConditions()
Correct boundary field.
Definition: GeometricField.C:940
deltaBoundary.H
Foam::Pstream::listCombineGather
static void listCombineGather(const List< commsStruct > &comms, List< T > &Value, const CombineOp &cop, const int tag, const label comm)
Definition: combineGatherScatter.C:290
Time.H
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
Foam::NURBS3DVolume::volumeDerivativeV
vector volumeDerivativeV(const scalar u, const scalar v, const scalar w) const
Volume point derivative wrt v at point u,v,w.
Definition: NURBS3DVolume.C:817
Foam::Vector::y
const Cmpt & y() const
Access to the vector y component.
Definition: VectorI.H:79
Foam::IOstreamOption::ASCII
"ascii" (normal default)
Definition: IOstreamOption.H:72
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::NURBS3DVolume::setControlPoints
void setControlPoints(const vectorField &newCps)
Set new control points.
Definition: NURBS3DVolume.C:1481
Foam::NURBS3DVolume::dndbBasedSensitivities
tmp< tensorField > dndbBasedSensitivities(const label patchI, const label cpI, bool DimensionedNormalSens=true)
Definition: NURBS3DVolume.C:1152
Foam::GeometricField::Boundary
The boundary fields.
Definition: GeometricField.H:115
DebugInfo
#define DebugInfo
Report an information message using Foam::Info.
Definition: messageStream.H:382
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::BitOps::count
unsigned int count(const UList< bool > &bools, const bool val=true)
Count number of 'true' entries.
Definition: BitOps.H:77
Foam::GeometricField::boundaryFieldRef
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
Definition: GeometricField.C:783
Foam::Vector< scalar >
Foam::NURBS3DVolume::getReverseMap
const labelList & getReverseMap()
Definition: NURBS3DVolume.C:1593
Foam::NURBS3DVolume::boundControlPointMovement
void boundControlPointMovement(vectorField &controlPointsMovement)
Definition: NURBS3DVolume.C:1495
Foam::List< label >
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::NURBS3DVolume::JacobianUVW
tensor JacobianUVW(const vector &u) const
Jacobian matrix wrt to the volume parametric coordinates.
Definition: NURBS3DVolume.C:893
points
const pointField & points
Definition: gmvOutputHeader.H:1
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
Foam::face::points
pointField points(const UList< point > &pts) const
Return the points corresponding to this face.
Definition: faceI.H:87
Foam::NURBS3DVolume::coordinates
tmp< vectorField > coordinates(const vectorField &uVector) const
Definition: NURBS3DVolume.C:1318
Foam::tmp::New
static tmp< T > New(Args &&... args)
Construct tmp of T with forwarding arguments.
NURBS3DVolume.H
Foam::NURBS3DVolume::NURBS3DVolume
NURBS3DVolume(const dictionary &dict, const fvMesh &mesh, bool computeParamCoors=true)
Construct from dictionary.
Definition: NURBS3DVolume.C:628
Foam::boolVector::x
bool x() const
The x component.
Definition: boolVectorI.H:110
Foam::NURBS3DVolume::computeMaxBoundaryDisplacement
scalar computeMaxBoundaryDisplacement(const vectorField &controlPointsMovement, const labelList &patchesToBeMoved)
Compute max. displacement at the boundary.
Definition: NURBS3DVolume.C:1518
Foam::NURBS3DVolume::patchDxDbFace
tmp< tensorField > patchDxDbFace(const label patchI, const label cpI)
Get patch dx/db.
Definition: NURBS3DVolume.C:1260
Foam::NURBS3DVolume::writeCps
void writeCps(const fileName &="cpsFile", const bool transform=true) const
Definition: NURBS3DVolume.C:1780
minValue
scalar minValue
Definition: LISASMDCalcMethod2.H:12
Foam::NURBS3DVolume::New
static autoPtr< NURBS3DVolume > New(const dictionary &dict, const fvMesh &mesh, bool computeParamCoors=true)
Return a reference to the selected NURBS model.
Definition: NURBS3DVolume.C:750
Foam::boolVector::y
bool y() const
The y component.
Definition: boolVectorI.H:111
Foam::NURBS3DVolume::getDxDb
tmp< pointTensorField > getDxDb(const label cpI)
Get dxCartesiandb for a certain control point.
Definition: NURBS3DVolume.C:1622
Foam::VectorSpace< Vector< scalar >, scalar, 3 >::zero
static const Vector< scalar > zero
Definition: VectorSpace.H:115
Foam::dictionary::add
entry * add(entry *entryPtr, bool mergeEntry=false)
Add a new entry.
Definition: dictionary.C:640
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:72
Foam::pointCells
Smooth ATC in cells having a point to a set of patches supplied by type.
Definition: pointCells.H:56
Foam::plusEqOp
Definition: ops.H:72
Foam::NURBS3DVolume::getDxCellsDb
tmp< volTensorField > getDxCellsDb(const label cpI)
Get dxCartesiandb for a certain control point on cells.
Definition: NURBS3DVolume.C:1666
Foam::Pstream::listCombineScatter
static void listCombineScatter(const List< commsStruct > &comms, List< T > &Value, const int tag, const label comm)
Scatter data. Reverse of combineGather.
Definition: combineGatherScatter.C:432
Foam::point
vector point
Point is a vector.
Definition: point.H:43
Foam::dictionary::getOrDefault
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:148
Foam::GeometricField< vector, pointPatchField, pointMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:188
Foam::mkDir
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:507
Foam::NURBS3DVolume::patchDxDb
tmp< tensorField > patchDxDb(const label patchI, const label cpI)
Get patch dx/db.
Definition: NURBS3DVolume.C:1221
Foam::NURBS3DVolume::computeNewPoints
tmp< vectorField > computeNewPoints(const vectorField &controlPointsMovement)
Mesh movement based on given control point movement.
Definition: NURBS3DVolume.C:1376
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:328
Foam::dimensionedTensor
dimensioned< tensor > dimensionedTensor
Dimensioned tensor obtained from generic dimensioned type.
Definition: dimensionedTensor.H:52
Foam::gMax
Type gMax(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:592
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Definition: GeometricFieldI.H:62
Foam::dimless
const dimensionSet dimless
Dimensionless.
Definition: dimensionSets.C:189
maxValue
scalar maxValue
Definition: LISASMDCalcMethod1.H:5
Foam::dictionary::dictName
word dictName() const
The local dictionary name (final part of scoped name)
Definition: dictionaryI.H:60