NURBS3DVolume.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2007-2023 PCOpt/NTUA
9  Copyright (C) 2013-2023 FOSS GP
10  Copyright (C) 2019-2020 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 Class
29  Foam::NURBS3DVolume
30 
31 Description
32  NURBS3DVolume morpher. Includes support functions for gradient computations
33  Base class providing support for different coordinate systems
34 
35  Reference:
36  \verbatim
37  For a short introduction to a volumetric B-Splines morpher and its use
38  in shape optimisation
39 
40  Papoutsis-Kiachagias, E., Magoulas, N., Mueller, J.,
41  Othmer, C., & Giannakoglou, K. (2015).
42  Noise reduction in car aerodynamics using a surrogate objective
43  function and the continuous adjoint method with wall functions.
44  Computers & Fluids, 122, 223-232.
45  http://doi.org/10.1016/j.compfluid.2015.09.002
46  \endverbatim
47 
48 SourceFiles
49  NURBS3DVolume.C
50  NURBS3DVolumeI.H
51 
52 \*---------------------------------------------------------------------------*/
53 
54 #ifndef NURBS3DVolume_H
55 #define NURBS3DVolume_H
56 
57 #include "NURBSbasis.H"
58 #include "boolVector.H"
59 #include "localIOdictionary.H"
60 #include "pointMesh.H"
61 #include "pointPatchField.H"
62 #include "pointPatchFieldsFwd.H"
63 #include "boundaryFieldsFwd.H"
64 
65 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
66 
67 namespace Foam
68 {
69 
70 /*---------------------------------------------------------------------------*\
71  Class NURBS3DVolume Declaration
72 \*---------------------------------------------------------------------------*/
73 
74 class NURBS3DVolume
75 :
76  public localIOdictionary
77 {
78 protected:
79 
80  // Protected Data
81 
83 
84  const fvMesh& mesh_;
87  //- NURBS basis functions
91 
92  //- Max iterations for Newton-Raphson
93  label maxIter_;
94 
95  //- Tolerance for Newton-Raphson
96  scalar tolerance_;
97 
98  //- How many times to bound parametric coordinates until deciding it is
99  //- outside the box
100  label nMaxBound_;
101 
102  //- The volumetric B-Splines control points
104 
105  //- Map of points-in-box to mesh points
107 
108  //- Map of mesh points to points-in-box
109  // Return -1 if point not within the box
111 
112  //- Parametric coordinates of pointsInBox
114 
115  //- Coordinates in the local system for which CPs are defined
117 
118  //- Confine movement in certain directions and control points. Refers
119  //- to the local system
121 
122  label confineVMovement_;
123 
124  label confineWMovement_;
127 
128  //- Which movement components to freeze in each plane
132 
134 
141  //- Which of the cps are moved in an optimisation
143 
144  //- Which design variables are changed in an optimisation
146 
147  //- Folder to store control points
148  string cpsFolder_;
150  //- Read parametric coordinates from file if present in the folder
152 
154  // Protected Member Functions
156  //- Find points within control points box
157  void findPointsInBox(const vectorField& meshPoints);
158 
159  //- Compute parametric coordinates in order to match a given set
160  //- of coordinates, based on the cps of the class
161  // Uses a Newton-Raphson loop.
162  // Argument is the points residing in the box
164 
166 
167  //- Bound components to certain limits
168  bool bound
169  (
170  vector& vec,
171  scalar minValue = 1e-7,
172  scalar maxValue = 0.999999
173  ) const;
174 
175  //- Create lists with active design variables and control points
178  //- Confine movement in boundary control points if necessary
180 
181  //- Confine control point movement to maintain user-defined continuity
183 
184  //- Confine movement in all control points for user-defined directions
186 
187  //- Deactivate control points if they affect no mesh point
189 
190  //- Confine all three movements for a prescribed control point
191  void confineControlPoint(const label cpI);
192 
193  //- Confine specific movements for a prescribed control point
194  void confineControlPoint(const label cpI, const boolVector&);
195 
196  //- Create folders to store cps and derivatives
197  void makeFolders();
198 
199  //- Transform a point from its coordinate system to a cartesian system
201  (
202  const vector& localCoordinates
203  ) const = 0;
204 
205  //- Transformation tensor for dxdb, from local coordinate system to
206  //- cartesian
207  virtual tensor transformationTensorDxDb(label globalPointIndex) = 0;
208 
209  //- Update coordinates in the local system based on the cartesian points
210  virtual void updateLocalCoordinateSystem
211  (
212  const vectorField& cartesianPoints
213  ) = 0;
214 
215 
216 public:
217 
218  //- Runtime type information
219  TypeName("NURBS3DVolume");
220 
221 
222  // Declare run-time constructor selection table
223 
225  (
226  autoPtr,
228  dictionary,
229  (
230  const dictionary& dict,
231  const fvMesh& mesh,
232  bool computeParamCoors
233  ),
234  (dict, mesh, computeParamCoors)
235  );
236 
237  // Constructors
238 
239  //- Construct from dictionary
241  (
242  const dictionary& dict,
243  const fvMesh& mesh,
244  bool computeParamCoors = true
245  );
246 
247  //- Construct as copy
249 
250 
251  // Selectors
252 
253  //- Return a reference to the selected NURBS model
255  (
256  const dictionary& dict,
257  const fvMesh& mesh,
258  bool computeParamCoors = true
259  );
260 
261 
262  //- Destructor
263  virtual ~NURBS3DVolume() = default;
264 
265 
266  // Member Functions
267 
268  //- Compute parametric coordinates for a given set of points
269  // (not necessarily the mesh ones)
271  (
272  const vectorField& points
273  ) const;
274 
275  // Derivatives wrt parametric coordinates
276 
277  //- Volume point derivative wrt u at point u,v,w
279  (
280  const scalar u,
281  const scalar v,
282  const scalar w
283  ) const;
284 
285  //- Volume point derivative wrt v at point u,v,w
287  (
288  const scalar u,
289  const scalar v,
290  const scalar w
291  ) const;
292 
293  //- Volume point derivative wrt w at point u,v,w
295  (
296  const scalar u,
297  const scalar v,
298  const scalar w
299  ) const;
300 
301  //- Jacobian matrix wrt to the volume parametric coordinates
302  tensor JacobianUVW(const vector& u) const;
303 
304 
305  // Derivatives wrt to control points
306 
307  //- Volume point derivative wrt to control point cpI at point u,v,w
308  // Scalar since in the local system!
309  scalar volumeDerivativeCP(const vector& u, const label cpI) const;
310 
311  //- Control point sensitivities computed using point-based surface
312  //- sensitivities
314  (
315  const pointVectorField& pointSens,
316  const labelList& sensitivityPatchIDs
317  );
318 
319  //- Control point sensitivities computed using face-based surface
320  //- sensitivities
322  (
323  const volVectorField& faceSens,
324  const labelList& sensitivityPatchIDs
325  );
326 
327  //- Control point sensitivities computed using face-based surface
328  //- sensitivities
330  (
331  const boundaryVectorField& faceSens,
332  const labelList& sensitivityPatchIDs
333  );
334 
335  //- Control point sensitivities computed using face-based surface
336  //- sensitivities
338  (
339  const vectorField& faceSens,
340  const label patchI,
341  const label cpI
342  );
343 
344  //- Part of control point sensitivities related to the face normal
345  //- variations
347  (
348  const label patchI,
349  const label cpI,
350  bool DimensionedNormalSens = true
351  );
352 
353  //- Get patch dx/db
355  (
356  const label patchI,
357  const label cpI
358  );
359 
360  //- Get patch dx/db
362  (
363  const label patchI,
364  const label cpI
365  );
366 
367 
368  // Cartesian coordinates
369 
370  //- Compute cartesian coordinates based on control points
371  //- and parametric coordinates
372  tmp<vectorField> coordinates(const vectorField& uVector) const;
373 
374  //- The same, for a specific point
375  vector coordinates(const vector& uVector) const;
376 
377  //- Mesh movement based on given control point movement
379  (
380  const vectorField& controlPointsMovement
381  );
382 
383  //- Boundary mesh movement based on given control point movement
385  (
386  const vectorField& controlPointsMovement,
387  const labelList& patchesToBeMoved,
388  const bool moveCPs = true
389  );
390 
391 
392  // Control points management
393 
394  //- Get control point ID from its I-J-K coordinates
395  label getCPID(const label i, const label j, const label k) const;
396 
397  //- Get I-J-K coordinates from control point ID
398  void getIJK(label& i, label& j, label& k, const label cpID) const;
399 
400  //- Set new control points
401  // New values should be on the coordinates system original CPs
402  // were defined
403  void setControlPoints(const vectorField& newCps);
404 
405 
406  //- Bound control points movement in the boundary control points
407  //- and in certain directions if needed
409  (
410  vectorField& controlPointsMovement
411  ) const;
412 
413  //- Compute max. displacement at the boundary
415  (
416  const vectorField& controlPointsMovement,
417  const labelList& patchesToBeMoved
418  );
419 
420 
421  // Access Functions
422 
423  // Functions triggering calculations
424 
425  //- Get mesh points that reside within the control points box
427 
428  //- Get map of points in box to mesh points
429  const labelList& getMap();
430 
431  //- Get map of mesh points to points in box.
432  //- Return -1 if point is outside the box
433  const labelList& getReverseMap();
434 
435  //- Get parametric coordinates
437 
438  //- Get dxCartesiandb for a certain control point
439  tmp<pointTensorField> getDxDb(const label cpI);
440 
441  //- Get dxCartesiandb for a certain control point on cells
442  tmp<volTensorField> getDxCellsDb(const label cpI);
443 
444  //- Get number of variables if CPs are moved symmetrically in U
445  label nUSymmetry() const;
446 
447  //- Get number of variables if CPs are moved symmetrically in V
448  label nVSymmetry() const;
449 
450  //- Get number of variables if CPs are moved symmetrically in W
451  label nWSymmetry() const;
452 
453  //- Get number of variables per direction,
454  //- if CPs are moved symmetrically
455  Vector<label> nSymmetry() const;
456 
457 
458  // Inline access functions
459 
460  //- Get box name
461  inline const word& name() const;
462 
463  //- Which control points are active?
464  // A control point is active if at least one component can move
465  inline const boolList& getActiveCPs() const;
466 
467  //- Which design variables are active?
468  // Numbering is (X1,Y1,Z1), (X2,Y2,Z2) ...
469  inline const boolList& getActiveDesignVariables() const;
470 
471  //- Get control points
472  inline const vectorField& getControlPoints() const;
473 
474  inline vectorField& getControlPoints();
475 
476  //- Get confine movements
477  inline bool confineUMovement() const;
478 
479  inline bool confineVMovement() const;
480 
481  inline bool confineWMovement() const;
482 
483  //- Get basis functions
484  inline const NURBSbasis& basisU() const;
485 
486  inline const NURBSbasis& basisV() const;
487 
488  inline const NURBSbasis& basisW() const;
489 
490  //- Get number of control points per direction
491  inline Vector<label> nCPsPerDirection() const;
492 
493  //- Get mesh
494  inline const fvMesh& mesh() const;
495 
496  //- Get dictionary
497  inline const dictionary& dict() const;
498 
499 
500  // Write Functions
501 
502  //- Write control points on a cartesian coordinates system for
503  //- visualization
504  void writeCps
505  (
506  const fileName& = "cpsFile",
507  const bool transform = true
508  ) const;
509 
510  //- Write parametric coordinates
511  void writeParamCoordinates() const;
512 
513  //- Write the control points to support restart
514  virtual bool writeData(Ostream& os) const;
515 };
516 
517 
518 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
519 
520 } // End namespace Foam
521 
522 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
523 
524 #include "NURBS3DVolumeI.H"
525 
526 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
527 
528 #endif
529 
530 // ************************************************************************* //
void continuityRealatedConfinement()
Confine control point movement to maintain user-defined continuity.
label maxIter_
Max iterations for Newton-Raphson.
Definition: NURBS3DVolume.H:92
vector volumeDerivativeU(const scalar u, const scalar v, const scalar w) const
Volume point derivative wrt u at point u,v,w.
label nWSymmetry() const
Get number of variables if CPs are moved symmetrically in W.
declareRunTimeSelectionTable(autoPtr, NURBS3DVolume, dictionary,(const dictionary &dict, const fvMesh &mesh, bool computeParamCoors),(dict, mesh, computeParamCoors))
tmp< vectorField > getPointsInBox()
Get mesh points that reside within the control points box.
A class for handling file names.
Definition: fileName.H:72
vectorField localSystemCoordinates_
Coordinates in the local system for which CPs are defined.
string cpsFolder_
Folder to store control points.
autoPtr< labelList > reverseMapPtr_
Map of mesh points to points-in-box.
tmp< tensorField > patchDxDb(const label patchI, const label cpI)
Get patch dx/db.
boolVectorList confineWMinCPs_
bool confineWMovement() const
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
boolVectorList confineWMaxCPs_
void confineControlPoint(const label cpI)
Confine all three movements for a prescribed control point.
scalar minValue
boolVectorList confineVMinCPs_
const fvMesh & mesh() const
Get mesh.
vectorField cps_
The volumetric B-Splines control points.
const NURBSbasis & basisW() const
bool bound(vector &vec, scalar minValue=1e-7, scalar maxValue=0.999999) const
Bound components to certain limits.
void getIJK(label &i, label &j, label &k, const label cpID) const
Get I-J-K coordinates from control point ID.
tmp< tensorField > dndbBasedSensitivities(const label patchI, const label cpI, bool DimensionedNormalSens=true)
Part of control point sensitivities related to the face normal variations.
autoPtr< labelList > mapPtr_
Map of points-in-box to mesh points.
const NURBSbasis & basisU() const
Get basis functions.
tmp< tensorField > patchDxDbFace(const label patchI, const label cpI)
Get patch dx/db.
bool readStoredData_
Read parametric coordinates from file if present in the folder.
label k
Boltzmann constant.
label nUSymmetry() const
Get number of variables if CPs are moved symmetrically in U.
const pointVectorField & getParametricCoordinates()
Get parametric coordinates.
NURBS3DVolume morpher. Includes support functions for gradient computations Base class providing supp...
Definition: NURBS3DVolume.H:69
const boolList & getActiveDesignVariables() const
Which design variables are active?
static autoPtr< NURBS3DVolume > New(const dictionary &dict, const fvMesh &mesh, bool computeParamCoors=true)
Return a reference to the selected NURBS model.
virtual tensor transformationTensorDxDb(label globalPointIndex)=0
Transformation tensor for dxdb, from local coordinate system to cartesian.
Vector< label > nCPsPerDirection() const
Get number of control points per direction.
void confineInertControlPoints()
Deactivate control points if they affect no mesh point.
localIOdictionary is derived from IOdictionary but excludes parallel master reading.
vectorField computeControlPointSensitivities(const pointVectorField &pointSens, const labelList &sensitivityPatchIDs)
Control point sensitivities computed using point-based surface sensitivities.
void determineActiveDesignVariablesAndPoints()
Create lists with active design variables and control points.
void makeFolders()
Create folders to store cps and derivatives.
scalar maxValue
const fvMesh & mesh_
Definition: NURBS3DVolume.H:79
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
tensor JacobianUVW(const vector &u) const
Jacobian matrix wrt to the volume parametric coordinates.
const pointField & points
void writeParamCoordinates() const
Write parametric coordinates.
A class for handling words, derived from Foam::string.
Definition: word.H:63
boolVectorList confineUMinCPs_
Which movement components to freeze in each plane.
virtual bool writeData(Ostream &os) const
Write the control points to support restart.
void setControlPoints(const vectorField &newCps)
Set new control points.
tmp< vectorField > computeNewBoundaryPoints(const vectorField &controlPointsMovement, const labelList &patchesToBeMoved, const bool moveCPs=true)
Boundary mesh movement based on given control point movement.
void confineControlPointsDirections()
Confine movement in all control points for user-defined directions.
NURBSbasis function. Used to construct NURBS curves, surfaces and volumes.
Definition: NURBSbasis.H:51
vector volumeDerivativeV(const scalar u, const scalar v, const scalar w) const
Volume point derivative wrt v at point u,v,w.
scalar tolerance_
Tolerance for Newton-Raphson.
Definition: NURBS3DVolume.H:97
void findPointsInBox(const vectorField &meshPoints)
Find points within control points box.
Definition: NURBS3DVolume.C:43
label getCPID(const label i, const label j, const label k) const
Get control point ID from its I-J-K coordinates.
void boundControlPointMovement(vectorField &controlPointsMovement) const
Bound control points movement in the boundary control points and in certain directions if needed...
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
label nVSymmetry() const
Get number of variables if CPs are moved symmetrically in V.
void writeCps(const fileName &="cpsFile", const bool transform=true) const
Write control points on a cartesian coordinates system for visualization.
Vector< label > nSymmetry() const
Get number of variables per direction, if CPs are moved symmetrically.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
Specialized bundling of boolean values as a vector of 3 components, element access using x()...
Definition: boolVector.H:52
OBJstream os(runTime.globalPath()/outputName)
tmp< vectorField > coordinates(const vectorField &uVector) const
Compute cartesian coordinates based on control points and parametric coordinates. ...
boolVectorList confineVMaxCPs_
void computeParametricCoordinates(const vectorField &points)
Compute parametric coordinates in order to match a given set of coordinates, based on the cps of the ...
Definition: NURBS3DVolume.C:98
label confineBoundaryControlPoints_
virtual void updateLocalCoordinateSystem(const vectorField &cartesianPoints)=0
Update coordinates in the local system based on the cartesian points.
label nMaxBound_
How many times to bound parametric coordinates until deciding it is outside the box.
tmp< vectorField > computeNewPoints(const vectorField &controlPointsMovement)
Mesh movement based on given control point movement.
virtual ~NURBS3DVolume()=default
Destructor.
Useful typenames for fields defined only at the boundaries.
const labelList & getMap()
Get map of points in box to mesh points.
void confineBoundaryControlPoints()
Confine movement in boundary control points if necessary.
TypeName("NURBS3DVolume")
Runtime type information.
autoPtr< pointVectorField > parametricCoordinatesPtr_
Parametric coordinates of pointsInBox.
const NURBSbasis & basisV() const
boolList activeDesignVariables_
Which design variables are changed in an optimisation.
vector volumeDerivativeW(const scalar u, const scalar v, const scalar w) const
Volume point derivative wrt w at point u,v,w.
List< boolVector > boolVectorList
Definition: NURBS3DVolume.H:77
NURBS3DVolume(const dictionary &dict, const fvMesh &mesh, bool computeParamCoors=true)
Construct from dictionary.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
const word & name() const
Get box name.
scalar computeMaxBoundaryDisplacement(const vectorField &controlPointsMovement, const labelList &patchesToBeMoved)
Compute max. displacement at the boundary.
bool confineUMovement() const
Get confine movements.
const boolList & getActiveCPs() const
Which control points are active?
const labelList & getReverseMap()
Get map of mesh points to points in box. Return -1 if point is outside the box.
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
bool confineVMovement() const
label confineUMovement_
Confine movement in certain directions and control points. Refers to the local system.
tmp< pointTensorField > getDxDb(const label cpI)
Get dxCartesiandb for a certain control point.
tmp< volTensorField > getDxCellsDb(const label cpI)
Get dxCartesiandb for a certain control point on cells.
A class for managing temporary objects.
Definition: HashPtrTable.H:50
boolVectorList confineUMaxCPs_
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:521
Tensor of scalars, i.e. Tensor<scalar>.
NURBSbasis basisU_
NURBS basis functions.
Definition: NURBS3DVolume.H:85
const vectorField & getControlPoints() const
Get control points.
const dictionary & dict() const
Get dictionary.
virtual vector transformPointToCartesian(const vector &localCoordinates) const =0
Transform a point from its coordinate system to a cartesian system.
boolList activeControlPoints_
Which of the cps are moved in an optimisation.
Namespace for OpenFOAM.
scalar volumeDerivativeCP(const vector &u, const label cpI) const
Volume point derivative wrt to control point cpI at point u,v,w.