sensitivityBezierFIIncompressible.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-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 \*---------------------------------------------------------------------------*/
29 
32 #include "IOmanip.H"
33 
34 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 
39 namespace incompressible
40 {
41 
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 
44 defineTypeNameAndDebug(sensitivityBezierFI, 0);
46 (
47  adjointSensitivity, sensitivityBezierFI, dictionary
48 );
49 
50 
51 // * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * * //
52 
54 {
55  // Laplace solution controls
56  const dictionary dxdbDict = dict_.subOrEmptyDict("dxdbSolver");
57  meshMovementIters_ = dxdbDict.getOrDefault<label>("iters", 1000);
59  dxdbDict.getOrDefault<scalar>("tolerance", 1.e-07);
60 
61  // Read variables related to the adjoint eikonal solver
62  FIBase::read();
63 }
64 
65 
67 (
68  const label iCP,
69  const label idir
70 )
71 {
72  read();
73  tmp<volVectorField> tm(new volVectorField("m", dxdb_));
74  volVectorField& m = tm.ref();
75 
76  // SOLVE FOR DXDB
77  //~~~~~~~~~~~~~~~~
78  // set boundary conditions
79  for (const label patchI : sensitivityPatchIDs_)
80  {
81  // interpolate parameterization info to faces
82  tmp<vectorField> tdxidXjFace = Bezier_.dxdbFace(patchI, iCP, idir);
83  const vectorField& dxidXjFace = tdxidXjFace();
84 
85  m.boundaryFieldRef()[patchI] == dxidXjFace;
86  }
87 
88  // iterate the adjoint to the eikonal equation
89  for (label iter = 0; iter < meshMovementIters_; iter++)
90  {
91  Info<< "Mesh Movement Propagation(direction, CP), ("
92  << idir << ", " << iCP << "), Iteration : "<< iter << endl;
93 
94  fvVectorMatrix mEqn
95  (
97  );
98 
99  // Scalar residual = max(mEqn.solve().initialResidual());
100  scalar residual = mag(mEqn.solve().initialResidual());
101 
102  Info<< "Max dxdb " << gMax(mag(m)()) << endl;
103 
104  mesh_.time().printExecutionTime(Info);
105 
106  // Check convergence
107  if (residual < meshMovementResidualLimit_)
108  {
109  Info<< "\n***Reached dxdb convergence limit, iteration " << iter
110  << "***\n\n";
111  break;
112  }
113  }
114 
115  return tm;
116 }
117 
118 
119 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
120 
121 sensitivityBezierFI::sensitivityBezierFI
122 (
123  const fvMesh& mesh,
124  const dictionary& dict,
125  incompressibleVars& primalVars,
126  incompressibleAdjointVars& adjointVars,
128 )
129 :
130  FIBase
131  (
132  mesh,
133  dict,
134  primalVars,
135  adjointVars,
137  ),
138  //Bezier_(mesh, dict), // AJH Read locally?
139  Bezier_(mesh, mesh.lookupObject<IOdictionary>("optimisationDict")),
140  flowSens_(3*Bezier_.nBezier(), Zero),
141  dSdbSens_(3*Bezier_.nBezier(), Zero),
142  dndbSens_(3*Bezier_.nBezier(), Zero),
143  dxdbDirectSens_(3*Bezier_.nBezier(), Zero),
144  dVdbSens_(3*Bezier_.nBezier(), Zero),
145  distanceSens_(3*Bezier_.nBezier(), Zero),
146  optionsSens_(3*Bezier_.nBezier(), Zero),
147  bcSens_(3*Bezier_.nBezier(), Zero),
148 
149  derivativesFolder_("optimisation"/type() + "Derivatives"),
150 
151  meshMovementIters_(-1),
152  meshMovementResidualLimit_(1.e-7),
153  dxdb_
154  (
156  (
157  mesh,
158  "mTilda",
160  )
161  )
162 {
163  read();
164 
165  derivatives_ = scalarField(3*Bezier_.nBezier(), Zero),
166  // Create folder to store sensitivities
167  mkDir(derivativesFolder_);
168 }
169 
170 
171 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
172 
174 {
175  // Adjoint to the eikonal equation
176  autoPtr<volTensorField> distanceSensPtr(nullptr);
177  if (includeDistance_)
178  {
179  // Solver equation
180  eikonalSolver_->solve();
181 
182  // Allocate memory and compute grad(dxdb) multiplier
183  distanceSensPtr.reset
184  (
185  createZeroFieldPtr<tensor>
186  (
187  mesh_,
188  "distanceSensPtr",
189  dimensionSet(0, 2, -3, 0, 0, 0, 0)
190  )
191  );
192  distanceSensPtr() = eikonalSolver_->getFISensitivityTerm()().T();
193  }
194 
195  const label nBezier = Bezier_.nBezier();
196  const label nDVs = 3*nBezier;
197  for (label iDV = 0; iDV < nDVs; iDV++)
198  {
199  label iCP = iDV%nBezier;
200  label idir = iDV/nBezier;
201  if
202  (
203  (idir == 0 && Bezier_.confineXmovement()[iCP])
204  || (idir == 1 && Bezier_.confineYmovement()[iCP])
205  || (idir == 2 && Bezier_.confineZmovement()[iCP])
206  )
207  {
208  continue;
209  }
210  else
211  {
212  // Flow term
213  // ~~~~~~~~~~~
214  // compute dxdb and its grad
216  const volVectorField& m = tm();
217  volTensorField gradDxDb(fvc::grad(m, "grad(dxdb)"));
218 
219  flowSens_[iDV] =
220  gSum
221  (
223  * mesh_.V()
224  );
225 
226  for (const label patchI : sensitivityPatchIDs_)
227  {
228  // Contribution from objective function
229  // term from delta(n dS)/delta b and
230  // term from delta(n)/delta b
231  //~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
232  tmp<vectorField> tdSdb =
233  Bezier_.dndbBasedSensitivities(patchI, iCP, idir);
234  const vectorField& dSdb = tdSdb();
235  tmp<vectorField> tdndb =
236  Bezier_.dndbBasedSensitivities(patchI, iCP, idir, false);
237  const vectorField& dndb = tdndb();
238  dSdbSens_[iDV] += gSum(dSfdbMult_()[patchI] & dSdb);
239  dndbSens_[iDV] += gSum(dnfdbMult_()[patchI] & dndb);
240 
241  // Contribution from objective function
242  // term from delta( x ) / delta b
243  // Only for objectives directly including
244  // x, like moments
245  //~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
246  tmp<vectorField> tdxdbFace =
247  Bezier_.dxdbFace(patchI, iCP, idir);
248  const vectorField& dxdbFace = tdxdbFace();
249  dxdbDirectSens_[iDV] +=
250  gSum(dxdbDirectMult_()[patchI] & dxdbFace);
251 
252  // Contribution from boundary conditions
253  bcSens_[iDV] += gSum(bcDxDbMult_()[patchI] & dxdbFace);
254  }
255 
256  // Contribution from delta (V) / delta b
257  // For objectives defined as volume integrals only
258  //~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
259  dVdbSens_[iDV] =
260  gSum
261  (
263  * fvc::div(m)().primitiveField()
264  * mesh_.V()
265  );
266 
267  // Distance dependent term
268  //~~~~~~~~~~~~~~~~~~~~~~~~~
269  if (includeDistance_)
270  {
271  distanceSens_[iDV] =
272  gSum
273  (
274  (
275  distanceSensPtr().primitiveField()
276  && gradDxDb.primitiveField()
277  )
278  *mesh_.V()
279  );
280  }
281 
282  // Terms from fvOptions
283  optionsSens_[iDV] +=
285  }
286 
287  // Sum contributions
288  derivatives_ =
289  flowSens_
290  + dSdbSens_
291  + dndbSens_
293  + dVdbSens_
294  + distanceSens_
295  + optionsSens_
296  + bcSens_;
297  }
298 }
299 
300 
302 {
303  flowSens_ = Zero;
304  dSdbSens_ = Zero;
305  dndbSens_ = Zero;
307  dVdbSens_ = Zero;
309  optionsSens_ = Zero;
310  bcSens_ = Zero;
311 
313 }
314 
315 
316 void sensitivityBezierFI::write(const word& baseName)
317 {
318  Info<< "Writing control point sensitivities to file" << endl;
319  if (Pstream::master())
320  {
321  OFstream derivFile
322  (
324  baseName + adjointVars_.solverName() + mesh_.time().timeName()
325  );
326  unsigned int widthDV = max(int(name(flowSens_.size()).size()), int(3));
327  unsigned int width = IOstream::defaultPrecision() + 7;
328  derivFile
329  << setw(widthDV) << "#dv" << " "
330  << setw(width) << "total" << " "
331  << setw(width) << "flow" << " "
332  << setw(width) << "dSdb" << " "
333  << setw(width) << "dndb" << " "
334  << setw(width) << "dxdbDirect" << " "
335  << setw(width) << "dVdb" << " "
336  << setw(width) << "distance" << " "
337  << setw(width) << "options" << " "
338  << setw(width) << "dvdb" << endl;
339  const label nDVs = derivatives_.size();
340  const label nBezier = Bezier_.nBezier();
341  const boolListList& confineMovement = Bezier_.confineMovement();
342  label lastActive(-1);
343 
344  for (label iDV = 0; iDV < nDVs; iDV++)
345  {
346  const label iCP(iDV%nBezier);
347  const label idir(iDV/nBezier);
348  if (!confineMovement[idir][iCP])
349  {
350  if (iDV!=lastActive + 1)
351  {
352  derivFile << "\n";
353  }
354  lastActive = iDV;
355  derivFile
356  << setw(widthDV) << iDV << " "
357  << setw(width) << derivatives_[iDV] << " "
358  << setw(width) << flowSens_[iDV] << " "
359  << setw(width) << dSdbSens_[iDV] << " "
360  << setw(width) << dndbSens_[iDV] << " "
361  << setw(width) << dxdbDirectSens_[iDV] << " "
362  << setw(width) << dVdbSens_[iDV] << " "
363  << setw(width) << distanceSens_[iDV] << " "
364  << setw(width) << optionsSens_[iDV] << " "
365  << setw(width) << bcSens_[iDV] << endl;
366  }
367  }
368  }
369 }
370 
371 
372 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
373 
374 } // End namespace incompressible
375 } // End namespace Foam
376 
377 // ************************************************************************* //
Foam::incompressible::FIBase::optionsDxDbMult_
vectorField optionsDxDbMult_
dx/db multiplier coming from fvOptions
Definition: FIBaseIncompressible.H:72
Foam::incompressible::sensitivityBezierFI::clearSensitivities
virtual void clearSensitivities()
Zero sensitivity fields and their constituents.
Definition: sensitivityBezierFIIncompressible.C:301
Foam::incompressible::adjointSensitivity::derivatives_
scalarField derivatives_
Definition: adjointSensitivityIncompressible.H:83
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:54
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::objectiveManager
class for managing incompressible objective functions.
Definition: objectiveManager.H:54
Foam::variablesSet::solverName
const word & solverName() const
Return solver name.
Definition: variablesSet.C:84
Foam::incompressible::sensitivityBezierFI::optionsSens_
scalarField optionsSens_
Term depending on fvOptions.
Definition: sensitivityBezierFIIncompressible.H:93
Foam::incompressible::sensitivityBezierFI::meshMovementIters_
label meshMovementIters_
Definition: sensitivityBezierFIIncompressible.H:100
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
sensitivityBezierFIIncompressible.H
Foam::incompressible::shapeSensitivities::bcDxDbMult_
autoPtr< boundaryVectorField > bcDxDbMult_
Definition: shapeSensitivitiesIncompressible.H:68
Foam::Bezier::confineXmovement
const boolList & confineXmovement() const
Confine x movement.
Definition: Bezier.C:139
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::incompressible::FIBase::clearSensitivities
virtual void clearSensitivities()
Zero sensitivity fields and their constituents.
Definition: FIBaseIncompressible.C:170
Foam::incompressible::defineTypeNameAndDebug
defineTypeNameAndDebug(adjointEikonalSolver, 0)
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::incompressible::sensitivityBezierFI::dVdbSens_
scalarField dVdbSens_
Term depending on delta(V)/delta b.
Definition: sensitivityBezierFIIncompressible.H:87
Foam::GeometricField::primitiveField
const Internal::FieldType & primitiveField() const
Return a const-reference to the internal field.
Definition: GeometricFieldI.H:53
Foam::read
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:108
Foam::incompressible::FIBase::gradDxDbMult_
volTensorField gradDxDbMult_
grad(dx/db) multiplier
Definition: FIBaseIncompressible.H:66
Foam::Time::timeName
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
Foam::incompressible::sensitivityBezierFI::dxdbDirectSens_
scalarField dxdbDirectSens_
Definition: sensitivityBezierFIIncompressible.H:84
Foam::incompressible::FIBase::read
void read()
Read options and update solver pointers if necessary.
Definition: FIBaseIncompressible.C:49
Foam::incompressible::sensitivityBezierFI::distanceSens_
scalarField distanceSens_
Term depending on distance differentiation.
Definition: sensitivityBezierFIIncompressible.H:90
Foam::UPstream::master
static bool master(const label communicator=worldComm)
Am I the master process.
Definition: UPstream.H:457
Foam::fvc::div
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:49
Foam::incompressible::FIBase::divDxDbMult_
scalarField divDxDbMult_
div(dx/db) multiplier
Definition: FIBaseIncompressible.H:69
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:594
Foam::Bezier::confineZmovement
const boolList & confineZmovement() const
Confine z movement.
Definition: Bezier.C:151
Foam::dimensionSet
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:108
Foam::incompressible::sensitivityBezierFI::assembleSensitivities
virtual void assembleSensitivities()
Assemble sensitivities.
Definition: sensitivityBezierFIIncompressible.C:173
Foam::incompressibleAdjointVars
Class including all adjoint fields for incompressible flows.
Definition: incompressibleAdjointVars.H:52
Foam::incompressible::sensitivityBezierFI::read
void read()
Definition: sensitivityBezierFIIncompressible.C:53
Foam::incompressible::FIBase::eikonalSolver_
autoPtr< adjointEikonalSolver > eikonalSolver_
Adjoint eikonal equation solver.
Definition: FIBaseIncompressible.H:78
Foam::incompressible::sensitivityBezierFI::dSdbSens_
scalarField dSdbSens_
Term depending on delta(n dS)/delta b.
Definition: sensitivityBezierFIIncompressible.H:77
Foam::Bezier::nBezier
label nBezier() const
Number of Bezier control points.
Definition: Bezier.C:127
Foam::incompressible::sensitivityBezierFI::solveMeshMovementEqn
tmp< volVectorField > solveMeshMovementEqn(const label iCP, const label idir)
Definition: sensitivityBezierFIIncompressible.C:67
Foam::incompressible::shapeSensitivities::dSfdbMult_
autoPtr< boundaryVectorField > dSfdbMult_
Fields related to direct sensitivities.
Definition: shapeSensitivitiesIncompressible.H:65
Foam::incompressible::addToRunTimeSelectionTable
addToRunTimeSelectionTable(adjointSensitivity, sensitivityBezier, dictionary)
Foam::incompressible::sensitivityBezierFI::flowSens_
scalarField flowSens_
Flow related term.
Definition: sensitivityBezierFIIncompressible.H:74
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:227
Foam::sensitivity::mesh_
const fvMesh & mesh_
Definition: sensitivity.H:69
Foam::Field< vector >
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::fvm::laplacian
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:48
Foam::incompressible::sensitivityBezierFI::dndbSens_
scalarField dndbSens_
Term depending on delta(n)/delta b.
Definition: sensitivityBezierFIIncompressible.H:80
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:58
IOmanip.H
Istream and Ostream manipulators taking arguments.
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::Bezier::dxdbFace
tmp< tensorField > dxdbFace(const label patchI, const label cpI, bool useChainRule=true) const
dxdb tensor for a Bezier parameterized patch
Definition: Bezier.C:281
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::variablesSet::autoCreateMeshMovementField
static tmp< volVectorField > autoCreateMeshMovementField(const fvMesh &mesh, const word &name, const dimensionSet &dims)
Auto create variable for mesh movement.
Definition: variablesSet.C:173
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::incompressible::shapeSensitivities::dxdbDirectMult_
autoPtr< boundaryVectorField > dxdbDirectMult_
Definition: shapeSensitivitiesIncompressible.H:67
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:85
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::setw
Omanip< int > setw(const int i)
Definition: IOmanip.H:199
Foam::dictionary::subOrEmptyDict
dictionary subOrEmptyDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX, const bool mandatory=false) const
Definition: dictionary.C:540
Foam::incompressible::sensitivityBezierFI::meshMovementResidualLimit_
scalar meshMovementResidualLimit_
Definition: sensitivityBezierFIIncompressible.H:101
Foam::OFstream
Output to file stream, using an OSstream.
Definition: OFstream.H:53
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:62
Foam::autoPtr::reset
void reset(autoPtr< T > &&other) noexcept
Delete managed object and set to new given pointer.
Definition: autoPtrI.H:117
Foam::IOstream::defaultPrecision
static unsigned int defaultPrecision() noexcept
Return the default precision.
Definition: IOstream.H:342
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
Foam::sensitivity::dict_
dictionary dict_
Definition: sensitivity.H:70
Foam::incompressible::FIBase
Base class for Field Integral-based sensitivity derivatives.
Definition: FIBaseIncompressible.H:57
Foam::incompressible::FIBase::includeDistance_
bool includeDistance_
Include distance variation in sens computation.
Definition: FIBaseIncompressible.H:75
Foam::GeometricField::boundaryFieldRef
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
Definition: GeometricField.C:783
Foam::List< List< bool > >
Foam::incompressible::sensitivityBezierFI::bcSens_
scalarField bcSens_
Term depending on the differenation of boundary conditions.
Definition: sensitivityBezierFIIncompressible.H:96
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:68
Foam::fvMatrix::solve
SolverPerformance< Type > solve(const dictionary &)
Solve returning the solution statistics.
Definition: fvMatrixSolve.C:319
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:280
Foam::Bezier::confineYmovement
const boolList & confineYmovement() const
Confine y movement.
Definition: Bezier.C:145
Foam::Bezier::confineMovement
const boolListList & confineMovement() const
Info about confining movement in all directions.
Definition: Bezier.C:157
Foam::incompressible::sensitivityBezierFI::write
virtual void write(const word &baseName=word::null)
Write sensitivities to file.
Definition: sensitivityBezierFIIncompressible.C:316
Foam::incompressible::sensitivityBezierFI::Bezier_
Bezier Bezier_
Definition: sensitivityBezierFIIncompressible.H:71
Foam::incompressible::shapeSensitivities::dnfdbMult_
autoPtr< boundaryVectorField > dnfdbMult_
Definition: shapeSensitivitiesIncompressible.H:66
Foam::dictionary::getOrDefault
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:148
Foam::GeometricField< vector, fvPatchField, volMesh >
Foam::mkDir
bool mkDir(const fileName &pathName, mode_t mode=0777)
Make a directory and return an error if it could not be created.
Definition: MSwindows.C:507
Foam::gMax
Type gMax(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:592
Foam::incompressibleVars
Base class for solution control classes.
Definition: incompressibleVars.H:54
Foam::incompressible::sensitivityBezierFI::derivativesFolder_
fileName derivativesFolder_
Definition: sensitivityBezierFIIncompressible.H:98
Foam::fvMesh::V
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
Definition: fvMeshGeometry.C:179
Foam::Bezier::dndbBasedSensitivities
tmp< tensorField > dndbBasedSensitivities(const label patchI, const label cpI, bool returnDimensionedNormalSens=true) const
Definition: Bezier.C:164