Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
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
28Class
29 Foam::NURBS3DVolume
30
31Description
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
48SourceFiles
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
67namespace Foam
68{
70/*---------------------------------------------------------------------------*\
71 Class NURBS3DVolume Declaration
72\*---------------------------------------------------------------------------*/
73
74class NURBS3DVolume
75:
78protected:
80 // Protected Data
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_;
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
216public:
217
218 //- Runtime type information
219 TypeName("NURBS3DVolume");
220
221
222 // Declare run-time constructor selection table
223
225 (
226 autoPtr,
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 );
261
262 //- Destructor
263 virtual ~NURBS3DVolume() = default;
264
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
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
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// ************************************************************************* //
scalar maxValue
scalar minValue
label k
Useful typenames for fields defined only at the boundaries.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
bool readStoredData_
Read parametric coordinates from file if present in the folder.
void getIJK(label &i, label &j, label &k, const label cpID) const
Get I-J-K coordinates from control point ID.
label maxIter_
Max iterations for Newton-Raphson.
virtual void updateLocalCoordinateSystem(const vectorField &cartesianPoints)=0
Update coordinates in the local system based on the cartesian points.
boolVectorList confineWMinCPs_
static autoPtr< NURBS3DVolume > New(const dictionary &dict, const fvMesh &mesh, bool computeParamCoors=true)
Return a reference to the selected NURBS model.
boolVectorList confineVMinCPs_
tmp< tensorField > patchDxDbFace(const label patchI, const label cpI)
Get patch dx/db.
void computeParametricCoordinates(const vectorField &points)
Compute parametric coordinates in order to match a given set of coordinates, based on the cps of the ...
label nUSymmetry() const
Get number of variables if CPs are moved symmetrically in U.
void confineInertControlPoints()
Deactivate control points if they affect no mesh point.
const fvMesh & mesh_
autoPtr< labelList > mapPtr_
Map of points-in-box to mesh points.
declareRunTimeSelectionTable(autoPtr, NURBS3DVolume, dictionary,(const dictionary &dict, const fvMesh &mesh, bool computeParamCoors),(dict, mesh, computeParamCoors))
vector volumeDerivativeU(const scalar u, const scalar v, const scalar w) const
Volume point derivative wrt u at point u,v,w.
void confineBoundaryControlPoints()
Confine movement in boundary control points if necessary.
const labelList & getMap()
Get map of points in box to mesh points.
scalar tolerance_
Tolerance for Newton-Raphson.
virtual bool writeData(Ostream &os) const
Write the control points to support restart.
scalar computeMaxBoundaryDisplacement(const vectorField &controlPointsMovement, const labelList &patchesToBeMoved)
Compute max. displacement at the boundary.
boolVectorList confineVMaxCPs_
void findPointsInBox(const vectorField &meshPoints)
Find points within control points box.
tmp< vectorField > computeNewBoundaryPoints(const vectorField &controlPointsMovement, const labelList &patchesToBeMoved, const bool moveCPs=true)
Boundary mesh movement based on given control point movement.
boolVectorList confineWMaxCPs_
const boolList & getActiveDesignVariables() const
Which design variables are active?
vectorField localSystemCoordinates_
Coordinates in the local system for which CPs are defined.
boolVectorList confineUMaxCPs_
vectorField computeControlPointSensitivities(const pointVectorField &pointSens, const labelList &sensitivityPatchIDs)
Control point sensitivities computed using point-based surface sensitivities.
boolVectorList confineUMinCPs_
Which movement components to freeze in each plane.
const pointVectorField & getParametricCoordinates()
Get parametric coordinates.
tmp< volTensorField > getDxCellsDb(const label cpI)
Get dxCartesiandb for a certain control point on cells.
virtual tensor transformationTensorDxDb(label globalPointIndex)=0
Transformation tensor for dxdb, from local coordinate system to cartesian.
autoPtr< labelList > reverseMapPtr_
Map of mesh points to points-in-box.
label nVSymmetry() const
Get number of variables if CPs are moved symmetrically in V.
tmp< pointTensorField > getDxDb(const label cpI)
Get dxCartesiandb for a certain control point.
NURBS3DVolume(const NURBS3DVolume &)
Construct as copy.
virtual ~NURBS3DVolume()=default
Destructor.
bool confineVMovement() const
bool bound(vector &vec, scalar minValue=1e-7, scalar maxValue=0.999999) const
Bound components to certain limits.
void makeFolders()
Create folders to store cps and derivatives.
void writeCps(const fileName &="cpsFile", const bool transform=true) const
Write control points on a cartesian coordinates system for visualization.
NURBS3DVolume(const dictionary &dict, const fvMesh &mesh, bool computeParamCoors=true)
Construct from dictionary.
boolList activeDesignVariables_
Which design variables are changed in an optimisation.
tmp< tensorField > patchDxDb(const label patchI, const label cpI)
Get patch dx/db.
tmp< vectorField > getPointsInBox()
Get mesh points that reside within the control points box.
vector volumeDerivativeV(const scalar u, const scalar v, const scalar w) const
Volume point derivative wrt v at point u,v,w.
virtual vector transformPointToCartesian(const vector &localCoordinates) const =0
Transform a point from its coordinate system to a cartesian system.
tmp< tensorField > dndbBasedSensitivities(const label patchI, const label cpI, bool DimensionedNormalSens=true)
Part of control point sensitivities related to the face normal variations.
string cpsFolder_
Folder to store control points.
NURBSbasis basisU_
NURBS basis functions.
label confineUMovement_
Confine movement in certain directions and control points. Refers to the local system.
void confineControlPoint(const label cpI)
Confine all three movements for a prescribed control point.
bool confineWMovement() const
const dictionary & dict() const
Get dictionary.
Vector< label > nSymmetry() const
Get number of variables per direction, if CPs are moved symmetrically.
vectorField cps_
The volumetric B-Splines control points.
const NURBSbasis & basisV() const
void boundControlPointMovement(vectorField &controlPointsMovement) const
Bound control points movement in the boundary control points and in certain directions if needed.
void continuityRealatedConfinement()
Confine control point movement to maintain user-defined continuity.
label getCPID(const label i, const label j, const label k) const
Get control point ID from its I-J-K coordinates.
autoPtr< pointVectorField > parametricCoordinatesPtr_
Parametric coordinates of pointsInBox.
scalar volumeDerivativeCP(const vector &u, const label cpI) const
Volume point derivative wrt to control point cpI at point u,v,w.
label confineBoundaryControlPoints_
label nMaxBound_
How many times to bound parametric coordinates until deciding it is outside the box.
const NURBSbasis & basisW() const
const fvMesh & mesh() const
Get mesh.
label nWSymmetry() const
Get number of variables if CPs are moved symmetrically in W.
void writeParamCoordinates() const
Write parametric coordinates.
void determineActiveDesignVariablesAndPoints()
Create lists with active design variables and control points.
Vector< label > nCPsPerDirection() const
Get number of control points per direction.
tmp< vectorField > coordinates(const vectorField &uVector) const
Compute cartesian coordinates based on control points and parametric coordinates.
tmp< vectorField > computeNewPoints(const vectorField &controlPointsMovement)
Mesh movement based on given control point movement.
const word & name() const
Get box name.
List< boolVector > boolVectorList
bool confineUMovement() const
Get confine movements.
const boolList & getActiveCPs() const
Which control points are active?
tensor JacobianUVW(const vector &u) const
Jacobian matrix wrt to the volume parametric coordinates.
void setControlPoints(const vectorField &newCps)
Set new control points.
boolList activeControlPoints_
Which of the cps are moved in an optimisation.
const labelList & getReverseMap()
Get map of mesh points to points in box. Return -1 if point is outside the box.
TypeName("NURBS3DVolume")
Runtime type information.
const NURBSbasis & basisU() const
Get basis functions.
const vectorField & getControlPoints() const
Get control points.
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.
NURBSbasis function. Used to construct NURBS curves, surfaces and volumes.
Definition NURBSbasis.H:52
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
Templated 3D Vector derived from VectorSpace adding construction from 3 components,...
Definition Vector.H:61
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
Specialized bundling of boolean values as a vector of 3 components, element access using x(),...
Definition boolVector.H:55
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
A class for handling file names.
Definition fileName.H:75
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
localIOdictionary(const IOobject &io, const dictionary *fallback=nullptr)
Construct given an IOobject and optional fallback dictionary content.
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
OBJstream os(runTime.globalPath()/outputName)
const pointField & points
Namespace for OpenFOAM.
GeometricField< vector, fvPatchField, volMesh > volVectorField
List< label > labelList
A List of labels.
Definition List.H:62
refinementData transform(const tensor &, const refinementData val)
No-op rotational transform for base types.
GeometricField< vector, pointPatchField, pointMesh > pointVectorField
Tensor< scalar > tensor
Definition symmTensor.H:57
Field< vector > vectorField
Specialisation of Field<T> for vector.
List< bool > boolList
A List of bools.
Definition List.H:60
Vector< scalar > vector
Definition vector.H:57
volVectorField::Boundary boundaryVectorField
#define declareRunTimeSelectionTable(ptrWrapper, baseType, argNames, argList, parList)
Declare a run-time selection (variables and adder classes).
volScalarField & e
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition typeInfo.H:68