Loading...
Searching...
No Matches
volumetricBSplinesDesignVariables.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2007-2023 PCOpt/NTUA
9 Copyright (C) 2013-2023 FOSS GP
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
33
34// * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * * //
35
36namespace Foam
37{
40 (
44 );
45}
46
47
48// * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * * //
49
51{
52 return 3*volBSplinesBase_.getTotalControlPointsNumber();
53}
54
55
58{
59 return volBSplinesBase_.getActiveDesignVariables();
60}
61
62
64{
65 // Active design variables pertaining to the CPs numbering
66 labelList activeVarsInCPs = volBSplinesBase_.getActiveDesignVariables();
67
68 // Convert the aforementioned list to the numbering of the actual design
69 // variables
71 constraint_().computeActiveDesignVariables(activeVarsInCPs);
72}
73
74
76{
77 const scalarField& dvs = *this;
78
79 // Convert design variables to CPs coordinates, stored in a scalarField
80 scalarField cpsScalar(constraint_().designVariablesToControlPoints(dvs));
81
82 // Convert the scalarField to vectorFields and transfer to morphing boxes
83 PtrList<NURBS3DVolume>& boxes = volBSplinesBase_.boxesRef();
84 label varID(0);
85 for (NURBS3DVolume& boxI : boxes)
86 {
87 vectorField cps(boxI.getControlPoints().size(), Zero);
88 for (vector& cpI : cps)
89 {
90 cpI.x() = cpsScalar[varID++];
91 cpI.y() = cpsScalar[varID++];
92 cpI.z() = cpsScalar[varID++];
93 }
94 boxI.setControlPoints(cps);
95 }
96}
97
98
100{
101 // Store CP coordinates to a scalarField
102 scalarField cpsScalar(3*volBSplinesBase_.getTotalControlPointsNumber());
103 const PtrList<NURBS3DVolume>& boxes = volBSplinesBase_.boxes();
104 label varID(0);
105 for (const NURBS3DVolume& boxI : boxes)
106 {
107 const vectorField& cps = boxI.getControlPoints();
108 for (const vector& cpI : cps)
109 {
110 cpsScalar[varID++] = cpI.x();
111 cpsScalar[varID++] = cpI.y();
112 cpsScalar[varID++] = cpI.z();
113 }
114 }
116 // Convert this scalarField to the design variables
117 scalarField::operator=
118 (constraint_().controlPointsToDesignVariables(cpsScalar));
119}
120
121
123(
124 const vectorField& controlPoints
125)
126{
127 // Store CP coordinates to a scalarField
128 scalarField cpsScalar(3*volBSplinesBase_.getTotalControlPointsNumber());
129 const PtrList<NURBS3DVolume>& boxes = volBSplinesBase_.boxes();
130 label varID(0);
131 for (const NURBS3DVolume& boxI : boxes)
132 {
133 const label nCPs(boxI.getControlPoints().size());
134 for (label cpI = 0; cpI < nCPs; ++cpI)
135 {
136 cpsScalar[varID++] = controlPoints[cpI].x();
137 cpsScalar[varID++] = controlPoints[cpI].y();
138 cpsScalar[varID++] = controlPoints[cpI].z();
139 }
140 }
142 // Convert this scalarField to the design variables
143 scalarField::operator=
144 (constraint_().controlPointsToDesignVariables(cpsScalar));
145}
146
147
149(
150 autoPtr<scalar> lowerBoundPtr,
151 autoPtr<scalar> upperBoundPtr
152)
153{
154 designVariables::readBounds(lowerBoundPtr, upperBoundPtr);
155 readBounds(lowerBounds_, "lower", -1);
157
158 // Update bounds based on the constraints - WIP
159 constraint_().computeBounds(lowerBounds_, upperBounds_);
160}
161
162
164(
166 const word& boundsName,
167 const label sign
168)
169{
170 // Read global bounds for the control points
171 if (dict_.found(boundsName + "CPBounds"))
172 {
173 bounds.reset(new scalarField(getVars().size()));
174
175 vector CPBounds(dict_.get<vector>(boundsName + "CPBounds"));
176 const PtrList<NURBS3DVolume>& boxes = volBSplinesBase_.boxesRef();
177 label varID(0);
178 for (const NURBS3DVolume& boxI : boxes)
179 {
180 const label nCPs(boxI.getControlPoints().size());
181 for (label iCP = 0; iCP < nCPs; ++iCP)
182 {
183 bounds()[varID++] = CPBounds.x();
184 bounds()[varID++] = CPBounds.y();
185 bounds()[varID++] = CPBounds.z();
186 }
187 }
188 }
189 // Read in bounds from the designVariables dictionary if present.
190 // If nonOverlappingCPs is used, the current CPs are used to determine the
191 // bounds of the CPs. If we continue from a previous solution, the current
192 // CPs are different from the initial ones and, hence, different bounds
193 // will be computed for the continuation run. Instead, read the bounds
194 // from the designVariables dict, if present
195 else if (localIOdictionary::found(boundsName + "Bounds"))
196 {
198 << "Reading " << boundsName << "Bounds from dict " << endl;
199 bounds.reset
200 (new scalarField(boundsName + "Bounds", *this, getVars().size()));
201
202 }
203 else if (nonOverlappingCPs_)
204 {
206 << "Setting " << boundsName << "Bounds from nonOverlappingCPs"
207 << endl;
208 bounds.reset(new scalarField(getVars().size()));
209 const PtrList<NURBS3DVolume>& boxes = volBSplinesBase_.boxesRef();
210 label varID(0);
211 for (const NURBS3DVolume& boxI : boxes)
212 {
213 const vectorField& cps = boxI.getControlPoints();
214 const Vector<label> nCPsDir = boxI.nCPsPerDirection();
215 vector dists(Zero);
216 for (label idir = 0; idir < 3; ++idir)
217 {
218 dists[idir] =
219 (max(cps.component(idir)) - min(cps.component(idir)))
220 /scalar(nCPsDir[idir] - 1);
221 }
222 const label nCPs(boxI.getControlPoints().size());
223 for (label iCP = 0; iCP < nCPs; ++iCP)
224 {
225 const vector& cp = cps[iCP];
226 bounds()[varID++] = cp.x() + sign*0.5*dists.x();
227 bounds()[varID++] = cp.y() + sign*0.5*dists.y();
228 bounds()[varID++] = cp.z() + sign*0.5*dists.z();
229 }
230 }
231 }
232}
233
234
236(
237 const scalarField& bounds,
238 const word& name
239) const
240{
241 if (Pstream::master())
242 {
243 const PtrList<NURBS3DVolume>& boxes = volBSplinesBase_.boxesRef();
244 label passed(0);
245 for (const NURBS3DVolume& boxI : boxes)
246 {
247 OFstream file
248 (
249 word("optimisation")/word("controlPoints")/boxI.name()
250 + name + mesh_.time().timeName() + ".csv"
251 );
252 // Write header
253 file<< "\"Points : 0\", \"Points : 1\", \"Points : 2\","
254 << "\"i\", \"j\", \"k\""<< endl;
255
256 const vectorField& cps = boxI.getControlPoints();
257 const label nCPsU = boxI.basisU().nCPs();
258 const label nCPsV = boxI.basisV().nCPs();
259 forAll(cps, cpI)
260 {
261 const label k = cpI/label(nCPsU*nCPsV);
262 const label j = (cpI - k*nCPsU*nCPsV)/nCPsU;
263 const label i = (cpI - k*nCPsU*nCPsV - j*nCPsU);
264
265 file<< bounds[3*cpI + passed] << ", "
266 << bounds[3*cpI + 1 + passed] << ", "
267 << bounds[3*cpI + 2 + passed] << ", "
268 << i << ", "
269 << j << ", "
270 << k << endl;
272 passed += 3*cps.size();
273 }
274 }
275}
276
277
279(
280 const vectorField& cpMovement
281)
282{
283 displacementMethod& dm = displMethodPtr_.ref();
284 // Are volumetric B-Splines also used to move the mesh ?
286 {
287 // Communicate the control points movement to the displacement method
288 displMethodPtr_->setControlField(cpMovement);
289 }
290 else
291 {
292 // This will also update the control point positions
293 tmp<vectorField> tnewPoints =
294 volBSplinesBase_.computeBoundaryDisplacement
295 (
296 cpMovement,
297 parametertisedPatches_.toc()
298 );
299 const vectorField& newPoints = tnewPoints();
300
302 (
304 (
305 "dx",
306 mesh_.time().timeName(),
307 mesh_,
311 ),
312 pointMesh::New(mesh_),
314 );
315
316 for (const label pI : parametertisedPatches_)
317 {
318 dx.boundaryField()[pI].setInInternalField
319 (
320 dx.primitiveFieldRef(),
321 vectorField(newPoints, mesh_.boundaryMesh()[pI].meshPoints())
322 );
323 }
324
325 // Set boundary movement of motion solver
326 displMethodPtr_->setMotionField(dx);
327 }
328}
329
330
331// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
332
333Foam::volumetricBSplinesDesignVariables::volumetricBSplinesDesignVariables
334(
335 fvMesh& mesh,
336 const dictionary& dict
337)
338:
339 shapeDesignVariables(mesh, dict),
340 localIOdictionary
341 (
342 IOobject
343 (
344 "volumetricBSplinesDesignVariables",
345 mesh.time().timeName(),
346 fileName("uniform"),
347 mesh,
348 IOobject::READ_IF_PRESENT,
349 IOobject::AUTO_WRITE
350 ),
351 word::null
352 ),
353 volBSplinesBase_(const_cast<volBSplinesBase&>(volBSplinesBase::New(mesh))),
354 nonOverlappingCPs_(dict_.getOrDefault<bool>("nonOverlappingCPs", false)),
355 updateBounds_(dict_.getOrDefault<bool>("updateBounds", true)),
356 constraint_(morphingBoxConstraint::New(mesh, dict, *this))
357{
358 // Read in design variables if present or initialise them
359 if (localIOdictionary::found("designVariables"))
360 {
361 scalarField::operator=
362 (scalarField("designVariables", *this, scalarField::size()));
363 }
364 else if (constraint_().initialiseVars())
365 {
367 }
368
369 // Set the active design variables
371
372 // Read bounds for design variables, if present
374}
375
376
377// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * //
378
381(
383)
384{
385 auto tcpMovement
386 (
388 (
389 volBSplinesBase_.getTotalControlPointsNumber(),
390 Zero
391 )
392 );
393 vectorField& cpMovement = tcpMovement.ref();
394
395 // Convert the correction pertaining to the design variables to a
396 // scalarField correction for the control points
397 const scalarField correctionCPs(constraint_().correctionCPs(correction));
398
399 // scalarField to vectorField conversion
400 forAll(cpMovement, iCP)
401 {
402 cpMovement[iCP].x() = correctionCPs[3*iCP];
403 cpMovement[iCP].y() = correctionCPs[3*iCP + 1];
404 cpMovement[iCP].z() = correctionCPs[3*iCP + 2];
406 volBSplinesBase_.boundControlPointMovement(cpMovement);
407
408 return tcpMovement;
409}
410
411
413{
414 // Get controlPoint movement from correction
415 tmp<vectorField> tcpMovement = controlPointMovement(correction);
416 const vectorField& cpMovement = tcpMovement();
417
418 // Set the field driving the displacement method
419 setDisplacement(cpMovement);
420
421 // Do the actual mesh movement
422 // Updates also the control point positions
424
425 // Update the design variables
427}
428
429
435
436
441}
442
443
456 return
457 constraint_().postProcessSens
460 adjointSens.getAdjointSolver().solverName()
461 );
462}
463
466{
467 constraint_().updateBounds(lowerBounds_, upperBounds_);
468}
469
470
472{
473 scalarField::writeEntry("designVariables", os);
474 if (lowerBounds_)
475 {
476 lowerBounds_().writeEntry("lowerBounds", os);
477 writeBounds(lowerBounds_(), "lowerBounds");
478 }
479 if (upperBounds_)
480 {
481 upperBounds_().writeEntry("upperBounds", os);
482 writeBounds(upperBounds_(), "upperBounds");
483 }
484 return constraint_().writeData(os);
485}
486
487
489(
490 const label varID
491) const
492{
493 const displacementMethod& dm = displMethodPtr_();
495 {
496 Vector<label> decomposed = volBSplinesBase_.decomposeDV(varID);
497 const label boxI = decomposed.x();
498 const label cpILocal = decomposed.y();
499 const label dir = decomposed.z();
500
501 pointTensorField dxdb(volBSplinesBase_.boxRef(boxI).getDxDb(cpILocal));
502 return unzipCol(dxdb, dir);
503 }
504 return tmp<vectorField>::New(0);
505}
506
507
509(
510 const label patchI,
511 const label varID
512) const
513{
514 Vector<label> decomposed = volBSplinesBase_.decomposeDV(varID);
515 const label boxI = decomposed.x();
516 const label cpILocal = decomposed.y();
517 const label dir = decomposed.z();
519 tensorField dxdb
520 (volBSplinesBase_.boxRef(boxI).patchDxDbFace(patchI, cpILocal));
521 return unzipCol(dxdb, dir);
522}
523
524
526(
527 const label patchI,
528 const label varID
529) const
530{
531 Vector<label> decomposed = volBSplinesBase_.decomposeDV(varID);
532 const label boxI = decomposed.x();
533 const label cpILocal = decomposed.y();
534 const label dir = decomposed.z();
535
536 tensorField dndb
537 (
538 volBSplinesBase_.boxRef(boxI).
539 dndbBasedSensitivities(patchI, cpILocal, false)
540 );
541 return unzipCol(dndb, dir);
542}
543
544
546(
547 const label patchI,
548 const label varID
549) const
550{
551 Vector<label> decomposed = volBSplinesBase_.decomposeDV(varID);
552 const label boxI = decomposed.x();
553 const label cpILocal = decomposed.y();
554 const label dir = decomposed.z();
555
556 tensorField dndb
558 volBSplinesBase_.boxRef(boxI).dndbBasedSensitivities(patchI, cpILocal)
559 );
560 return unzipCol(dndb, dir);
561}
562
563
565(
566 const label varID
567) const
568{
569 Vector<label> decomposed = volBSplinesBase_.decomposeDV(varID);
570 const label boxI = decomposed.x();
571 const label cpILocal = decomposed.y();
572 const label dir = decomposed.z();
573 NURBS3DVolume& box = volBSplinesBase_.boxRef(boxI);
574 pointVolInterpolation volPointInter(pointMesh::New(mesh_), mesh_);
575 // WIP: we compute the entire dxdb tensor corresponding to the contol point
576 // and then extract the desired direction. This is quite expensive.
577 // Specific functions returning what we want should be implemented in
578 // NURBS3DVolume to reduce the cost
579 tmp<volTensorField> dxdb = volPointInter.interpolate(box.getDxDb(cpILocal));
580 auto tdxdbDir =
582 (
584 (
585 "dxdbDir",
586 mesh_.time().timeName(),
587 mesh_,
590 ),
591 mesh_,
593 );
594 volVectorField& dxdbDir = tdxdbDir.ref();
595 unzipCol(dxdb(), vector::components(dir), dxdbDir);
596 return tdxdbDir;
597}
598
599
600// ************************************************************************* //
label k
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
void operator+=(const UList< scalar > &)
void writeEntry(const word &keyword, Ostream &os) const
tmp< Field< cmptType > > component(const direction) const
Return a component field of the field.
Definition Field.C:607
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field values.
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
@ NO_REGISTER
Do not request registration (bool: false).
@ NO_READ
Nothing to be read.
@ READ_IF_PRESENT
Reading is optional [identical to LAZY_READ].
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const Time & time() const noexcept
Return Time associated with the objectRegistry.
Definition IOobject.C:456
static FOAM_NO_DANGLING_REFERENCE const pointMesh & New(const polyMesh &mesh, Args &&... args)
NURBS3DVolume morpher. Includes support functions for gradient computations Base class providing supp...
tmp< pointTensorField > getDxDb(const label cpI)
Get dxCartesiandb for a certain control point.
Output to file stream as an OSstream, normally using std::ofstream for the actual output.
Definition OFstream.H:75
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
label size() const noexcept
Definition UList.H:706
void size(const label n)
Definition UList.H:118
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition UPstream.H:1714
Templated 3D Vector derived from VectorSpace adding construction from 3 components,...
Definition Vector.H:61
const Cmpt & x() const noexcept
Access to the vector x component.
Definition Vector.H:135
const Cmpt & z() const noexcept
Access to the vector z component.
Definition Vector.H:145
const Cmpt & y() const noexcept
Access to the vector y component.
Definition Vector.H:140
Abstract base class for adjoint-based sensitivities.
const adjointSolver & getAdjointSolver() const
Const access to adjoint solver.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
const word & name() const
Name function is needed to disambiguate those inherited from regIOobject and dictionary.
autoPtr< scalarField > upperBounds_
Upper bounds of the design variables.
autoPtr< scalar > maxInitChange_
Maximum design variables' change in the first optimisation cycle.
void readBounds(autoPtr< scalar > lowerBoundPtr=nullptr, autoPtr< scalar > upperBoundPtr=nullptr)
Read bounds for design variables, if present.
virtual const scalarField & getVars() const
Get the design variables.
autoPtr< scalarField > lowerBounds_
Lower bounds of the design variables.
labelList activeDesignVariables_
Which of the design variables will be updated.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry (const access) with the given keyword.
static const dictionary null
An empty dictionary, which is also the parent for all dictionaries.
Definition dictionary.H:487
Abstract base class for displacement methods, which are a set or wrapper classes allowing to change t...
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.
Abstract base class for defining constraints for the control points of volumetric B-Splines morphing ...
void interpolate(const GeometricField< Type, pointPatchField, pointMesh > &, GeometricField< Type, fvPatchField, volMesh > &) const
Interpolate from pointField to volField.
Abstract base class for defining design variables for shape optimisation.
virtual void resetDesignVariables()
Reset to starting point of line search.
virtual void moveMesh()
Move mesh based on displacementMethod.
autoPtr< displacementMethod > displMethodPtr_
Mesh movement mechanism.
static autoPtr< shapeDesignVariables > New(fvMesh &mesh, const dictionary &dict)
Construct and return the selected shapeDesignVariables.
labelHashSet parametertisedPatches_
Patches to be moved by the design variables.
virtual tmp< scalarField > assembleSensitivities(adjointSensitivity &adjointSens)
Add part of sensitivity derivatives related to geometry variations.
const word & solverName() const
Return the solver name.
Definition solverI.H:30
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
Class constructing a number of volumetric B-Splines boxes, read from dynamicMeshDict....
Volumetric B-Splines design variables for shape optimisation.
virtual void update(scalarField &correction)
Update design variables based on a given correction.
virtual void evolveNumber()
For design variables with a dynamic character (i.e. changing.
virtual tmp< vectorField > dSdb(const label patchI, const label varID) const
Get dSdb for given design variable and patch.
virtual scalar computeEta(scalarField &correction)
Compute eta if not set in the first step.
virtual void resetDesignVariables()
Reset to starting point of line search.
void controlPointsToDesignVariables()
Set the design variables based on the current control points.
void setDisplacement(const vectorField &cpMovement)
Set the field driving the displacement method.
virtual bool writeData(Ostream &os) const
Write fields to support continuation.
void readBounds(autoPtr< scalar > lowerBoundPtr=nullptr, autoPtr< scalar > upperBoundPtr=nullptr)
Read bounds for design variables, if present.
void designVariablesToControlPoints()
Set control points based on current design variables values.
void setActiveDesignVariables()
Set IDs of active design variables.
virtual bool globalSum() const
Whether to use global sum when computing matrix-vector products in update methods.
void writeBounds(const scalarField &bounds, const word &name) const
Write current bounds to file.
bool updateBounds_
Should the bounds of the non-overlapping control points be updated in each optimisation cycle.
tmp< vectorField > controlPointMovement(const scalarField &correction)
Transform correction of design variables to control points movement.
virtual tmp< vectorField > dxdbVol(const label varID) const
Get dxdb for all mesh points.
virtual const labelList & activeSensitivities() const
Components of the active control points.
bool nonOverlappingCPs_
Should the control points be non-overlapping.
volBSplinesBase & volBSplinesBase_
Reference to underlaying volumetric B-Splines morpher.
virtual label sensSize() const
Size of the active control points.
autoPtr< morphingBoxConstraint > constraint_
Constraint imposed on the movement of the design variables.
virtual tmp< vectorField > dxdbFace(const label patchI, const label varID) const
Get dxdb for given design variable and patch.
virtual tmp< vectorField > dndb(const label patchI, const label varID) const
Get dndb for given design variable and patch.
virtual tmp< volVectorField > dCdb(const label varID) const
Get dCdb for given design variable.
virtual tmp< scalarField > assembleSensitivities(adjointSensitivity &adjointSens)
Assemble the sensitivity derivatives, by also applying possible constraints.
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
dynamicFvMesh & mesh
OBJstream os(runTime.globalPath()/outputName)
word timeName
Definition getTimeIndex.H:3
#define DebugInfo
Report an information message using Foam::Info.
Namespace for bounding specifications. At the moment, mostly for tables.
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimless
Dimensionless.
GeometricField< tensor, pointPatchField, pointMesh > pointTensorField
dimensionedScalar sign(const dimensionedScalar &ds)
List< label > labelList
A List of labels.
Definition List.H:62
GeometricField< vector, pointPatchField, pointMesh > pointVectorField
void unzipCol(const FieldField< Field, SymmTensor< Cmpt > > &input, const direction idx, FieldField< Field, Vector< Cmpt > > &result)
Extract a symmTensor field field column (x,y,z) == (0,1,2).
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
Field< vector > vectorField
Specialisation of Field<T> for vector.
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
Field< tensor > tensorField
Specialisation of Field<T> for tensor.
bool cp(const fileName &src, const fileName &dst, const bool followLink=true)
Copy the source to the destination (recursively if necessary).
Definition POSIX.C:1065
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix by subtracting the matrix multiplied by the current fi...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
Vector< scalar > vector
Definition vector.H:57
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299