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 "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  bcSens_(0),
87 
88  derivativesFolder_("optimisation"/type() + "Derivatives")
89 {
90  // No boundary field pointers need to be allocated
91 
92  label nCPs = volBSplinesBase_.getTotalControlPointsNumber();
93  derivatives_ = scalarField(3*nCPs, Zero);
94  flowSens_ = vectorField(nCPs, Zero);
95  dSdbSens_ = vectorField(nCPs, Zero);
96  dndbSens_ = vectorField(nCPs, Zero);
97  dxdbDirectSens_ = vectorField(nCPs, Zero);
98  dVdbSens_ = vectorField(nCPs, Zero);
99  distanceSens_ = vectorField(nCPs, Zero);
100  optionsSens_ = vectorField(nCPs, Zero);
101  bcSens_ = vectorField(nCPs, Zero);
102 
103  // Create folder to store sensitivities
104  mkDir(derivativesFolder_);
105 }
106 
107 
108 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
109 
111 {
112  read();
113 
114  // Interpolation engine
116 
117  // Adjoint to the eikonal equation
118  autoPtr<volTensorField> distanceSensPtr(nullptr);
119  if (includeDistance_)
120  {
121  // Solver equation
122  eikonalSolver_->solve();
123 
124  // Allocate memory and compute grad(dxdb) multiplier
125  distanceSensPtr.reset
126  (
127  createZeroFieldPtr<tensor>
128  (
129  mesh_,
130  "distanceSensPtr",
131  dimensionSet(0, 2, -3, 0, 0, 0, 0)
132  )
133  );
134  distanceSensPtr() = eikonalSolver_->getFISensitivityTerm()().T();
135  }
136 
137  // Integration
138  label passedCPs(0);
140  forAll(boxes, iNURB)
141  {
142  const label nb(boxes[iNURB].getControlPoints().size());
143  vectorField boxSensitivities(nb, Zero);
144 
145  vectorField dxdbSens = boxes[iNURB].computeControlPointSensitivities
146  (
147  dxdbDirectMult_(),
148  sensitivityPatchIDs_.toc()
149  );
150 
151  vectorField bcSens = boxes[iNURB].computeControlPointSensitivities
152  (
153  bcDxDbMult_(),
154  sensitivityPatchIDs_.toc()
155  );
156 
157  for (label cpI = 0; cpI < nb; cpI++)
158  {
159  label globalCP = passedCPs + cpI;
160 
161  // Parameterization info
162  tmp<pointTensorField> dxdbI(boxes[iNURB].getDxDb(cpI));
163  tmp<volTensorField> tvolDxDbI(volPointInter.interpolate(dxdbI));
164  const volTensorField& volDxDbI = tvolDxDbI();
165 
166  // Chain rule used to get dx/db at cells
167  // Gives practically the same results at a much higher CPU cost
168  /*
169  tmp<volTensorField> tvolDxDbI(boxes[iNURB].getDxCellsDb(cpI));
170  volTensorField& volDxDbI = tvolDxDbI.ref();
171  */
172 
173  // Gradient of parameterization info
174  volVectorField temp
175  (
176  IOobject
177  (
178  "dxdb",
179  mesh_.time().timeName(),
180  mesh_,
183  ),
184  mesh_,
186  );
187 
188  temp.replace(0, volDxDbI.component(0));
189  temp.replace(1, volDxDbI.component(3));
190  temp.replace(2, volDxDbI.component(6));
191  volTensorField gradDxDb1(fvc::grad(temp));
192 
193  temp.replace(0, volDxDbI.component(1));
194  temp.replace(1, volDxDbI.component(4));
195  temp.replace(2, volDxDbI.component(7));
196  volTensorField gradDxDb2(fvc::grad(temp));
197 
198  temp.replace(0, volDxDbI.component(2));
199  temp.replace(1, volDxDbI.component(5));
200  temp.replace(2, volDxDbI.component(8));
201  volTensorField gradDxDb3(fvc::grad(temp));
202 
203 
204  // Volume integral terms
205  flowSens_[globalCP].x() = gSum
206  (
207  (gradDxDbMult_.primitiveField() && gradDxDb1.primitiveField())
208  *mesh_.V()
209  );
210  flowSens_[globalCP].y() = gSum
211  (
212  (gradDxDbMult_.primitiveField() && gradDxDb2.primitiveField())
213  *mesh_.V()
214  );
215  flowSens_[globalCP].z() = gSum
216  (
217  (gradDxDbMult_.primitiveField() && gradDxDb3.primitiveField())
218  *mesh_.V()
219  );
220 
221  // Contribution from objective function term from
222  // delta( n dS ) / delta b and
223  // delta ( x ) / delta b
224  // for objectives directly depending on x
225  for (const label patchI : sensitivityPatchIDs_)
226  {
227  tensorField dSdb
228  (
229  boxes[iNURB].dndbBasedSensitivities(patchI, cpI)
230  );
231  dSdbSens_[globalCP] += gSum(dSfdbMult_()[patchI] & dSdb);
232  tensorField dndb
233  (
234  boxes[iNURB].dndbBasedSensitivities(patchI, cpI, false)
235  );
236  dndbSens_[globalCP] += gSum((dnfdbMult_()[patchI] & dndb));
237  }
238 
239  // Contribution from delta (V) / delta b
240  // For objectives defined as volume integrals only
241  dVdbSens_[globalCP] +=
242  gSum
243  (
245  *fvc::div(T(volDxDbI))().primitiveField()
246  *mesh_.V()
247  );
248 
249  // Distance dependent term
250  if (includeDistance_)
251  {
252  const tensorField& distSensInt =
253  distanceSensPtr().primitiveField();
254  distanceSens_[globalCP].x() =
255  gSum
256  (
257  (distSensInt && gradDxDb1.primitiveField())*mesh_.V()
258  );
259  distanceSens_[globalCP].y() =
260  gSum
261  (
262  (distSensInt && gradDxDb2.primitiveField())*mesh_.V()
263  );
264  distanceSens_[globalCP].z() =
265  gSum
266  (
267  (distSensInt && gradDxDb3.primitiveField()) *mesh_.V()
268  );
269  }
270 
271  // Terms from fvOptions
272  optionsSens_[globalCP] +=
273  gSum((optionsDxDbMult_ & volDxDbI.primitiveField())*mesh_.V());
274 
275  // dxdbSens storage
276  dxdbDirectSens_[globalCP] = dxdbSens[cpI];
277 
278  // bcSens storage
279  bcSens_[globalCP] = bcSens[cpI];
280 
281  boxSensitivities[cpI] =
282  flowSens_[globalCP]
283  + dSdbSens_[globalCP]
284  + dndbSens_[globalCP]
285  + dVdbSens_[globalCP]
286  + distanceSens_[globalCP]
287  + dxdbDirectSens_[globalCP]
288  + optionsSens_[globalCP]
289  + bcSens_[globalCP];
290  }
291 
292  // Zero sensitivities in non-active design variables
293  boxes[iNURB].boundControlPointMovement(boxSensitivities);
294 
295  // Transfer sensitivities to global list
296  for (label cpI = 0; cpI < nb; cpI++)
297  {
298  label globalCP = passedCPs + cpI;
299  derivatives_[3*globalCP] = boxSensitivities[cpI].x();
300  derivatives_[3*globalCP + 1] = boxSensitivities[cpI].y();
301  derivatives_[3*globalCP + 2] = boxSensitivities[cpI].z();
302  }
303 
304  // Increment number of passed sensitivities
305  passedCPs += nb;
306  }
307 
308  // Zero non-active sensitivity components.
309  // For consistent output only, does not affect optimisation
318 }
319 
320 
322 {
331 
333 }
334 
335 
337 {
338  Info<< "Writing control point sensitivities to file" << endl;
339  if (Pstream::master())
340  {
341  OFstream derivFile
342  (
344  baseName + adjointVars_.solverName() + mesh_.time().timeName()
345  );
346  unsigned int widthDV
347  (
348  max(int(Foam::name(flowSens_.size()).size()), int(3))
349  );
350  unsigned int width = IOstream::defaultPrecision() + 7;
351  derivFile
352  << setw(widthDV) << "#cp" << " "
353  << setw(width) << "total::x" << " "
354  << setw(width) << "total::y" << " "
355  << setw(width) << "total::z" << " "
356  << setw(width) << "flow::x" << " "
357  << setw(width) << "flow::y" << " "
358  << setw(width) << "flow::z" << " "
359  << setw(width) << "dSdb::x" << " "
360  << setw(width) << "dSdb::y" << " "
361  << setw(width) << "dSdb::z" << " "
362  << setw(width) << "dndb::x" << " "
363  << setw(width) << "dndb::y" << " "
364  << setw(width) << "dndb::z" << " "
365  << setw(width) << "dxdbDirect::x" << " "
366  << setw(width) << "dxdbDirect::y" << " "
367  << setw(width) << "dxdbDirect::z" << " "
368  << setw(width) << "dVdb::x" << " "
369  << setw(width) << "dVdb::y" << " "
370  << setw(width) << "dVdb::z" << " "
371  << setw(width) << "distance::x" << " "
372  << setw(width) << "distance::y" << " "
373  << setw(width) << "distance::z" << " "
374  << setw(width) << "options::x" << " "
375  << setw(width) << "options::y" << " "
376  << setw(width) << "options::z" << " "
377  << setw(width) << "dvdb::x" << " "
378  << setw(width) << "dvdb::y" << " "
379  << setw(width) << "dvdb::z" << endl;
380 
381  label passedCPs(0);
382  label lastActive(-1);
384  forAll(boxes, iNURB)
385  {
386  label nb = boxes[iNURB].getControlPoints().size();
387  const boolList& activeCPs = boxes[iNURB].getActiveCPs();
388  for (label iCP = 0; iCP < nb; iCP++)
389  {
390  if (activeCPs[iCP])
391  {
392  label globalCP = passedCPs + iCP;
393  if (globalCP!=lastActive + 1) derivFile << "\n";
394  lastActive = globalCP;
395 
396  derivFile
397  << setw(widthDV) << globalCP << " "
398  << setw(width) << derivatives_[3*globalCP] << " "
399  << setw(width) << derivatives_[3*globalCP + 1] << " "
400  << setw(width) << derivatives_[3*globalCP + 2] << " "
401  << setw(width) << flowSens_[globalCP].x() << " "
402  << setw(width) << flowSens_[globalCP].y() << " "
403  << setw(width) << flowSens_[globalCP].z() << " "
404  << setw(width) << dSdbSens_[globalCP].x() << " "
405  << setw(width) << dSdbSens_[globalCP].y() << " "
406  << setw(width) << dSdbSens_[globalCP].z() << " "
407  << setw(width) << dndbSens_[globalCP].x() << " "
408  << setw(width) << dndbSens_[globalCP].y() << " "
409  << setw(width) << dndbSens_[globalCP].z() << " "
410  << setw(width) << dxdbDirectSens_[globalCP].x() << " "
411  << setw(width) << dxdbDirectSens_[globalCP].y() << " "
412  << setw(width) << dxdbDirectSens_[globalCP].z() << " "
413  << setw(width) << dVdbSens_[globalCP].x() << " "
414  << setw(width) << dVdbSens_[globalCP].y() << " "
415  << setw(width) << dVdbSens_[globalCP].z() << " "
416  << setw(width) << distanceSens_[globalCP].x() << " "
417  << setw(width) << distanceSens_[globalCP].y() << " "
418  << setw(width) << distanceSens_[globalCP].z() << " "
419  << setw(width) << optionsSens_[globalCP].x() << " "
420  << setw(width) << optionsSens_[globalCP].y() << " "
421  << setw(width) << optionsSens_[globalCP].z() << " "
422  << setw(width) << bcSens_[globalCP].x() << " "
423  << setw(width) << bcSens_[globalCP].y() << " "
424  << setw(width) << bcSens_[globalCP].z() << endl;
425  }
426  }
427  passedCPs += nb;
428  }
429  }
430 }
431 
432 
433 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
434 
435 } // End namespace incompressible
436 } // End namespace Foam
437 
438 // ************************************************************************* //
Foam::incompressible::FIBase::optionsDxDbMult_
vectorField optionsDxDbMult_
dx/db multiplier coming from fvOptions
Definition: FIBaseIncompressible.H:72
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:109
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::incompressible::shapeSensitivities::bcDxDbMult_
autoPtr< boundaryVectorField > bcDxDbMult_
Definition: shapeSensitivitiesIncompressible.H:68
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::sensitivityVolBSplinesFI::bcSens_
vectorField bcSens_
Term depending on the differentiation of boundary conditions.
Definition: sensitivityVolBSplinesFIIncompressible.H:91
Foam::incompressible::FIBase::clearSensitivities
virtual void clearSensitivities()
Zero sensitivity fields and their constituents.
Definition: FIBaseIncompressible.C:172
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 (0)
Definition: zero.H:131
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:66
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:785
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:69
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
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::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:296
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:78
Foam::incompressible::shapeSensitivities::dSfdbMult_
autoPtr< boundaryVectorField > dSfdbMult_
Fields related to direct sensitivities.
Definition: shapeSensitivitiesIncompressible.H:65
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::Field< vector >
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
Foam::volBSplinesBase::boundControlPointMovement
void boundControlPointMovement(vectorField &controlPointsMovement)
Bound control points movement.
Definition: volBSplinesBase.C:252
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:62
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
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: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::OFstream
Output to file stream, using an OSstream.
Definition: OFstream.H:87
Foam::UPstream::master
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:439
Foam::incompressible::sensitivityVolBSplinesFI::write
virtual void write(const word &baseName=word::null)
Write sensitivities to file.
Definition: sensitivityVolBSplinesFIIncompressible.C:336
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:333
Foam::incompressible::FIBase::includeDistance_
bool includeDistance_
Include distance variation in sens computation.
Definition: FIBaseIncompressible.H:75
Foam::incompressible::sensitivityVolBSplinesFI::assembleSensitivities
virtual void assembleSensitivities()
Assemble sensitivities.
Definition: sensitivityVolBSplinesFIIncompressible.C:110
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:321
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:248
Foam::VectorSpace< Vector< scalar >, scalar, 3 >::zero
static const Vector< scalar > zero
Definition: VectorSpace.H:115
Foam::incompressible::shapeSensitivities::dnfdbMult_
autoPtr< boundaryVectorField > dnfdbMult_
Definition: shapeSensitivitiesIncompressible.H:66
Foam::incompressible::sensitivityVolBSplinesFI::derivativesFolder_
fileName derivativesFolder_
Definition: sensitivityVolBSplinesFIIncompressible.H:93
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:179