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-2019 PCOpt/NTUA
9  Copyright (C) 2013-2019 FOSS GP
10  Copyright (C) 2019 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 "pointMesh.H"
59 #include "pointPatchField.H"
60 #include "pointPatchFieldsFwd.H"
61 #include "boundaryFieldsFwd.H"
62 
63 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
64 
65 namespace Foam
66 {
67 
68 /*---------------------------------------------------------------------------*\
69  Class NURBS3DVolume Declaration
70 \*---------------------------------------------------------------------------*/
71 
72 class NURBS3DVolume
73 {
74 protected:
75 
76  // Protected data
77 
79 
80  //- NURBS basis functions
81  const fvMesh& mesh_;
86 
87  //- Max iterations for Newton-Raphson
89 
90  //- Tolerance for Newton-Raphson
91  scalar tolerance_;
92 
93  //- How many times to bound parametric coordinates until deciding it is
94  //- outside the box
96 
97  //- The volumetric B-Splines control points
99 
100  //- Map of points-in-box to mesh points
102 
103  //- Map of mesh points to points-in-box
104  // Return -1 if point not whithin the box
106 
107  //- Parametric coordinates of pointsInBox
109 
110  //- Coordinates in the local system for which CPs are defined
112 
113  //- Confine movement in certain directions and control points. Refers
114  //- to the local system
116 
118 
120 
122 
123  //- Which movement components to freeze in each plane
125 
127 
129 
131 
133 
135 
136  //- Which of the cps are moved in an optimisation
138 
139  //- Which design variables are changed in an optimisation
141 
142  //- Folder to store control points
143  string cpsFolder_;
144 
145  //- Read parametric coordinates from file if present in the folder
146  bool readStoredData_;
147 
148 
149  // Protected Member Functions
150 
151  //- Get control point ID from its I-J-K coordinates
152  label getCPID(const label i, const label j, const label k) const;
153 
154  //- Find points within control points box
155  void findPointsInBox(const vectorField& meshPoints);
156 
157  //- Compute parametric coordinates in order to match a given set
158  //- of coordinates, based on the cps of the class
159  // Uses a Newton-Raphson loop.
160  // Argument is the points residing in the box
162 
164 
165  //- Bound components to certain limits
166  bool bound
167  (
168  vector& vec,
169  scalar minValue = 1e-7,
170  scalar maxValue = 0.999999
171  );
172 
173  //- Create lists with active design variables and control points
175 
176  //- Confine movement in boundary control points if necessary
178 
179  //- Confine control point movement to maintain user-defined continuity
181 
182  //- Confine movement in all control points for user-defined directions
184 
185  //- Confine all three movements for a prescribed control point
186  void confineControlPoint(const label cpI);
187 
188  //- Confine specific movements for a prescribed control point
189  void confineControlPoint(const label cpI, const FixedList<bool, 3>&);
190 
191  //- Create folders to store cps and derivatives
192  void makeFolders();
193 
194  //- Transform a point from its coordinate system to a cartesian system
196  (
197  const vector& localCoordinates
198  ) const = 0;
199 
200  //- Transformation tensor for dxdb, from local coordinate system to
201  //- cartesian
202  virtual tensor transformationTensorDxDb(label globalPointIndex) = 0;
203 
204  //- Update coordinates in the local system based on the cartesian points
205  virtual void updateLocalCoordinateSystem
206  (
207  const vectorField& cartesianPoints
208  ) = 0;
209 
210 
211 public:
212 
213  //- Runtime type information
214  TypeName("NURBS3DVolume");
215 
216 
217  // Declare run-time constructor selection table
218 
220  (
221  autoPtr,
223  dictionary,
224  (
225  const dictionary& dict,
226  const fvMesh& mesh,
227  bool computeParamCoors
228  ),
229  (dict, mesh, computeParamCoors)
230  );
231 
232  // Constructors
233 
234  //- Construct from dictionary
236  (
237  const dictionary& dict,
238  const fvMesh& mesh,
239  bool computeParamCoors = true
240  );
241 
242  //- Construct as copy
244 
245 
246  // Selectors
247 
248  //- Return a reference to the selected NURBS model
250  (
251  const dictionary& dict,
252  const fvMesh& mesh,
253  bool computeParamCoors = true
254  );
255 
256 
257  //- Destructor
258  virtual ~NURBS3DVolume() = default;
259 
260 
261  // Member Functions
262 
263  // Derivatives wrt parametric coordinates
264 
265  //- Volume point derivative wrt u at point u,v,w
267  (
268  const scalar u,
269  const scalar v,
270  const scalar w
271  ) const;
272 
273  //- Volume point derivative wrt v at point u,v,w
275  (
276  const scalar u,
277  const scalar v,
278  const scalar w
279  ) const;
280 
281  //- Volume point derivative wrt w at point u,v,w
283  (
284  const scalar u,
285  const scalar v,
286  const scalar w
287  ) const;
288 
289  //- Jacobian matrix wrt to the volume parametric coordinates
290  tensor JacobianUVW(const vector& u) const;
291 
292 
293  // Derivatives wrt to control points
294 
295  //- Volume point derivative wrt to control point cpI at point u,v,w
296  // Scalar since in the local system!
297  scalar volumeDerivativeCP(const vector& u, const label cpI) const;
298 
299  //- Control point sensitivities computed using point-based surface
300  //- sensitivities
302  (
303  const pointVectorField& pointSens,
304  const labelList& sensitivityPatchIDs
305  );
306 
307  //- Control point sensitivities computed using face-based surface
308  //- sensitivities
310  (
311  const volVectorField& faceSens,
312  const labelList& sensitivityPatchIDs
313  );
314 
315  //- Control point sensitivities computed using face-based surface
316  //- sensitivities
318  (
319  const boundaryVectorField& faceSens,
320  const labelList& sensitivityPatchIDs
321  );
322 
323  //- Control point sensitivities computed using face-based surface
324  //- sensitivities
326  (
327  const vectorField& faceSens,
328  const label patchI,
329  const label cpI
330  );
331 
332  //- Part of control point sensitivities related to the face normal
333  //- variations
335  (
336  const label patchI,
337  const label cpI,
338  bool DimensionedNormalSens = true
339  );
340 
341  //- Get patch dx/db
343  (
344  const label patchI,
345  const label cpI
346  );
347 
348  //- Get patch dx/db
350  (
351  const label patchI,
352  const label cpI
353  );
354 
355 
356  // Cartesian coordinates
357 
358  //- Compute cartesian coordinates based on control points
359  //- and parametric coordinates
360  tmp<vectorField> coordinates(const vectorField& uVector) const;
361 
362  //- The same, for a specific point
363  vector coordinates(const vector& uVector) const;
364 
365  //- Mesh movement based on given control point movement
367  (
368  const vectorField& controlPointsMovement
369  );
370 
371  //- Boundary mesh movement based on given control point movement
373  (
374  const vectorField& controlPointsMovement,
375  const labelList& patchesToBeMoved
376  );
377 
378 
379  // Control points management
380 
381  //- Set new control points
382  // New values should be on the coordinates system original CPs
383  // were defined
384  void setControlPoints(const vectorField& newCps);
385 
386 
387  //- Bound control points movement in the boundary control points
388  //- and in certain directions if needed
390  (
391  vectorField& controlPointsMovement
392  );
393 
394  //- Compute max. displacement at the boundary
396  (
397  const vectorField& controlPointsMovement,
398  const labelList& patchesToBeMoved
399  );
400 
401 
402  // Access Functions
403 
404  // Functions triggering calculations
405 
406  //- Get points of mesh which recide within the control points
407  //- box
409 
410  //- Get map of points in box to mesh points
411  const labelList& getMap();
412 
413  //- Get map of mesh points to points in box.
414  //- Return -1 if point is outside the box
415  const labelList& getReverseMap();
416 
417  //- Get parametric coordinates
419 
420  //- Get dxCartesiandb for a certain control point
422 
423  //- Get dxCartesiandb for a certain control point on cells
425 
426  //- Get number of variables if CPs are moved symmetrically in U
427  label nUSymmetry() const;
428 
429  //- Get number of variables if CPs are moved symmetrically in V
430  label nVSymmetry() const;
431 
432  //- Get number of variables if CPs are moved symmetrically in W
433  label nWSymmetry() const;
434 
435 
436  // Inline access functions
437 
438  //- Get box name
439  inline const word& name() const;
440 
441  //- Which control points are active?
442  // A control point is active if at least one component can move
443  inline const boolList& getActiveCPs() const;
444 
445  //- Which design variables are active?
446  // Numbering is (X1,Y1,Z1), (X2,Y2,Z2) ...
447  inline const boolList& getActiveDesignVariables() const;
448 
449  //- Get control points
450  inline const vectorField& getControlPoints() const;
451 
452  //- Get confine movements
453  inline bool confineX1movement() const;
454 
455  inline bool confineX2movement() const;
456 
457  inline bool confineX3movement() const;
458 
459  //- Get basis functions
460  inline const NURBSbasis& basisU() const;
461 
462  inline const NURBSbasis& basisV() const;
463 
464  inline const NURBSbasis& basisW() const;
465 
466 
467  // Write Functions
468 
469  //- Write control points on a cartesian coordinates system for
470  //- visualization
471  void writeCps(const string="cpsFile") const;
472 
473  //- Write control points on the local coordinate system.
474  // For continuation
475  void writeCpsInDict() const;
476 
477  //- Write parametric coordinates
478  void write() const;
479 };
480 
481 
482 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
483 
484 } // End namespace Foam
485 
486 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
487 
488 #include "NURBS3DVolumeI.H"
489 
490 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
491 
492 #endif
493 
494 // ************************************************************************* //
Foam::NURBS3DVolume::computeParametricCoordinates
void computeParametricCoordinates(const vectorField &points)
Definition: NURBS3DVolume.C:118
Foam::Tensor< scalar >
Foam::NURBS3DVolume::name_
word name_
Definition: NURBS3DVolume.H:81
Foam::NURBS3DVolume::maxIter_
label maxIter_
Max iterations for Newton-Raphson.
Definition: NURBS3DVolume.H:87
Foam::NURBS3DVolume::nWSymmetry
label nWSymmetry() const
Get number of variables if CPs are moved symmetrically in W.
Definition: NURBS3DVolume.C:1819
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::NURBS3DVolume::boundVMinCPs_
boolListList3 boundVMinCPs_
Definition: NURBS3DVolume.H:127
Foam::NURBS3DVolume::basisV
const NURBSbasis & basisV() const
Definition: NURBS3DVolumeI.H:83
Foam::NURBS3DVolume::nUSymmetry
label nUSymmetry() const
Get number of variables if CPs are moved symmetrically in U.
Definition: NURBS3DVolume.C:1789
Foam::NURBS3DVolume::getMap
const labelList & getMap()
Get map of points in box to mesh points.
Definition: NURBS3DVolume.C:1637
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:832
Foam::NURBS3DVolume::write
void write() const
Write parametric coordinates.
Definition: NURBS3DVolume.C:1911
Foam::NURBS3DVolume::activeDesignVariables_
boolList activeDesignVariables_
Which design variables are changed in an optimisation.
Definition: NURBS3DVolume.H:139
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
pointPatchField.H
Foam::NURBS3DVolume::activeControlPoints_
boolList activeControlPoints_
Which of the cps are moved in an optimisation.
Definition: NURBS3DVolume.H:136
Foam::NURBS3DVolume
NURBS3DVolume morpher. Includes support functions for gradient computations Base class providing supp...
Definition: NURBS3DVolume.H:71
Foam::NURBSbasis
NURBSbasis fuction. Used to construct NURBS curves, surfaces and volumes.
Definition: NURBSbasis.H:55
Foam::NURBS3DVolume::boundUMaxCPs_
boolListList3 boundUMaxCPs_
Definition: NURBS3DVolume.H:125
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:904
Foam::NURBS3DVolume::parametricCoordinatesPtr_
autoPtr< pointVectorField > parametricCoordinatesPtr_
Parametric coordinates of pointsInBox.
Definition: NURBS3DVolume.H:107
Foam::NURBS3DVolume::~NURBS3DVolume
virtual ~NURBS3DVolume()=default
Destructor.
Foam::NURBS3DVolume::cpsFolder_
string cpsFolder_
Folder to store control points.
Definition: NURBS3DVolume.H:142
Foam::NURBS3DVolume::localSystemCoordinates_
vectorField localSystemCoordinates_
Coordinates in the local system for which CPs are defined.
Definition: NURBS3DVolume.H:110
Foam::NURBS3DVolume::computeNewBoundaryPoints
tmp< vectorField > computeNewBoundaryPoints(const vectorField &controlPointsMovement, const labelList &patchesToBeMoved)
Boundary mesh movement based on given control point movement.
Definition: NURBS3DVolume.C:1483
Foam::NURBS3DVolume::makeFolders
void makeFolders()
Create folders to store cps and derivatives.
Definition: NURBS3DVolume.C:405
Foam::NURBS3DVolume::bound
bool bound(vector &vec, scalar minValue=1e-7, scalar maxValue=0.999999)
Bound components to certain limits.
Definition: NURBS3DVolume.C:362
Foam::NURBS3DVolume::confineControlPoint
void confineControlPoint(const label cpI)
Confine all three movements for a prescribed control point.
Definition: NURBS3DVolume.C:598
Foam::NURBS3DVolume::cps_
vectorField cps_
The volumetric B-Splines control points.
Definition: NURBS3DVolume.H:97
Foam::NURBS3DVolume::basisW_
NURBSbasis basisW_
Definition: NURBS3DVolume.H:84
pointPatchFieldsFwd.H
Foam::NURBS3DVolume::updateLocalCoordinateSystem
virtual void updateLocalCoordinateSystem(const vectorField &cartesianPoints)=0
Update coordinates in the local system based on the cartesian points.
Foam::NURBS3DVolume::confineControlPointsDirections
void confineControlPointsDirections()
Confine movement in all control points for user-defined directions.
Definition: NURBS3DVolume.C:587
Foam::NURBS3DVolume::continuityRealatedConfinement
void continuityRealatedConfinement()
Confine control point movement to maintain user-defined continuity.
Definition: NURBS3DVolume.C:492
Foam::NURBS3DVolume::computeControlPointSensitivities
vectorField computeControlPointSensitivities(const pointVectorField &pointSens, const labelList &sensitivityPatchIDs)
Definition: NURBS3DVolume.C:1002
Foam::NURBS3DVolume::reverseMapPtr_
autoPtr< labelList > reverseMapPtr_
Map of mesh points to points-in-box.
Definition: NURBS3DVolume.H:104
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:50
Foam::NURBS3DVolume::getParametricCoordinates
const pointVectorField & getParametricCoordinates()
Get parametric coordinates.
Definition: NURBS3DVolume.C:1659
Foam::NURBS3DVolume::confineBoundaryControlPoints
void confineBoundaryControlPoints()
Confine movement in boundary control points if necessary.
Definition: NURBS3DVolume.C:446
Foam::NURBS3DVolume::TypeName
TypeName("NURBS3DVolume")
Runtime type information.
Foam::NURBS3DVolume::mesh_
const fvMesh & mesh_
NURBS basis functions.
Definition: NURBS3DVolume.H:80
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:62
Foam::Field< vector >
Foam::NURBS3DVolume::findPointsInBox
void findPointsInBox(const vectorField &meshPoints)
Find points within control points box.
Definition: NURBS3DVolume.C:63
Foam::NURBS3DVolume::nMaxBound_
label nMaxBound_
Definition: NURBS3DVolume.H:94
Foam::NURBS3DVolume::basisU_
NURBSbasis basisU_
Definition: NURBS3DVolume.H:82
Foam::NURBS3DVolume::getActiveDesignVariables
const boolList & getActiveDesignVariables() const
Which design variables are active?
Definition: NURBS3DVolumeI.H:47
Foam::NURBS3DVolume::determineActiveDesignVariablesAndPoints
void determineActiveDesignVariablesAndPoints()
Create lists with active design variables and control points.
Definition: NURBS3DVolume.C:414
Foam::NURBS3DVolume::mapPtr_
autoPtr< labelList > mapPtr_
Map of points-in-box to mesh points.
Definition: NURBS3DVolume.H:100
Foam::NURBS3DVolume::writeCps
void writeCps(const string="cpsFile") const
Definition: NURBS3DVolume.C:1834
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:969
Foam::NURBS3DVolume::nVSymmetry
label nVSymmetry() const
Get number of variables if CPs are moved symmetrically in V.
Definition: NURBS3DVolume.C:1804
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::NURBS3DVolume::confineX3movement_
label confineX3movement_
Definition: NURBS3DVolume.H:118
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::NURBS3DVolume::writeCpsInDict
void writeCpsInDict() const
Write control points on the local coordinate system.
Definition: NURBS3DVolume.C:1876
Foam::NURBS3DVolume::confineX3movement
bool confineX3movement() const
Definition: NURBS3DVolumeI.H:71
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:84
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::NURBS3DVolume::confineBoundaryControlPoints_
label confineBoundaryControlPoints_
Definition: NURBS3DVolume.H:120
Foam::NURBS3DVolume::confineX2movement
bool confineX2movement() const
Definition: NURBS3DVolumeI.H:65
Foam::NURBS3DVolume::getPointsInBox
tmp< vectorField > getPointsInBox()
Definition: NURBS3DVolume.C:1622
Foam::NURBS3DVolume::basisW
const NURBSbasis & basisW() const
Definition: NURBS3DVolumeI.H:89
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:868
Foam::NURBS3DVolume::transformationTensorDxDb
virtual tensor transformationTensorDxDb(label globalPointIndex)=0
Foam::NURBS3DVolume::confineX2movement_
label confineX2movement_
Definition: NURBS3DVolume.H:116
NURBSbasis.H
Foam::NURBS3DVolume::setControlPoints
void setControlPoints(const vectorField &newCps)
Set new control points.
Definition: NURBS3DVolume.C:1536
Foam::NURBS3DVolume::dndbBasedSensitivities
tmp< tensorField > dndbBasedSensitivities(const label patchI, const label cpI, bool DimensionedNormalSens=true)
Definition: NURBS3DVolume.C:1199
Foam::GeometricField::Boundary
Definition: GeometricField.H:114
boundaryFieldsFwd.H
Useful typenames for fields defined only at the boundaries.
Foam::NURBS3DVolume::getActiveCPs
const boolList & getActiveCPs() const
Which control points are active?
Definition: NURBS3DVolumeI.H:40
NURBS3DVolumeI.H
Foam::NURBS3DVolume::boundVMaxCPs_
boolListList3 boundVMaxCPs_
Definition: NURBS3DVolume.H:129
Foam::Vector< scalar >
Foam::NURBS3DVolume::transformPointToCartesian
virtual vector transformPointToCartesian(const vector &localCoordinates) const =0
Transform a point from its coordinate system to a cartesian system.
Foam::NURBS3DVolume::getReverseMap
const labelList & getReverseMap()
Definition: NURBS3DVolume.C:1648
Foam::NURBS3DVolume::boundControlPointMovement
void boundControlPointMovement(vectorField &controlPointsMovement)
Definition: NURBS3DVolume.C:1550
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:102
Foam::NURBS3DVolume::tolerance_
scalar tolerance_
Tolerance for Newton-Raphson.
Definition: NURBS3DVolume.H:90
Foam::NURBS3DVolume::JacobianUVW
tensor JacobianUVW(const vector &u) const
Jacobian matrix wrt to the volume parametric coordinates.
Definition: NURBS3DVolume.C:940
Foam::FixedList
A 1D vector of objects of type <T> with a fixed length <N>.
Definition: HashTable.H:104
Foam::NURBS3DVolume::boundUMinCPs_
boolListList3 boundUMinCPs_
Which movement components to freeze in each plane.
Definition: NURBS3DVolume.H:123
points
const pointField & points
Definition: gmvOutputHeader.H:1
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::NURBS3DVolume::coordinates
tmp< vectorField > coordinates(const vectorField &uVector) const
Definition: NURBS3DVolume.C:1365
Foam::NURBS3DVolume::NURBS3DVolume
NURBS3DVolume(const dictionary &dict, const fvMesh &mesh, bool computeParamCoors=true)
Construct from dictionary.
Definition: NURBS3DVolume.C:641
Foam::NURBS3DVolume::computeMaxBoundaryDisplacement
scalar computeMaxBoundaryDisplacement(const vectorField &controlPointsMovement, const labelList &patchesToBeMoved)
Compute max. displacement at the boundary.
Definition: NURBS3DVolume.C:1573
Foam::NURBS3DVolume::name
const word & name() const
Get box name.
Definition: NURBS3DVolumeI.H:34
Foam::NURBS3DVolume::confineX1movement
bool confineX1movement() const
Get confine movements.
Definition: NURBS3DVolumeI.H:59
Foam::NURBS3DVolume::patchDxDbFace
tmp< tensorField > patchDxDbFace(const label patchI, const label cpI)
Get patch dx/db.
Definition: NURBS3DVolume.C:1307
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:802
Foam::NURBS3DVolume::getDxDb
tmp< pointTensorField > getDxDb(const label cpI)
Get dxCartesiandb for a certain control point.
Definition: NURBS3DVolume.C:1677
Foam::NURBS3DVolume::basisV_
NURBSbasis basisV_
Definition: NURBS3DVolume.H:83
Foam::NURBS3DVolume::getDxCellsDb
tmp< volTensorField > getDxCellsDb(const label cpI)
Get dxCartesiandb for a certain control point on cells.
Definition: NURBS3DVolume.C:1721
Foam::NURBS3DVolume::basisU
const NURBSbasis & basisU() const
Get basis functions.
Definition: NURBS3DVolumeI.H:77
Foam::NURBS3DVolume::confineX1movement_
label confineX1movement_
Definition: NURBS3DVolume.H:114
Foam::NURBS3DVolume::getControlPoints
const vectorField & getControlPoints() const
Get control points.
Definition: NURBS3DVolumeI.H:53
Foam::GeometricField< vector, pointPatchField, pointMesh >
Foam::NURBS3DVolume::boundWMaxCPs_
boolListList3 boundWMaxCPs_
Definition: NURBS3DVolume.H:133
Foam::NURBS3DVolume::patchDxDb
tmp< tensorField > patchDxDb(const label patchI, const label cpI)
Get patch dx/db.
Definition: NURBS3DVolume.C:1268
Foam::NURBS3DVolume::computeNewPoints
tmp< vectorField > computeNewPoints(const vectorField &controlPointsMovement)
Mesh movement based on given control point movement.
Definition: NURBS3DVolume.C:1445
Foam::NURBS3DVolume::boundWMinCPs_
boolListList3 boundWMinCPs_
Definition: NURBS3DVolume.H:131
pointMesh.H
Foam::NURBS3DVolume::readStoredData_
bool readStoredData_
Read parametric coordinates from file if present in the folder.
Definition: NURBS3DVolume.H:145
maxValue
scalar maxValue
Definition: LISASMDCalcMethod1.H:5
Foam::NURBS3DVolume::boolListList3
List< FixedList< bool, 3 > > boolListList3
Definition: NURBS3DVolume.H:77
Foam::NURBS3DVolume::declareRunTimeSelectionTable
declareRunTimeSelectionTable(autoPtr, NURBS3DVolume, dictionary,(const dictionary &dict, const fvMesh &mesh, bool computeParamCoors),(dict, mesh, computeParamCoors))