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-2019 PCOpt/NTUA
9  Copyright (C) 2013-2019 FOSS GP
10  Copyright (C) 2019 OpenCFD Ltd.
11 -------------------------------------------------------------------------------
12 License
13  This file is part of OpenFOAM.
14 
15  OpenFOAM is free software: you can redistribute it and/or modify it
16  under the terms of the GNU General Public License as published by
17  the Free Software Foundation, either version 3 of the License, or
18  (at your option) any later version.
19 
20  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
21  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
22  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
23  for more details.
24 
25  You should have received a copy of the GNU General Public License
26  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
27 
28 \*---------------------------------------------------------------------------*/
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.lookupOrDefault<label>("iters", 1000);
59  dxdbDict.lookupOrDefault<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,
129 )
130 :
131  FIBase
132  (
133  mesh,
134  dict,
135  primalVars,
136  adjointVars,
139  ),
140  //Bezier_(mesh, dict), // AJH Read locally?
141  Bezier_(mesh, mesh.lookupObject<IOdictionary>("optimisationDict")),
142  flowSens_(3*Bezier_.nBezier(), Zero),
143  dSdbSens_(3*Bezier_.nBezier(), Zero),
144  dndbSens_(3*Bezier_.nBezier(), Zero),
145  dxdbDirectSens_(3*Bezier_.nBezier(), Zero),
146  dVdbSens_(3*Bezier_.nBezier(), Zero),
147  distanceSens_(3*Bezier_.nBezier(), Zero),
148  optionsSens_(3*Bezier_.nBezier(), Zero),
149 
150  derivativesFolder_("optimisation"/type() + "Derivatives"),
151 
152  meshMovementIters_(-1),
153  meshMovementResidualLimit_(1.e-7),
154  dxdb_
155  (
157  (
158  mesh,
159  "mTilda",
161  )
162  )
163 {
164  read();
165 
166  derivatives_ = scalarField(3*Bezier_.nBezier(), Zero),
167  // Create folder to store sensitivities
168  mkDir(derivativesFolder_);
169 }
170 
171 
172 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
173 
175 {
176  // Adjoint to the eikonal equation
177  autoPtr<volTensorField> distanceSensPtr(nullptr);
178  if (includeDistance_)
179  {
180  // Solver equation
181  eikonalSolver_->solve();
182 
183  // Allocate memory and compute grad(dxdb) multiplier
184  distanceSensPtr.reset
185  (
186  createZeroFieldPtr<tensor>
187  (
188  mesh_,
189  "distanceSensPtr",
190  dimensionSet(0, 2, -3, 0, 0, 0, 0)
191  )
192  );
193  distanceSensPtr() = eikonalSolver_->getFISensitivityTerm()().T();
194  }
195 
196  const label nBezier = Bezier_.nBezier();
197  const label nDVs = 3*nBezier;
198  for (label iDV = 0; iDV < nDVs; iDV++)
199  {
200  label iCP = iDV%nBezier;
201  label idir = iDV/nBezier;
202  if
203  (
204  (idir == 0 && Bezier_.confineXmovement()[iCP])
205  || (idir == 1 && Bezier_.confineYmovement()[iCP])
206  || (idir == 2 && Bezier_.confineZmovement()[iCP])
207  )
208  {
209  continue;
210  }
211  else
212  {
213  // Flow term
214  // ~~~~~~~~~~~
215  // compute dxdb and its grad
217  const volVectorField& m = tm();
218  volTensorField gradDxDb(fvc::grad(m, "grad(dxdb)"));
219 
220  flowSens_[iDV] =
221  gSum
222  (
224  * mesh_.V()
225  );
226 
227  for (const label patchI : sensitivityPatchIDs_)
228  {
229  // Contribution from objective function
230  // term from delta(n dS)/delta b and
231  // term from delta(n)/delta b
232  //~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
233  tmp<vectorField> tdSdb =
234  Bezier_.dndbBasedSensitivities(patchI, iCP, idir);
235  const vectorField& dSdb = tdSdb();
236  tmp<vectorField> tdndb =
237  Bezier_.dndbBasedSensitivities(patchI, iCP, idir, false);
238  const vectorField& dndb = tdndb();
239  dSdbSens_[iDV] += gSum(dSfdbMult_()[patchI] & dSdb);
240  dndbSens_[iDV] += gSum(dnfdbMult_()[patchI] & dndb);
241 
242  // Contribution from objective function
243  // term from delta( x ) / delta b
244  // Only for objectives directly including
245  // x, like moments
246  //~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
247  tmp<vectorField> tdxdbFace =
248  Bezier_.dxdbFace(patchI, iCP, idir);
249  const vectorField& dxdbFace = tdxdbFace();
250  dxdbDirectSens_[iDV] +=
251  gSum(dxdbDirectMult_()[patchI] & dxdbFace);
252  }
253 
254  // Contribution from delta (V) / delta b
255  // For objectives defined as volume integrals only
256  //~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
257  dVdbSens_[iDV] =
258  gSum
259  (
261  * fvc::div(m)().primitiveField()
262  * mesh_.V()
263  );
264 
265  // Distance dependent term
266  //~~~~~~~~~~~~~~~~~~~~~~~~~
267  if (includeDistance_)
268  {
269  distanceSens_[iDV] =
270  gSum
271  (
272  (
273  distanceSensPtr().primitiveField()
274  && gradDxDb.primitiveField()
275  )
276  *mesh_.V()
277  );
278  }
279 
280  // Terms from fvOptions
281  optionsSens_[iDV] +=
283  }
284 
285  // Sum contributions
286  derivatives_ =
287  flowSens_
288  + dSdbSens_
289  + dndbSens_
291  + dVdbSens_
292  + distanceSens_
293  + optionsSens_;
294  }
295 }
296 
297 
299 {
300  flowSens_ = Zero;
301  dSdbSens_ = Zero;
302  dndbSens_ = Zero;
304  dVdbSens_ = Zero;
306  optionsSens_ = Zero;
307 
309 }
310 
311 
312 void sensitivityBezierFI::write(const word& baseName)
313 {
314  Info<< "Writing control point sensitivities to file" << endl;
315  if (Pstream::master())
316  {
317  OFstream derivFile
318  (
320  baseName + adjointVars_.solverName() + mesh_.time().timeName()
321  );
322  unsigned int widthDV = max(int(name(flowSens_.size()).size()), int(3));
323  unsigned int width = IOstream::defaultPrecision() + 7;
324  derivFile
325  << setw(widthDV) << "#dv" << " "
326  << setw(width) << "total" << " "
327  << setw(width) << "flow" << " "
328  << setw(width) << "dSdb" << " "
329  << setw(width) << "dndb" << " "
330  << setw(width) << "dxdbDirect" << " "
331  << setw(width) << "dVdb" << " "
332  << setw(width) << "distance" << endl;
333  const label nDVs = derivatives_.size();
334  const label nBezier = Bezier_.nBezier();
335  const boolListList& confineMovement = Bezier_.confineMovement();
336  label lastActive(-1);
337 
338  for (label iDV = 0; iDV < nDVs; iDV++)
339  {
340  label iCP = iDV%nBezier;
341  label idir = iDV/nBezier;
342  if (!confineMovement[idir][iCP])
343  {
344  if (iDV!=lastActive + 1) derivFile << "\n";
345  lastActive = iDV;
346  derivFile
347  << setw(widthDV) << iDV << " "
348  << setw(width) << derivatives_[iDV] << " "
349  << setw(width) << flowSens_[iDV] << " "
350  << setw(width) << dSdbSens_[iDV] << " "
351  << setw(width) << dndbSens_[iDV] << " "
352  << setw(width) << dxdbDirectSens_[iDV] << " "
353  << setw(width) << dVdbSens_[iDV] << " "
354  << setw(width) << distanceSens_[iDV] << endl;
355  }
356  }
357  }
358 }
359 
360 
361 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
362 
363 } // End namespace incompressible
364 } // End namespace Foam
365 
366 // ************************************************************************* //
Foam::incompressible::FIBase::optionsDxDbMult_
vectorField optionsDxDbMult_
dx/db multiplier coming from fvOptions
Definition: FIBaseIncompressible.H:73
Foam::incompressible::sensitivityBezierFI::clearSensitivities
virtual void clearSensitivities()
Zero sensitivity fields and their constituents.
Definition: sensitivityBezierFIIncompressible.C:298
Foam::incompressible::adjointSensitivity::derivatives_
scalarField derivatives_
Definition: adjointSensitivityIncompressible.H:84
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:54
Foam::autoPtr::reset
void reset(T *p=nullptr) noexcept
Delete managed object and set to new given pointer.
Definition: autoPtrI.H:158
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:97
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::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:62
Foam::incompressible::FIBase::clearSensitivities
virtual void clearSensitivities()
Zero sensitivity fields and their constituents.
Definition: FIBaseIncompressible.C:183
Foam::incompressible::defineTypeNameAndDebug
defineTypeNameAndDebug(adjointEikonalSolver, 0)
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:53
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
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:67
Foam::Time::timeName
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:764
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:48
Foam::incompressible::sensitivityBezierFI::distanceSens_
scalarField distanceSens_
Term depending on distance differentiation.
Definition: sensitivityBezierFIIncompressible.H:90
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:70
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
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.
Definition: dimensionSet.H:65
Foam::incompressible::FIBase::dSfdbMult_
autoPtr< boundaryVectorField > dSfdbMult_
Fields related to direct sensitivities.
Definition: FIBaseIncompressible.H:76
Foam::dictionary::lookupOrDefault
T lookupOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionary.H:1241
Foam::incompressible::sensitivityBezierFI::assembleSensitivities
virtual void assembleSensitivities()
Assemble sensitivities.
Definition: sensitivityBezierFIIncompressible.C:174
Foam::incompressibleAdjointVars
Class including all adjoint fields for incompressible flows.
Definition: incompressibleAdjointVars.H:52
Foam::fv::optionAdjointList
Definition: fvOptionAdjointList.H:59
Foam::incompressible::sensitivityBezierFI::read
void read()
Definition: sensitivityBezierFIIncompressible.C:53
Foam::incompressible::FIBase::eikonalSolver_
autoPtr< adjointEikonalSolver > eikonalSolver_
Adjoint eikonal equation solver.
Definition: FIBaseIncompressible.H:84
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::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:258
Foam::sensitivity::mesh_
const fvMesh & mesh_
Definition: sensitivity.H:69
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:62
Foam::Field< vector >
Foam::incompressible::FIBase::dnfdbMult_
autoPtr< boundaryVectorField > dnfdbMult_
Definition: FIBaseIncompressible.H:77
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::fvm::laplacian
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:48
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
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:121
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:86
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:84
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::setw
Omanip< int > setw(const int i)
Definition: IOmanip.H:199
Foam::incompressible::FIBase::dxdbDirectMult_
autoPtr< boundaryVectorField > dxdbDirectMult_
Definition: FIBaseIncompressible.H:78
Foam::dictionary::subOrEmptyDict
dictionary subOrEmptyDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX, const bool mandatory=false) const
Definition: dictionary.C:603
Foam::incompressible::sensitivityBezierFI::meshMovementResidualLimit_
scalar meshMovementResidualLimit_
Definition: sensitivityBezierFIIncompressible.H:98
Foam::OFstream
Output to file stream, using an OSstream.
Definition: OFstream.H:99
Foam::UPstream::master
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:438
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:60
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::IOstream::defaultPrecision
static unsigned int defaultPrecision()
Return the default precision.
Definition: IOstream.H:325
Foam::incompressible::FIBase::includeDistance_
bool includeDistance_
Include distance variation in sens computation.
Definition: FIBaseIncompressible.H:81
Foam::GeometricField::boundaryFieldRef
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
Definition: GeometricField.C:752
Foam::List< List< bool > >
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:76
Foam::fvMatrix::solve
SolverPerformance< Type > solve(const dictionary &)
Solve returning the solution statistics.
Definition: fvMatrixSolve.C:298
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:246
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:312
Foam::incompressible::sensitivityBezierFI::Bezier_
Bezier Bezier_
Definition: sensitivityBezierFIIncompressible.H:71
Foam::GeometricField< vector, fvPatchField, volMesh >
fvOptionsAdjoint
fv::IOoptionListAdjoint fvOptionsAdjoint(mesh)
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:95
Foam::fvMesh::V
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
Definition: fvMeshGeometry.C:190
Foam::Bezier::dndbBasedSensitivities
tmp< tensorField > dndbBasedSensitivities(const label patchI, const label cpI, bool returnDimensionedNormalSens=true) const
Definition: Bezier.C:164