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-2020 PCOpt/NTUA
9  Copyright (C) 2013-2020 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 "pointMesh.H"
60 #include "pointPatchField.H"
61 #include "pointPatchFieldsFwd.H"
62 #include "boundaryFieldsFwd.H"
63 
64 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
65 
66 namespace Foam
67 {
68 
69 /*---------------------------------------------------------------------------*\
70  Class NURBS3DVolume Declaration
71 \*---------------------------------------------------------------------------*/
72 
73 class NURBS3DVolume
74 {
75 protected:
76 
77  // Protected Data
78 
80 
81  const fvMesh& mesh_;
83  word name_;
84  //- NURBS basis functions
88 
89  //- Max iterations for Newton-Raphson
90  label maxIter_;
91 
92  //- Tolerance for Newton-Raphson
93  scalar tolerance_;
94 
95  //- How many times to bound parametric coordinates until deciding it is
96  //- outside the box
97  label nMaxBound_;
98 
99  //- The volumetric B-Splines control points
101 
102  //- Map of points-in-box to mesh points
104 
105  //- Map of mesh points to points-in-box
106  // Return -1 if point not within the box
108 
109  //- Parametric coordinates of pointsInBox
111 
112  //- Coordinates in the local system for which CPs are defined
114 
115  //- Confine movement in certain directions and control points. Refers
116  //- to the local system
117  label confineUMovement_;
118 
119  label confineVMovement_;
120 
121  label confineWMovement_;
122 
124 
125  //- Which movement components to freeze in each plane
127 
129 
131 
133 
135 
137 
138  //- Which of the cps are moved in an optimisation
140 
141  //- Which design variables are changed in an optimisation
143 
144  //- Folder to store control points
145  string cpsFolder_;
146 
147  //- Read parametric coordinates from file if present in the folder
148  bool readStoredData_;
149 
150 
151  // Protected Member Functions
152 
153  //- Find points within control points box
154  void findPointsInBox(const vectorField& meshPoints);
155 
156  //- Compute parametric coordinates in order to match a given set
157  //- of coordinates, based on the cps of the class
158  // Uses a Newton-Raphson loop.
159  // Argument is the points residing in the box
161 
163 
164  //- Bound components to certain limits
165  bool bound
166  (
167  vector& vec,
168  scalar minValue = 1e-7,
169  scalar maxValue = 0.999999
170  );
171 
172  //- Create lists with active design variables and control points
174 
175  //- Confine movement in boundary control points if necessary
177 
178  //- Confine control point movement to maintain user-defined continuity
180 
181  //- Confine movement in all control points for user-defined directions
183 
184  //- Confine all three movements for a prescribed control point
185  void confineControlPoint(const label cpI);
186 
187  //- Confine specific movements for a prescribed control point
188  void confineControlPoint(const label cpI, const boolVector&);
189 
190  //- Create folders to store cps and derivatives
191  void makeFolders();
192 
193  //- Transform a point from its coordinate system to a cartesian system
195  (
196  const vector& localCoordinates
197  ) const = 0;
198 
199  //- Transformation tensor for dxdb, from local coordinate system to
200  //- cartesian
201  virtual tensor transformationTensorDxDb(label globalPointIndex) = 0;
202 
203  //- Update coordinates in the local system based on the cartesian points
204  virtual void updateLocalCoordinateSystem
205  (
206  const vectorField& cartesianPoints
207  ) = 0;
208 
209 
210 public:
211 
212  //- Runtime type information
213  TypeName("NURBS3DVolume");
214 
215 
216  // Declare run-time constructor selection table
217 
219  (
220  autoPtr,
222  dictionary,
223  (
224  const dictionary& dict,
225  const fvMesh& mesh,
226  bool computeParamCoors
227  ),
228  (dict, mesh, computeParamCoors)
229  );
230 
231  // Constructors
232 
233  //- Construct from dictionary
235  (
236  const dictionary& dict,
237  const fvMesh& mesh,
238  bool computeParamCoors = true
239  );
240 
241  //- Construct as copy
243 
244 
245  // Selectors
246 
247  //- Return a reference to the selected NURBS model
249  (
250  const dictionary& dict,
251  const fvMesh& mesh,
252  bool computeParamCoors = true
253  );
254 
255 
256  //- Destructor
257  virtual ~NURBS3DVolume() = default;
258 
259 
260  // Member Functions
261 
262  // Derivatives wrt parametric coordinates
263 
264  //- Volume point derivative wrt u at point u,v,w
266  (
267  const scalar u,
268  const scalar v,
269  const scalar w
270  ) const;
271 
272  //- Volume point derivative wrt v at point u,v,w
274  (
275  const scalar u,
276  const scalar v,
277  const scalar w
278  ) const;
279 
280  //- Volume point derivative wrt w at point u,v,w
282  (
283  const scalar u,
284  const scalar v,
285  const scalar w
286  ) const;
287 
288  //- Jacobian matrix wrt to the volume parametric coordinates
289  tensor JacobianUVW(const vector& u) const;
290 
291 
292  // Derivatives wrt to control points
293 
294  //- Volume point derivative wrt to control point cpI at point u,v,w
295  // Scalar since in the local system!
296  scalar volumeDerivativeCP(const vector& u, const label cpI) const;
297 
298  //- Control point sensitivities computed using point-based surface
299  //- sensitivities
301  (
302  const pointVectorField& pointSens,
303  const labelList& sensitivityPatchIDs
304  );
305 
306  //- Control point sensitivities computed using face-based surface
307  //- sensitivities
309  (
310  const volVectorField& faceSens,
311  const labelList& sensitivityPatchIDs
312  );
313 
314  //- Control point sensitivities computed using face-based surface
315  //- sensitivities
317  (
318  const boundaryVectorField& faceSens,
319  const labelList& sensitivityPatchIDs
320  );
321 
322  //- Control point sensitivities computed using face-based surface
323  //- sensitivities
325  (
326  const vectorField& faceSens,
327  const label patchI,
328  const label cpI
329  );
330 
331  //- Part of control point sensitivities related to the face normal
332  //- variations
334  (
335  const label patchI,
336  const label cpI,
337  bool DimensionedNormalSens = true
338  );
339 
340  //- Get patch dx/db
342  (
343  const label patchI,
344  const label cpI
345  );
346 
347  //- Get patch dx/db
349  (
350  const label patchI,
351  const label cpI
352  );
353 
354 
355  // Cartesian coordinates
356 
357  //- Compute cartesian coordinates based on control points
358  //- and parametric coordinates
359  tmp<vectorField> coordinates(const vectorField& uVector) const;
360 
361  //- The same, for a specific point
362  vector coordinates(const vector& uVector) const;
363 
364  //- Mesh movement based on given control point movement
366  (
367  const vectorField& controlPointsMovement
368  );
369 
370  //- Boundary mesh movement based on given control point movement
372  (
373  const vectorField& controlPointsMovement,
374  const labelList& patchesToBeMoved
375  );
376 
377 
378  // Control points management
379 
380  //- Get control point ID from its I-J-K coordinates
381  label getCPID(const label i, const label j, const label k) const;
382 
383  //- Set new control points
384  // New values should be on the coordinates system original CPs
385  // were defined
386  void setControlPoints(const vectorField& newCps);
387 
388 
389  //- Bound control points movement in the boundary control points
390  //- and in certain directions if needed
392  (
393  vectorField& controlPointsMovement
394  );
395 
396  //- Compute max. displacement at the boundary
398  (
399  const vectorField& controlPointsMovement,
400  const labelList& patchesToBeMoved
401  );
402 
403 
404  // Access Functions
405 
406  // Functions triggering calculations
407 
408  //- Get mesh points that reside within the control points box
410 
411  //- Get map of points in box to mesh points
412  const labelList& getMap();
413 
414  //- Get map of mesh points to points in box.
415  //- Return -1 if point is outside the box
416  const labelList& getReverseMap();
417 
418  //- Get parametric coordinates
420 
421  //- Get dxCartesiandb for a certain control point
422  tmp<pointTensorField> getDxDb(const label cpI);
423 
424  //- Get dxCartesiandb for a certain control point on cells
425  tmp<volTensorField> getDxCellsDb(const label cpI);
426 
427  //- Get number of variables if CPs are moved symmetrically in U
428  label nUSymmetry() const;
429 
430  //- Get number of variables if CPs are moved symmetrically in V
431  label nVSymmetry() const;
432 
433  //- Get number of variables if CPs are moved symmetrically in W
434  label nWSymmetry() const;
435 
436 
437  // Inline access functions
438 
439  //- Get box name
440  inline const word& name() const;
441 
442  //- Which control points are active?
443  // A control point is active if at least one component can move
444  inline const boolList& getActiveCPs() const;
445 
446  //- Which design variables are active?
447  // Numbering is (X1,Y1,Z1), (X2,Y2,Z2) ...
448  inline const boolList& getActiveDesignVariables() const;
449 
450  //- Get control points
451  inline const vectorField& getControlPoints() const;
452 
453  inline vectorField& getControlPoints();
454 
455  //- Get confine movements
456  inline bool confineUMovement() const;
457 
458  inline bool confineVMovement() const;
459 
460  inline bool confineWMovement() const;
461 
462  //- Get basis functions
463  inline const NURBSbasis& basisU() const;
464 
465  inline const NURBSbasis& basisV() const;
466 
467  inline const NURBSbasis& basisW() const;
468 
469  //- Get mesh
470  inline const fvMesh& mesh() const;
471 
472  //- Get dictionary
473  inline const dictionary& dict() const;
474 
475 
476  // Write Functions
477 
478  //- Write control points on a cartesian coordinates system for
479  //- visualization
480  void writeCps
481  (
482  const fileName& = "cpsFile",
483  const bool transform = true
484  ) const;
485 
486  //- Write control points on the local coordinate system.
487  // For continuation
488  void writeCpsInDict() const;
489 
490  //- Write parametric coordinates
491  void write() const;
492 };
493 
494 
495 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
496 
497 } // End namespace Foam
498 
499 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
500 
501 #include "NURBS3DVolumeI.H"
502 
503 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
504 
505 #endif
506 
507 // ************************************************************************* //
Foam::NURBS3DVolume::computeParametricCoordinates
void computeParametricCoordinates(const vectorField &points)
Definition: NURBS3DVolume.C:105
Foam::NURBS3DVolume::dict
const dictionary & dict() const
Get dictionary.
Definition: NURBS3DVolumeI.H:106
Foam::Tensor< scalar >
Foam::NURBS3DVolume::name_
word name_
Definition: NURBS3DVolume.H:82
Foam::NURBS3DVolume::maxIter_
label maxIter_
Max iterations for Newton-Raphson.
Definition: NURBS3DVolume.H:89
Foam::NURBS3DVolume::nWSymmetry
label nWSymmetry() const
Get number of variables if CPs are moved symmetrically in W.
Definition: NURBS3DVolume.C:1764
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::NURBS3DVolume::basisV
const NURBSbasis & basisV() const
Definition: NURBS3DVolumeI.H:88
Foam::NURBS3DVolume::confineWMovement
bool confineWMovement() const
Definition: NURBS3DVolumeI.H:76
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::NURBS3DVolume::activeDesignVariables_
boolList activeDesignVariables_
Which design variables are changed in an optimisation.
Definition: NURBS3DVolume.H:141
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
pointPatchField.H
Foam::NURBS3DVolume::activeControlPoints_
boolList activeControlPoints_
Which of the cps are moved in an optimisation.
Definition: NURBS3DVolume.H:138
Foam::NURBS3DVolume
NURBS3DVolume morpher. Includes support functions for gradient computations Base class providing supp...
Definition: NURBS3DVolume.H:72
Foam::NURBS3DVolume::confineVMaxCPs_
boolVectorList confineVMaxCPs_
Definition: NURBS3DVolume.H:131
Foam::NURBSbasis
NURBSbasis function. Used to construct NURBS curves, surfaces and volumes.
Definition: NURBSbasis.H:55
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::NURBS3DVolume::parametricCoordinatesPtr_
autoPtr< pointVectorField > parametricCoordinatesPtr_
Parametric coordinates of pointsInBox.
Definition: NURBS3DVolume.H:109
Foam::NURBS3DVolume::~NURBS3DVolume
virtual ~NURBS3DVolume()=default
Destructor.
Foam::NURBS3DVolume::cpsFolder_
string cpsFolder_
Folder to store control points.
Definition: NURBS3DVolume.H:144
Foam::NURBS3DVolume::localSystemCoordinates_
vectorField localSystemCoordinates_
Coordinates in the local system for which CPs are defined.
Definition: NURBS3DVolume.H:112
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::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::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::NURBS3DVolume::basisW_
NURBSbasis basisW_
Definition: NURBS3DVolume.H:86
Foam::transform
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:521
Foam::NURBS3DVolume::dict_
dictionary dict_
Definition: NURBS3DVolume.H:81
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: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::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
Foam::NURBS3DVolume::getParametricCoordinates
const pointVectorField & getParametricCoordinates()
Get parametric coordinates.
Definition: NURBS3DVolume.C:1604
Foam::NURBS3DVolume::confineBoundaryControlPoints
void confineBoundaryControlPoints()
Confine movement in boundary control points if necessary.
Definition: NURBS3DVolume.C:433
Foam::NURBS3DVolume::TypeName
TypeName("NURBS3DVolume")
Runtime type information.
Foam::NURBS3DVolume::mesh_
const fvMesh & mesh_
Definition: NURBS3DVolume.H:80
Foam::NURBS3DVolume::boolVectorList
List< boolVector > boolVectorList
Definition: NURBS3DVolume.H:78
Foam::Field< vector >
Foam::NURBS3DVolume::confineVMovement
bool confineVMovement() const
Definition: NURBS3DVolumeI.H:70
Foam::NURBS3DVolume::findPointsInBox
void findPointsInBox(const vectorField &meshPoints)
Find points within control points box.
Definition: NURBS3DVolume.C:50
Foam::NURBS3DVolume::nMaxBound_
label nMaxBound_
Definition: NURBS3DVolume.H:96
Foam::NURBS3DVolume::basisU_
NURBSbasis basisU_
NURBS basis functions.
Definition: NURBS3DVolume.H:84
Foam::NURBS3DVolume::confineWMovement_
label confineWMovement_
Definition: NURBS3DVolume.H:120
Foam::NURBS3DVolume::getActiveDesignVariables
const boolList & getActiveDesignVariables() const
Which design variables are active?
Definition: NURBS3DVolumeI.H:46
Foam::NURBS3DVolume::confineWMaxCPs_
boolVectorList confineWMaxCPs_
Definition: NURBS3DVolume.H:135
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::NURBS3DVolume::nVSymmetry
label nVSymmetry() const
Get number of variables if CPs are moved symmetrically in V.
Definition: NURBS3DVolume.C:1749
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
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::NURBS3DVolume::confineBoundaryControlPoints_
label confineBoundaryControlPoints_
Definition: NURBS3DVolume.H:122
Foam::NURBS3DVolume::getPointsInBox
tmp< vectorField > getPointsInBox()
Get mesh points that reside within the control points box.
Definition: NURBS3DVolume.C:1567
boolVector.H
Foam::NURBS3DVolume::basisW
const NURBSbasis & basisW() const
Definition: NURBS3DVolumeI.H:94
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::NURBS3DVolume::confineVMinCPs_
boolVectorList confineVMinCPs_
Definition: NURBS3DVolume.H:129
Foam::NURBS3DVolume::transformationTensorDxDb
virtual tensor transformationTensorDxDb(label globalPointIndex)=0
Foam::NURBS3DVolume::confineUMovement_
label confineUMovement_
Definition: NURBS3DVolume.H:116
NURBSbasis.H
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
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:39
NURBS3DVolumeI.H
Foam::NURBS3DVolume::confineUMinCPs_
boolVectorList confineUMinCPs_
Which movement components to freeze in each plane.
Definition: NURBS3DVolume.H:125
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:1593
Foam::NURBS3DVolume::boundControlPointMovement
void boundControlPointMovement(vectorField &controlPointsMovement)
Definition: NURBS3DVolume.C:1495
Foam::List< boolVector >
Foam::NURBS3DVolume::tolerance_
scalar tolerance_
Tolerance for Newton-Raphson.
Definition: NURBS3DVolume.H:92
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::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:1318
Foam::NURBS3DVolume::confineVMovement_
label confineVMovement_
Definition: NURBS3DVolume.H:118
Foam::NURBS3DVolume::NURBS3DVolume
NURBS3DVolume(const dictionary &dict, const fvMesh &mesh, bool computeParamCoors=true)
Construct from dictionary.
Definition: NURBS3DVolume.C:628
Foam::NURBS3DVolume::computeMaxBoundaryDisplacement
scalar computeMaxBoundaryDisplacement(const vectorField &controlPointsMovement, const labelList &patchesToBeMoved)
Compute max. displacement at the boundary.
Definition: NURBS3DVolume.C:1518
Foam::NURBS3DVolume::name
const word & name() const
Get box name.
Definition: NURBS3DVolumeI.H:33
Foam::NURBS3DVolume::confineUMovement
bool confineUMovement() const
Get confine movements.
Definition: NURBS3DVolumeI.H:64
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::NURBS3DVolume::getDxDb
tmp< pointTensorField > getDxDb(const label cpI)
Get dxCartesiandb for a certain control point.
Definition: NURBS3DVolume.C:1622
Foam::NURBS3DVolume::basisV_
NURBSbasis basisV_
Definition: NURBS3DVolume.H:85
Foam::NURBS3DVolume::mesh
const fvMesh & mesh() const
Get mesh.
Definition: NURBS3DVolumeI.H:100
Foam::NURBS3DVolume::getDxCellsDb
tmp< volTensorField > getDxCellsDb(const label cpI)
Get dxCartesiandb for a certain control point on cells.
Definition: NURBS3DVolume.C:1666
Foam::NURBS3DVolume::basisU
const NURBSbasis & basisU() const
Get basis functions.
Definition: NURBS3DVolumeI.H:82
Foam::NURBS3DVolume::getControlPoints
const vectorField & getControlPoints() const
Get control points.
Definition: NURBS3DVolumeI.H:52
Foam::GeometricField< vector, pointPatchField, pointMesh >
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
pointMesh.H
Foam::NURBS3DVolume::readStoredData_
bool readStoredData_
Read parametric coordinates from file if present in the folder.
Definition: NURBS3DVolume.H:147
Foam::NURBS3DVolume::confineWMinCPs_
boolVectorList confineWMinCPs_
Definition: NURBS3DVolume.H:133
maxValue
scalar maxValue
Definition: LISASMDCalcMethod1.H:5
Foam::NURBS3DVolume::confineUMaxCPs_
boolVectorList confineUMaxCPs_
Definition: NURBS3DVolume.H:127
Foam::NURBS3DVolume::declareRunTimeSelectionTable
declareRunTimeSelectionTable(autoPtr, NURBS3DVolume, dictionary,(const dictionary &dict, const fvMesh &mesh, bool computeParamCoors),(dict, mesh, computeParamCoors))