sensitivityVolBSplinesFIIncompressible.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-2020 PCOpt/NTUA
9  Copyright (C) 2013-2020 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 \*---------------------------------------------------------------------------*/
29 
31 #include "pointVolInterpolation.H"
32 #include "IOmanip.H"
34 
35 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 
40 namespace incompressible
41 {
42 
43 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44 
45 defineTypeNameAndDebug(sensitivityVolBSplinesFI, 0);
47 (
48  adjointSensitivity,
51 );
52 
53 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
54 
55 sensitivityVolBSplinesFI::sensitivityVolBSplinesFI
56 (
57  const fvMesh& mesh,
58  const dictionary& dict,
59  incompressibleVars& primalVars,
60  incompressibleAdjointVars& adjointVars,
62 )
63 :
64  FIBase
65  (
66  mesh,
67  dict,
68  primalVars,
69  adjointVars,
71  ),
72  volBSplinesBase_
73  (
74  const_cast<volBSplinesBase&>(volBSplinesBase::New(mesh))
75  ),
76  flowSens_(0),
77  dSdbSens_(0),
78  dndbSens_(0),
79  dxdbDirectSens_(0),
80  dVdbSens_(0),
81  distanceSens_(0),
82  optionsSens_(0),
83  bcSens_(0),
84 
85  derivativesFolder_("optimisation"/type() + "Derivatives")
86 {
87  // No boundary field pointers need to be allocated
88 
90  derivatives_ = scalarField(3*nCPs, Zero);
91  flowSens_ = vectorField(nCPs, Zero);
92  dSdbSens_ = vectorField(nCPs, Zero);
93  dndbSens_ = vectorField(nCPs, Zero);
95  dVdbSens_ = vectorField(nCPs, Zero);
97  optionsSens_ = vectorField(nCPs, Zero);
98  bcSens_ = vectorField(nCPs, Zero);
99 
100  // Create folder to store sensitivities
102 }
103 
104 
105 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
106 
108 {
109  read();
110 
111  // Interpolation engine
113 
114  // Adjoint to the eikonal equation
115  autoPtr<volTensorField> distanceSensPtr(nullptr);
116  if (includeDistance_)
117  {
118  // Solver equation
119  eikonalSolver_->solve();
120 
121  // Allocate memory and compute grad(dxdb) multiplier
122  distanceSensPtr.reset
123  (
124  createZeroFieldPtr<tensor>
125  (
126  mesh_,
127  "distanceSensPtr",
128  dimensionSet(0, 2, -3, 0, 0, 0, 0)
129  )
130  );
131  distanceSensPtr() = eikonalSolver_->getFISensitivityTerm()().T();
132  }
133 
134  // Integration
135  label passedCPs(0);
137  forAll(boxes, iNURB)
138  {
139  const label nb(boxes[iNURB].getControlPoints().size());
140  vectorField boxSensitivities(nb, Zero);
141 
142  vectorField dxdbSens = boxes[iNURB].computeControlPointSensitivities
143  (
144  dxdbDirectMult_(),
145  sensitivityPatchIDs_.toc()
146  );
147 
148  vectorField bcSens = boxes[iNURB].computeControlPointSensitivities
149  (
150  bcDxDbMult_(),
151  sensitivityPatchIDs_.toc()
152  );
153 
154  for (label cpI = 0; cpI < nb; cpI++)
155  {
156  label globalCP = passedCPs + cpI;
157 
158  // Parameterization info
159  tmp<pointTensorField> dxdbI(boxes[iNURB].getDxDb(cpI));
160  tmp<volTensorField> tvolDxDbI(volPointInter.interpolate(dxdbI));
161  const volTensorField& volDxDbI = tvolDxDbI();
162 
163  // Chain rule used to get dx/db at cells
164  // Gives practically the same results at a much higher CPU cost
165  /*
166  tmp<volTensorField> tvolDxDbI(boxes[iNURB].getDxCellsDb(cpI));
167  volTensorField& volDxDbI = tvolDxDbI.ref();
168  */
169 
170  // Gradient of parameterization info
171  volVectorField temp
172  (
173  IOobject
174  (
175  "dxdb",
176  mesh_.time().timeName(),
177  mesh_,
180  ),
181  mesh_,
183  );
184 
185  temp.replace(0, volDxDbI.component(0));
186  temp.replace(1, volDxDbI.component(3));
187  temp.replace(2, volDxDbI.component(6));
188  volTensorField gradDxDb1(fvc::grad(temp));
189 
190  temp.replace(0, volDxDbI.component(1));
191  temp.replace(1, volDxDbI.component(4));
192  temp.replace(2, volDxDbI.component(7));
193  volTensorField gradDxDb2(fvc::grad(temp));
194 
195  temp.replace(0, volDxDbI.component(2));
196  temp.replace(1, volDxDbI.component(5));
197  temp.replace(2, volDxDbI.component(8));
198  volTensorField gradDxDb3(fvc::grad(temp));
199 
200 
201  // Volume integral terms
202  flowSens_[globalCP].x() = gSum
203  (
204  (gradDxDbMult_.primitiveField() && gradDxDb1.primitiveField())
205  *mesh_.V()
206  );
207  flowSens_[globalCP].y() = gSum
208  (
209  (gradDxDbMult_.primitiveField() && gradDxDb2.primitiveField())
210  *mesh_.V()
211  );
212  flowSens_[globalCP].z() = gSum
213  (
214  (gradDxDbMult_.primitiveField() && gradDxDb3.primitiveField())
215  *mesh_.V()
216  );
217 
218  // Contribution from objective function term from
219  // delta( n dS ) / delta b and
220  // delta ( x ) / delta b
221  // for objectives directly depending on x
222  for (const label patchI : sensitivityPatchIDs_)
223  {
224  tensorField dSdb
225  (
226  boxes[iNURB].dndbBasedSensitivities(patchI, cpI)
227  );
228  dSdbSens_[globalCP] += gSum(dSfdbMult_()[patchI] & dSdb);
229  tensorField dndb
230  (
231  boxes[iNURB].dndbBasedSensitivities(patchI, cpI, false)
232  );
233  dndbSens_[globalCP] += gSum((dnfdbMult_()[patchI] & dndb));
234  }
235 
236  // Contribution from delta (V) / delta b
237  // For objectives defined as volume integrals only
238  dVdbSens_[globalCP] +=
239  gSum
240  (
242  *fvc::div(T(volDxDbI))().primitiveField()
243  *mesh_.V()
244  );
245 
246  // Distance dependent term
247  if (includeDistance_)
248  {
249  const tensorField& distSensInt =
250  distanceSensPtr().primitiveField();
251  distanceSens_[globalCP].x() =
252  gSum
253  (
254  (distSensInt && gradDxDb1.primitiveField())*mesh_.V()
255  );
256  distanceSens_[globalCP].y() =
257  gSum
258  (
259  (distSensInt && gradDxDb2.primitiveField())*mesh_.V()
260  );
261  distanceSens_[globalCP].z() =
262  gSum
263  (
264  (distSensInt && gradDxDb3.primitiveField()) *mesh_.V()
265  );
266  }
267 
268  // Terms from fvOptions
269  optionsSens_[globalCP] +=
270  gSum((optionsDxDbMult_ & volDxDbI.primitiveField())*mesh_.V());
271 
272  // dxdbSens storage
273  dxdbDirectSens_[globalCP] = dxdbSens[cpI];
274 
275  // bcSens storage
276  bcSens_[globalCP] = bcSens[cpI];
277 
278  boxSensitivities[cpI] =
279  flowSens_[globalCP]
280  + dSdbSens_[globalCP]
281  + dndbSens_[globalCP]
282  + dVdbSens_[globalCP]
283  + distanceSens_[globalCP]
284  + dxdbDirectSens_[globalCP]
285  + optionsSens_[globalCP]
286  + bcSens_[globalCP];
287  }
288 
289  // Zero sensitivities in non-active design variables
290  boxes[iNURB].boundControlPointMovement(boxSensitivities);
291 
292  // Transfer sensitivities to global list
293  for (label cpI = 0; cpI < nb; cpI++)
294  {
295  label globalCP = passedCPs + cpI;
296  derivatives_[3*globalCP] = boxSensitivities[cpI].x();
297  derivatives_[3*globalCP + 1] = boxSensitivities[cpI].y();
298  derivatives_[3*globalCP + 2] = boxSensitivities[cpI].z();
299  }
300 
301  // Increment number of passed sensitivities
302  passedCPs += nb;
303  }
304 
305  // Zero non-active sensitivity components.
306  // For consistent output only, does not affect optimisation
315 }
316 
317 
319 {
328 
330 }
331 
332 
333 void sensitivityVolBSplinesFI::write(const word& baseName)
334 {
335  Info<< "Writing control point sensitivities to file" << endl;
336  if (Pstream::master())
337  {
338  OFstream derivFile
339  (
341  baseName + adjointVars_.solverName() + mesh_.time().timeName()
342  );
343  unsigned int widthDV
344  (
345  max(int(Foam::name(flowSens_.size()).size()), int(3))
346  );
347  unsigned int width = IOstream::defaultPrecision() + 7;
348  derivFile
349  << setw(widthDV) << "#cp" << " "
350  << setw(width) << "total::x" << " "
351  << setw(width) << "total::y" << " "
352  << setw(width) << "total::z" << " "
353  << setw(width) << "flow::x" << " "
354  << setw(width) << "flow::y" << " "
355  << setw(width) << "flow::z" << " "
356  << setw(width) << "dSdb::x" << " "
357  << setw(width) << "dSdb::y" << " "
358  << setw(width) << "dSdb::z" << " "
359  << setw(width) << "dndb::x" << " "
360  << setw(width) << "dndb::y" << " "
361  << setw(width) << "dndb::z" << " "
362  << setw(width) << "dxdbDirect::x" << " "
363  << setw(width) << "dxdbDirect::y" << " "
364  << setw(width) << "dxdbDirect::z" << " "
365  << setw(width) << "dVdb::x" << " "
366  << setw(width) << "dVdb::y" << " "
367  << setw(width) << "dVdb::z" << " "
368  << setw(width) << "distance::x" << " "
369  << setw(width) << "distance::y" << " "
370  << setw(width) << "distance::z" << " "
371  << setw(width) << "options::x" << " "
372  << setw(width) << "options::y" << " "
373  << setw(width) << "options::z" << " "
374  << setw(width) << "dvdb::x" << " "
375  << setw(width) << "dvdb::y" << " "
376  << setw(width) << "dvdb::z" << endl;
377 
378  label passedCPs(0);
379  label lastActive(-1);
381  forAll(boxes, iNURB)
382  {
383  label nb = boxes[iNURB].getControlPoints().size();
384  const boolList& activeCPs = boxes[iNURB].getActiveCPs();
385  for (label iCP = 0; iCP < nb; iCP++)
386  {
387  if (activeCPs[iCP])
388  {
389  label globalCP = passedCPs + iCP;
390  if (globalCP!=lastActive + 1) derivFile << "\n";
391  lastActive = globalCP;
392 
393  derivFile
394  << setw(widthDV) << globalCP << " "
395  << setw(width) << derivatives_[3*globalCP] << " "
396  << setw(width) << derivatives_[3*globalCP + 1] << " "
397  << setw(width) << derivatives_[3*globalCP + 2] << " "
398  << setw(width) << flowSens_[globalCP].x() << " "
399  << setw(width) << flowSens_[globalCP].y() << " "
400  << setw(width) << flowSens_[globalCP].z() << " "
401  << setw(width) << dSdbSens_[globalCP].x() << " "
402  << setw(width) << dSdbSens_[globalCP].y() << " "
403  << setw(width) << dSdbSens_[globalCP].z() << " "
404  << setw(width) << dndbSens_[globalCP].x() << " "
405  << setw(width) << dndbSens_[globalCP].y() << " "
406  << setw(width) << dndbSens_[globalCP].z() << " "
407  << setw(width) << dxdbDirectSens_[globalCP].x() << " "
408  << setw(width) << dxdbDirectSens_[globalCP].y() << " "
409  << setw(width) << dxdbDirectSens_[globalCP].z() << " "
410  << setw(width) << dVdbSens_[globalCP].x() << " "
411  << setw(width) << dVdbSens_[globalCP].y() << " "
412  << setw(width) << dVdbSens_[globalCP].z() << " "
413  << setw(width) << distanceSens_[globalCP].x() << " "
414  << setw(width) << distanceSens_[globalCP].y() << " "
415  << setw(width) << distanceSens_[globalCP].z() << " "
416  << setw(width) << optionsSens_[globalCP].x() << " "
417  << setw(width) << optionsSens_[globalCP].y() << " "
418  << setw(width) << optionsSens_[globalCP].z() << " "
419  << setw(width) << bcSens_[globalCP].x() << " "
420  << setw(width) << bcSens_[globalCP].y() << " "
421  << setw(width) << bcSens_[globalCP].z() << endl;
422  }
423  }
424  passedCPs += nb;
425  }
426  }
427 }
428 
429 
430 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
431 
432 } // End namespace incompressible
433 } // End namespace Foam
434 
435 // ************************************************************************* //
Foam::incompressible::FIBase::optionsDxDbMult_
vectorField optionsDxDbMult_
Definition: FIBaseIncompressible.H:68
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:191
Foam::incompressible::adjointSensitivity::derivatives_
scalarField derivatives_
Definition: adjointSensitivityIncompressible.H:79
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:46
Foam::objectiveManager
class for managing incompressible objective functions.
Definition: objectiveManager.H:50
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:165
Foam::GeometricField::component
tmp< GeometricField< cmptType, PatchField, GeoMesh > > component(const direction) const
Foam::variablesSet::solverName
const word & solverName() const
Definition: variablesSet.C:77
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:47
Foam::incompressible::shapeSensitivities::bcDxDbMult_
autoPtr< boundaryVectorField > bcDxDbMult_
Definition: shapeSensitivitiesIncompressible.H:64
Foam::GeometricField::replace
void replace(const direction d, const GeometricField< cmptType, PatchField, GeoMesh > &gcf)
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:63
Foam::incompressible::sensitivityVolBSplinesFI::bcSens_
vectorField bcSens_
Definition: sensitivityVolBSplinesFIIncompressible.H:87
Foam::incompressible::FIBase::clearSensitivities
virtual void clearSensitivities()
Definition: FIBaseIncompressible.C:163
Foam::incompressible::defineTypeNameAndDebug
defineTypeNameAndDebug(adjointEikonalSolver, 0)
Foam::volBSplinesBase
Class constructing a number of volumetric B-Splines boxes, read from dynamicMeshDict....
Definition: volBSplinesBase.H:55
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:57
Foam::incompressible::sensitivityVolBSplinesFI::dSdbSens_
vectorField dSdbSens_
Definition: sensitivityVolBSplinesFIIncompressible.H:68
Foam::Zero
static constexpr const zero Zero
Definition: zero.H:131
Foam::incompressible::sensitivityVolBSplinesFI::flowSens_
vectorField flowSens_
Definition: sensitivityVolBSplinesFIIncompressible.H:65
pointVolInterpolation.H
Foam::incompressible::sensitivityVolBSplinesFI::dndbSens_
vectorField dndbSens_
Definition: sensitivityVolBSplinesFIIncompressible.H:71
Foam::GeometricField::primitiveField
const Internal::FieldType & primitiveField() const
Definition: GeometricFieldI.H:46
Foam::incompressible::FIBase::gradDxDbMult_
volTensorField gradDxDbMult_
Definition: FIBaseIncompressible.H:62
Foam::MeshObject< fvMesh, UpdateableMeshObject, volBSplinesBase >::New
static const volBSplinesBase & New(const fvMesh &mesh, Args &&... args)
Definition: MeshObject.C:41
Foam::Time::timeName
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:773
Foam::incompressible::FIBase::read
void read()
Definition: FIBaseIncompressible.C:42
Foam::tensorField
Field< tensor > tensorField
Specialisation of Field<T> for tensor.
Definition: primitiveFieldsFwd.H:51
sensitivityVolBSplinesFIIncompressible.H
Foam::incompressible::New
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
Definition: VoFphaseTurbulentTransportModel.C:31
Foam::UPstream::master
static bool master(const label communicator=worldComm)
Definition: UPstream.H:453
Foam::fvc::div
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:42
Foam::incompressible::FIBase::divDxDbMult_
scalarField divDxDbMult_
Definition: FIBaseIncompressible.H:65
Foam::endl
Ostream & endl(Ostream &os)
Definition: Ostream.H:381
Foam::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:587
Foam::dimensionSet
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:104
Foam::incompressible::sensitivityVolBSplinesFI::dVdbSens_
vectorField dVdbSens_
Definition: sensitivityVolBSplinesFIIncompressible.H:78
Foam::incompressibleAdjointVars
Class including all adjoint fields for incompressible flows.
Definition: incompressibleAdjointVars.H:46
Foam::incompressible::sensitivityVolBSplinesFI::dxdbDirectSens_
vectorField dxdbDirectSens_
Definition: sensitivityVolBSplinesFIIncompressible.H:75
forAll
#define forAll(list, i)
Definition: stdFoam.H:349
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:48
Foam::incompressible::FIBase::eikonalSolver_
autoPtr< adjointEikonalSolver > eikonalSolver_
Definition: FIBaseIncompressible.H:74
Foam::incompressible::shapeSensitivities::dSfdbMult_
autoPtr< boundaryVectorField > dSfdbMult_
Definition: shapeSensitivitiesIncompressible.H:61
Foam::incompressible::addToRunTimeSelectionTable
addToRunTimeSelectionTable(adjointSensitivity, sensitivityBezier, dictionary)
Foam::pointVolInterpolation::interpolate
void interpolate(const GeometricField< Type, pointPatchField, pointMesh > &, GeometricField< Type, fvPatchField, volMesh > &) const
Definition: pointVolInterpolate.C:33
Foam::incompressible::sensitivityVolBSplinesFI::optionsSens_
vectorField optionsSens_
Definition: sensitivityVolBSplinesFIIncompressible.H:84
Foam::sensitivity::mesh_
const fvMesh & mesh_
Definition: sensitivity.H:65
Foam::Field
Generic templated field type.
Definition: Field.H:59
Foam::Info
messageStream Info
Foam::volBSplinesBase::getTotalControlPointsNumber
label getTotalControlPointsNumber() const
Definition: volBSplinesBase.C:150
Foam::pointVolInterpolation
Definition: pointVolInterpolation.H:58
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:51
Foam::volBSplinesBase::boundControlPointMovement
void boundControlPointMovement(vectorField &controlPointsMovement)
Definition: volBSplinesBase.C:245
IOmanip.H
Istream and Ostream manipulators taking arguments.
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:55
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Definition: hashSets.C:40
Foam::volBSplinesBase::boxesRef
PtrList< NURBS3DVolume > & boxesRef()
Definition: volBSplinesBase.C:114
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:119
Foam::incompressible::shapeSensitivities::dxdbDirectMult_
autoPtr< boundaryVectorField > dxdbDirectMult_
Definition: shapeSensitivitiesIncompressible.H:63
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::incompressible::adjointSensitivity::adjointVars_
incompressibleAdjointVars & adjointVars_
Definition: adjointSensitivityIncompressible.H:81
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:81
Foam
Definition: atmBoundaryLayer.C:26
Foam::setw
Omanip< int > setw(const int i)
Definition: IOmanip.H:199
Foam::OFstream
Output to file stream, using an OSstream.
Definition: OFstream.H:49
Foam::incompressible::sensitivityVolBSplinesFI::write
virtual void write(const word &baseName=word::null)
Definition: sensitivityVolBSplinesFIIncompressible.C:326
Foam::autoPtr::reset
void reset(autoPtr< T > &&other) noexcept
Foam::IOstream::defaultPrecision
static unsigned int defaultPrecision() noexcept
Definition: IOstream.H:338
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:49
Foam::incompressible::FIBase
Base class for Field Integral-based sensitivity derivatives.
Definition: FIBaseIncompressible.H:53
Foam::incompressible::sensitivityVolBSplinesFI
Calculation of adjoint based sensitivities at vol B-Splines control points using the FI approach.
Definition: sensitivityVolBSplinesFIIncompressible.H:53
Foam::incompressible::FIBase::includeDistance_
bool includeDistance_
Definition: FIBaseIncompressible.H:71
Foam::incompressible::sensitivityVolBSplinesFI::assembleSensitivities
virtual void assembleSensitivities()
Definition: sensitivityVolBSplinesFIIncompressible.C:100
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:58
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Definition: POSIX.C:717
Foam::incompressible::sensitivityVolBSplinesFI::clearSensitivities
virtual void clearSensitivities()
Definition: sensitivityVolBSplinesFIIncompressible.C:311
Foam::name
word name(const expressions::valueTypeCode typeCode)
Definition: exprTraits.C:52
Foam::fvMesh::time
const Time & time() const
Definition: fvMesh.H:276
Foam::VectorSpace< Vector< scalar >, scalar, 3 >::zero
static const Form zero
Definition: VectorSpace.H:111
Foam::incompressible::shapeSensitivities::dnfdbMult_
autoPtr< boundaryVectorField > dnfdbMult_
Definition: shapeSensitivitiesIncompressible.H:62
Foam::incompressible::sensitivityVolBSplinesFI::derivativesFolder_
fileName derivativesFolder_
Definition: sensitivityVolBSplinesFIIncompressible.H:89
Foam::GeometricField
Generic GeometricField class.
Definition: areaFieldsFwd.H:49
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:184
Foam::mkDir
bool mkDir(const fileName &pathName, mode_t mode=0777)
Definition: POSIX.C:533
Foam::incompressible::sensitivityVolBSplinesFI::volBSplinesBase_
volBSplinesBase & volBSplinesBase_
Definition: sensitivityVolBSplinesFIIncompressible.H:62
Foam::incompressible::sensitivityVolBSplinesFI::distanceSens_
vectorField distanceSens_
Definition: sensitivityVolBSplinesFIIncompressible.H:81
Foam::incompressibleVars
Base class for solution control classes.
Definition: incompressibleVars.H:48
Foam::fvMesh::V
const DimensionedField< scalar, volMesh > & V() const
Definition: fvMeshGeometry.C:172