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-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 
31 #include "pointVolInterpolation.H"
32 #include "fvOptionAdjoint.H"
33 #include "IOmanip.H"
35 
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40 
41 namespace incompressible
42 {
43 
44 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45 
46 defineTypeNameAndDebug(sensitivityVolBSplinesFI, 0);
48 (
49  adjointSensitivity,
50  sensitivityVolBSplinesFI,
51  dictionary
52 );
53 
54 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
55 
56 sensitivityVolBSplinesFI::sensitivityVolBSplinesFI
57 (
58  const fvMesh& mesh,
59  const dictionary& dict,
60  incompressibleVars& primalVars,
61  incompressibleAdjointVars& adjointVars,
64 )
65 :
66  FIBase
67  (
68  mesh,
69  dict,
70  primalVars,
71  adjointVars,
74  ),
75  volBSplinesBase_
76  (
78  ),
79  flowSens_(0),
80  dSdbSens_(0),
81  dndbSens_(0),
82  dxdbDirectSens_(0),
83  dVdbSens_(0),
84  distanceSens_(0),
85  optionsSens_(0),
86 
87  derivativesFolder_("optimisation"/type() + "Derivatives")
88 {
89  // No boundary field pointers need to be allocated
90 
91  label nCPs = volBSplinesBase_.getTotalControlPointsNumber();
92  derivatives_ = scalarField(3*nCPs, Zero);
93  flowSens_ = vectorField(nCPs, Zero);
94  dSdbSens_ = vectorField(nCPs, Zero);
95  dndbSens_ = vectorField(nCPs, Zero);
96  dxdbDirectSens_ = vectorField(nCPs, Zero);
97  dVdbSens_ = vectorField(nCPs, Zero);
98  distanceSens_ = vectorField(nCPs, Zero);
99  optionsSens_ = vectorField(nCPs, Zero);
100 
101  // Create folder to store sensitivities
102  mkDir(derivativesFolder_);
103 }
104 
105 
106 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
107 
109 {
110  read();
111 
112  // Interpolation engine
114 
115  // Adjoint to the eikonal equation
116  autoPtr<volTensorField> distanceSensPtr(nullptr);
117  if (includeDistance_)
118  {
119  // Solver equation
120  eikonalSolver_->solve();
121 
122  // Allocate memory and compute grad(dxdb) multiplier
123  distanceSensPtr.reset
124  (
125  createZeroFieldPtr<tensor>
126  (
127  mesh_,
128  "distanceSensPtr",
129  dimensionSet(0, 2, -3, 0, 0, 0, 0)
130  )
131  );
132  distanceSensPtr() = eikonalSolver_->getFISensitivityTerm()().T();
133  }
134 
135  // Integration
136  label passedCPs(0);
138  forAll(boxes, iNURB)
139  {
140  label nb = boxes[iNURB].getControlPoints().size();
141  vectorField boxSensitivities(nb, Zero);
142 
143  vectorField dxdbSens = boxes[iNURB].computeControlPointSensitivities
144  (
145  dxdbDirectMult_(),
146  sensitivityPatchIDs_.toc()
147  );
148 
149  for (label cpI = 0; cpI < nb; cpI++)
150  {
151  label globalCP = passedCPs + cpI;
152 
153  // Parameterization info
154  tmp<pointTensorField> dxdbI(boxes[iNURB].getDxDb(cpI));
155  tmp<volTensorField> tvolDxDbI(volPointInter.interpolate(dxdbI));
156  const volTensorField& volDxDbI = tvolDxDbI();
157 
158  // Chain rule used to get dx/db at cells
159  // Gives practically the same results at a much higher CPU cost
160  /*
161  tmp<volTensorField> tvolDxDbI(boxes[iNURB].getDxCellsDb(cpI));
162  volTensorField& volDxDbI = tvolDxDbI.ref();
163  */
164 
165  // Gradient of parameterization info
166  volVectorField temp
167  (
168  IOobject
169  (
170  "dxdb",
171  mesh_.time().timeName(),
172  mesh_,
175  ),
176  mesh_,
178  );
179 
180  temp.replace(0, volDxDbI.component(0));
181  temp.replace(1, volDxDbI.component(3));
182  temp.replace(2, volDxDbI.component(6));
183  volTensorField gradDxDb1(fvc::grad(temp));
184 
185  temp.replace(0, volDxDbI.component(1));
186  temp.replace(1, volDxDbI.component(4));
187  temp.replace(2, volDxDbI.component(7));
188  volTensorField gradDxDb2(fvc::grad(temp));
189 
190  temp.replace(0, volDxDbI.component(2));
191  temp.replace(1, volDxDbI.component(5));
192  temp.replace(2, volDxDbI.component(8));
193  volTensorField gradDxDb3(fvc::grad(temp));
194 
195 
196  // Volume integral terms
197  flowSens_[globalCP].x() = gSum
198  (
199  (gradDxDbMult_.primitiveField() && gradDxDb1.primitiveField())
200  *mesh_.V()
201  );
202  flowSens_[globalCP].y() = gSum
203  (
204  (gradDxDbMult_.primitiveField() && gradDxDb2.primitiveField())
205  *mesh_.V()
206  );
207  flowSens_[globalCP].z() = gSum
208  (
209  (gradDxDbMult_.primitiveField() && gradDxDb3.primitiveField())
210  *mesh_.V()
211  );
212 
213  // Contribution from objective function term from
214  // delta( n dS ) / delta b and
215  // delta ( x ) / delta b
216  // for objectives directly depending on x
217  for (const label patchI : sensitivityPatchIDs_)
218  {
219  tensorField dSdb
220  (
221  boxes[iNURB].dndbBasedSensitivities(patchI, cpI)
222  );
223  dSdbSens_[globalCP] += gSum(dSfdbMult_()[patchI] & dSdb);
224  tensorField dndb
225  (
226  boxes[iNURB].dndbBasedSensitivities(patchI, cpI, false)
227  );
228  dndbSens_[globalCP] += gSum((dnfdbMult_()[patchI] & dndb));
229  }
230 
231  // Contribution from delta (V) / delta b
232  // For objectives defined as volume integrals only
233  dVdbSens_[globalCP] +=
234  gSum
235  (
237  *fvc::div(T(volDxDbI))().primitiveField()
238  *mesh_.V()
239  );
240 
241  // Distance dependent term
242  if (includeDistance_)
243  {
244  const tensorField& distSensInt =
245  distanceSensPtr().primitiveField();
246  distanceSens_[globalCP].x() =
247  gSum
248  (
249  (distSensInt && gradDxDb1.primitiveField())*mesh_.V()
250  );
251  distanceSens_[globalCP].y() =
252  gSum
253  (
254  (distSensInt && gradDxDb2.primitiveField())*mesh_.V()
255  );
256  distanceSens_[globalCP].z() =
257  gSum
258  (
259  (distSensInt && gradDxDb3.primitiveField()) *mesh_.V()
260  );
261  }
262 
263  // Terms from fvOptions
264  optionsSens_[globalCP] +=
265  gSum((optionsDxDbMult_ & volDxDbI.primitiveField())*mesh_.V());
266 
267  // dxdbSens storage
268  dxdbDirectSens_[globalCP] = dxdbSens[cpI];
269 
270  boxSensitivities[cpI] =
271  flowSens_[globalCP]
272  + dSdbSens_[globalCP]
273  + dndbSens_[globalCP]
274  + dVdbSens_[globalCP]
275  + distanceSens_[globalCP]
276  + dxdbDirectSens_[globalCP]
277  + optionsSens_[globalCP];
278  }
279 
280  // Zero sensitivities in non-active design variables
281  boxes[iNURB].boundControlPointMovement(flowSens_);
282  boxes[iNURB].boundControlPointMovement(dSdbSens_);
283  boxes[iNURB].boundControlPointMovement(dndbSens_);
284  boxes[iNURB].boundControlPointMovement(dVdbSens_);
285  boxes[iNURB].boundControlPointMovement(distanceSens_);
286  boxes[iNURB].boundControlPointMovement(dxdbDirectSens_);
287  boxes[iNURB].boundControlPointMovement(optionsSens_);
288  boxes[iNURB].boundControlPointMovement(boxSensitivities);
289 
290  // Transfer sensitivities to global list
291  for (label cpI = 0; cpI < nb; cpI++)
292  {
293  label globalCP = passedCPs + cpI;
294  derivatives_[3*globalCP] = boxSensitivities[cpI].x();
295  derivatives_[3*globalCP + 1] = boxSensitivities[cpI].y();
296  derivatives_[3*globalCP + 2] = boxSensitivities[cpI].z();
297  }
298 
299  // Increment number of passed sensitivities
300  passedCPs += nb;
301  }
302 }
303 
304 
306 {
314 
316 }
317 
318 
320 {
321  Info<< "Writing control point sensitivities to file" << endl;
322  if (Pstream::master())
323  {
324  OFstream derivFile
325  (
327  baseName + adjointVars_.solverName() + mesh_.time().timeName()
328  );
329  unsigned int widthDV
330  (
331  max(int(Foam::name(flowSens_.size()).size()), int(3))
332  );
333  unsigned int width = IOstream::defaultPrecision() + 7;
334  derivFile
335  << setw(widthDV) << "#cp" << " "
336  << setw(width) << "total::x" << " "
337  << setw(width) << "total::y" << " "
338  << setw(width) << "total::z" << " "
339  << setw(width) << "flow::x" << " "
340  << setw(width) << "flow::y" << " "
341  << setw(width) << "flow::z" << " "
342  << setw(width) << "dSdb::x" << " "
343  << setw(width) << "dSdb::y" << " "
344  << setw(width) << "dSdb::z" << " "
345  << setw(width) << "dndb::x" << " "
346  << setw(width) << "dndb::y" << " "
347  << setw(width) << "dndb::z" << " "
348  << setw(width) << "dxdbDirect::x" << " "
349  << setw(width) << "dxdbDirect::y" << " "
350  << setw(width) << "dxdbDirect::z" << " "
351  << setw(width) << "dVdb::x" << " "
352  << setw(width) << "dVdb::y" << " "
353  << setw(width) << "dVdb::z" << " "
354  << setw(width) << "distance::x" << " "
355  << setw(width) << "distance::y" << " "
356  << setw(width) << "distance::z" << " "
357  << setw(width) << "options::x" << " "
358  << setw(width) << "options::y" << " "
359  << setw(width) << "options::z" << endl;
360 
361  label passedCPs(0);
362  label lastActive(-1);
364  forAll(boxes, iNURB)
365  {
366  label nb = boxes[iNURB].getControlPoints().size();
367  const boolList& activeCPs = boxes[iNURB].getActiveCPs();
368  for (label iCP = 0; iCP < nb; iCP++)
369  {
370  if (activeCPs[iCP])
371  {
372  label globalCP = passedCPs + iCP;
373  if (globalCP!=lastActive + 1) derivFile << "\n";
374  lastActive = globalCP;
375 
376  derivFile
377  << setw(widthDV) << globalCP << " "
378  << setw(width) << derivatives_[3*globalCP] << " "
379  << setw(width) << derivatives_[3*globalCP + 1] << " "
380  << setw(width) << derivatives_[3*globalCP + 2] << " "
381  << setw(width) << flowSens_[globalCP].x() << " "
382  << setw(width) << flowSens_[globalCP].y() << " "
383  << setw(width) << flowSens_[globalCP].z() << " "
384  << setw(width) << dSdbSens_[globalCP].x() << " "
385  << setw(width) << dSdbSens_[globalCP].y() << " "
386  << setw(width) << dSdbSens_[globalCP].z() << " "
387  << setw(width) << dndbSens_[globalCP].x() << " "
388  << setw(width) << dndbSens_[globalCP].y() << " "
389  << setw(width) << dndbSens_[globalCP].z() << " "
390  << setw(width) << dxdbDirectSens_[globalCP].x() << " "
391  << setw(width) << dxdbDirectSens_[globalCP].y() << " "
392  << setw(width) << dxdbDirectSens_[globalCP].z() << " "
393  << setw(width) << dVdbSens_[globalCP].x() << " "
394  << setw(width) << dVdbSens_[globalCP].y() << " "
395  << setw(width) << dVdbSens_[globalCP].z() << " "
396  << setw(width) << distanceSens_[globalCP].x() << " "
397  << setw(width) << distanceSens_[globalCP].y() << " "
398  << setw(width) << distanceSens_[globalCP].z() << " "
399  << setw(width) << optionsSens_[globalCP].x() << " "
400  << setw(width) << optionsSens_[globalCP].y() << " "
401  << setw(width) << optionsSens_[globalCP].z() << endl;
402  }
403  }
404  passedCPs += nb;
405  }
406  }
407 }
408 
409 
410 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
411 
412 } // End namespace incompressible
413 } // End namespace Foam
414 
415 // ************************************************************************* //
Foam::incompressible::FIBase::optionsDxDbMult_
vectorField optionsDxDbMult_
dx/db multiplier coming from fvOptions
Definition: FIBaseIncompressible.H:73
Foam::IOobject::NO_WRITE
Definition: IOobject.H:130
Foam::incompressible::adjointSensitivity::derivatives_
scalarField derivatives_
Definition: adjointSensitivityIncompressible.H:84
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::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::variablesSet::solverName
const word & solverName() const
Return solver name.
Definition: variablesSet.C:84
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
Foam::GeometricField::replace
void replace(const direction d, const GeometricField< cmptType, PatchField, GeoMesh > &gcf)
Replace specified field component with content from another field.
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::volBSplinesBase
Class constructing a number of volumetric B-Splines boxes, read from dynamicMeshDict....
Definition: volBSplinesBase.H:59
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
Foam::incompressible::sensitivityVolBSplinesFI::dSdbSens_
vectorField dSdbSens_
Term depending on delta(n dS)/delta b.
Definition: sensitivityVolBSplinesFIIncompressible.H:72
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::incompressible::sensitivityVolBSplinesFI::flowSens_
vectorField flowSens_
Flow related term.
Definition: sensitivityVolBSplinesFIIncompressible.H:69
pointVolInterpolation.H
Foam::incompressible::sensitivityVolBSplinesFI::dndbSens_
vectorField dndbSens_
Term depending on delta(n)/delta b.
Definition: sensitivityVolBSplinesFIIncompressible.H:75
Foam::GeometricField::primitiveField
const Internal::FieldType & primitiveField() const
Return a const-reference to the internal field.
Definition: GeometricFieldI.H:53
Foam::incompressible::FIBase::gradDxDbMult_
volTensorField gradDxDbMult_
grad(dx/db) multiplier
Definition: FIBaseIncompressible.H:67
Foam::MeshObject< fvMesh, UpdateableMeshObject, volBSplinesBase >::New
static const volBSplinesBase & New(const fvMesh &mesh, Args &&... args)
Get existing or create a new MeshObject.
Definition: MeshObject.C:48
Foam::Time::timeName
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:764
Foam::incompressible::FIBase::read
void read()
Read options and update solver pointers if necessary.
Definition: FIBaseIncompressible.C:48
sensitivityVolBSplinesFIIncompressible.H
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::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::incompressible::sensitivityVolBSplinesFI::dVdbSens_
vectorField dVdbSens_
Term depending on delta(V)/delta b.
Definition: sensitivityVolBSplinesFIIncompressible.H:82
Foam::incompressibleAdjointVars
Class including all adjoint fields for incompressible flows.
Definition: incompressibleAdjointVars.H:52
Foam::incompressible::sensitivityVolBSplinesFI::dxdbDirectSens_
vectorField dxdbDirectSens_
Definition: sensitivityVolBSplinesFIIncompressible.H:79
Foam::fv::optionAdjointList
Definition: fvOptionAdjointList.H:59
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:54
Foam::incompressible::FIBase::eikonalSolver_
autoPtr< adjointEikonalSolver > eikonalSolver_
Adjoint eikonal equation solver.
Definition: FIBaseIncompressible.H:84
Foam::incompressible::addToRunTimeSelectionTable
addToRunTimeSelectionTable(adjointSensitivity, sensitivityBezier, dictionary)
Foam::pointVolInterpolation::interpolate
void interpolate(const GeometricField< Type, pointPatchField, pointMesh > &, GeometricField< Type, fvPatchField, volMesh > &) const
Interpolate from pointField to volField.
Definition: pointVolInterpolate.C:40
Foam::incompressible::sensitivityVolBSplinesFI::optionsSens_
vectorField optionsSens_
Term depending on fvOptions.
Definition: sensitivityVolBSplinesFIIncompressible.H:88
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::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::pointVolInterpolation
Definition: pointVolInterpolation.H:62
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::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:65
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::volBSplinesBase::boxesRef
PtrList< NURBS3DVolume > & boxesRef()
Get non-const reference to the vol. B-splines boxes.
Definition: volBSplinesBase.C:121
dict
dictionary dict
Definition: searchingEngine.H:14
fvOptionAdjoint.H
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::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::incompressible::sensitivityVolBSplinesFI::write
virtual void write(const word &baseName=word::null)
Write sensitivities to file.
Definition: sensitivityVolBSplinesFIIncompressible.C:319
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
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::incompressible::sensitivityVolBSplinesFI::assembleSensitivities
virtual void assembleSensitivities()
Assemble sensitivities.
Definition: sensitivityVolBSplinesFIIncompressible.C:108
Foam::List< bool >
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::incompressible::sensitivityVolBSplinesFI::clearSensitivities
virtual void clearSensitivities()
Zero sensitivity fields and their constituents.
Definition: sensitivityVolBSplinesFIIncompressible.C:305
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:246
Foam::VectorSpace< Vector< scalar >, scalar, 3 >::zero
static const Vector< scalar > zero
Definition: VectorSpace.H:115
Foam::incompressible::sensitivityVolBSplinesFI::derivativesFolder_
fileName derivativesFolder_
Definition: sensitivityVolBSplinesFIIncompressible.H:90
Foam::GeometricField< tensor, fvPatchField, volMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:123
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::incompressible::sensitivityVolBSplinesFI::volBSplinesBase_
volBSplinesBase & volBSplinesBase_
Reference to underlaying volumetric B-Splines morpher.
Definition: sensitivityVolBSplinesFIIncompressible.H:66
Foam::incompressible::sensitivityVolBSplinesFI::distanceSens_
vectorField distanceSens_
Term depending on distance differentiation.
Definition: sensitivityVolBSplinesFIIncompressible.H:85
Foam::incompressibleVars
Base class for solution control classes.
Definition: incompressibleVars.H:54
Foam::fvMesh::V
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
Definition: fvMeshGeometry.C:190