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-2021 PCOpt/NTUA
9 Copyright (C) 2013-2021 FOSS GP
10 Copyright (C) 2019 OpenCFD Ltd.
11-------------------------------------------------------------------------------
12License
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"
34
35// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36
37namespace Foam
38{
39
40namespace incompressible
41{
42
43// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44
47(
51);
52
53// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
54
56(
57 const fvMesh& mesh,
58 const dictionary& dict,
60)
61:
63 volBSplinesBase_
64 (
66 ),
67 flowSens_(0),
68 dSdbSens_(0),
69 dndbSens_(0),
70 dxdbDirectSens_(0),
71 dVdbSens_(0),
72 distanceSens_(0),
73 optionsSens_(0),
74 bcSens_(0),
75
76 derivativesFolder_("optimisation"/type() + "Derivatives")
77{
78 // No boundary field pointers need to be allocated
79
81 derivatives_ = scalarField(3*nCPs, Zero);
82 flowSens_ = vectorField(nCPs, Zero);
83 dSdbSens_ = vectorField(nCPs, Zero);
84 dndbSens_ = vectorField(nCPs, Zero);
86 dVdbSens_ = vectorField(nCPs, Zero);
89 bcSens_ = vectorField(nCPs, Zero);
90
91 // Create folder to store sensitivities
93}
94
95
96// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
97
99{
100 /*
101 addProfiling
102 (
103 sensitivityVolBSplinesFI,
104 "sensitivityVolBSplinesFI::assembleSensitivities"
105 );
106 */
107 read();
108
109 // Interpolation engine
111
112 // Adjoint to the eikonal equation
113 autoPtr<volTensorField> distanceSensPtr(nullptr);
115 {
116 // Solver equation
117 eikonalSolver_->solve();
118
119 // Allocate memory and compute grad(dxdb) multiplier
120 distanceSensPtr.reset
121 (
122 createZeroFieldPtr<tensor>
123 (
124 mesh_,
125 "distanceSensPtr",
126 dimensionSet(0, 2, -3, 0, 0, 0, 0)
127 )
128 );
129 distanceSensPtr() = eikonalSolver_->getFISensitivityTerm()().T();
130 }
131
132 // Integration
133 label passedCPs(0);
135 forAll(boxes, iNURB)
136 {
137 const label nb(boxes[iNURB].getControlPoints().size());
138 vectorField boxSensitivities(nb, Zero);
139
140 vectorField dxdbSens = boxes[iNURB].computeControlPointSensitivities
141 (
143 sensitivityPatchIDs_.toc()
144 );
145
146 vectorField bcSens = boxes[iNURB].computeControlPointSensitivities
147 (
148 bcDxDbMult_(),
149 sensitivityPatchIDs_.toc()
150 );
151
152 for (label cpI = 0; cpI < nb; cpI++)
153 {
154 label globalCP = passedCPs + cpI;
155
156 // Parameterization info
157 tmp<volTensorField> tvolDxDbI
158 (
159 volPointInter.interpolate(boxes[iNURB].getDxDb(cpI))
160 );
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 const tensorField& gradDxDbMultInt = gradDxDbMult_.primitiveField();
171 for (label idir = 0; idir < pTraits<vector>::nComponents; ++idir)
172 {
173 // Gradient of parameterization info
174 auto ttemp =
176 (
178 (
179 "dxdb",
180 mesh_.time().timeName(),
181 mesh_,
184 ),
185 mesh_,
187 );
188 volVectorField& temp = ttemp.ref();
189 unzipCol(volDxDbI, vector::components(idir), temp);
190
191 volTensorField gradDxDb(fvc::grad(temp));
192 // Volume integral terms
193 flowSens_[globalCP].component(idir) = gSum
194 (
195 (gradDxDbMultInt && gradDxDb.primitiveField())
196 *mesh_.V()
197 );
198
200 {
201 const tensorField& distSensInt =
202 distanceSensPtr().primitiveField();
203 distanceSens_[globalCP].component(idir) =
204 gSum
205 (
206 (distSensInt && gradDxDb.primitiveField())
207 *mesh_.V()
208 );
209 }
210 }
211
212 // Contribution from objective function term from
213 // delta( n dS ) / delta b and
214 // delta ( x ) / delta b
215 // for objectives directly depending on x
216 for (const label patchI : sensitivityPatchIDs_)
217 {
218 tensorField dSdb
219 (
220 boxes[iNURB].dndbBasedSensitivities(patchI, cpI)
221 );
222 dSdbSens_[globalCP] += gSum(dSfdbMult_()[patchI] & dSdb);
223 tensorField dndb
224 (
225 boxes[iNURB].dndbBasedSensitivities(patchI, cpI, false)
226 );
227 dndbSens_[globalCP] += gSum((dnfdbMult_()[patchI] & dndb));
228 }
229
230 // Contribution from delta (V) / delta b
231 // For objectives defined as volume integrals only
232 dVdbSens_[globalCP] +=
233 gSum
234 (
236 *fvc::div(T(volDxDbI))().primitiveField()
237 *mesh_.V()
238 );
239
240 // Terms from fvOptions
241 optionsSens_[globalCP] +=
242 gSum((optionsDxDbMult_ & volDxDbI.primitiveField())*mesh_.V());
243
244 // dxdbSens storage
245 dxdbDirectSens_[globalCP] = dxdbSens[cpI];
246
247 // bcSens storage
248 bcSens_[globalCP] = bcSens[cpI];
249
250 boxSensitivities[cpI] =
251 flowSens_[globalCP]
252 + dSdbSens_[globalCP]
253 + dndbSens_[globalCP]
254 + dVdbSens_[globalCP]
255 + distanceSens_[globalCP]
256 + dxdbDirectSens_[globalCP]
257 + optionsSens_[globalCP]
258 + bcSens_[globalCP];
259 }
260
261 // Zero sensitivities in non-active design variables
262 boxes[iNURB].boundControlPointMovement(boxSensitivities);
263
264 // Transfer sensitivities to global list
265 for (label cpI = 0; cpI < nb; cpI++)
266 {
267 label globalCP = passedCPs + cpI;
268 derivatives_[3*globalCP] = boxSensitivities[cpI].x();
269 derivatives_[3*globalCP + 1] = boxSensitivities[cpI].y();
270 derivatives_[3*globalCP + 2] = boxSensitivities[cpI].z();
271 }
272
273 // Increment number of passed sensitivities
274 passedCPs += nb;
275 }
276
277 // Zero non-active sensitivity components.
278 // For consistent output only, does not affect optimisation
287
288 //profiling::writeNow();
289}
290
291
293{
302
304}
305
306
308{
309 Info<< "Writing control point sensitivities to file" << endl;
310 if (Pstream::master())
311 {
312 OFstream derivFile
313 (
315 baseName + adjointVars_.solverName() + mesh_.time().timeName()
316 );
317 unsigned int widthDV
318 (
319 max(int(Foam::name(flowSens_.size()).size()), int(3))
320 );
321 unsigned int width = IOstream::defaultPrecision() + 7;
322 derivFile
323 << setw(widthDV) << "#cp" << " "
324 << setw(width) << "total::x" << " "
325 << setw(width) << "total::y" << " "
326 << setw(width) << "total::z" << " "
327 << setw(width) << "flow::x" << " "
328 << setw(width) << "flow::y" << " "
329 << setw(width) << "flow::z" << " "
330 << setw(width) << "dSdb::x" << " "
331 << setw(width) << "dSdb::y" << " "
332 << setw(width) << "dSdb::z" << " "
333 << setw(width) << "dndb::x" << " "
334 << setw(width) << "dndb::y" << " "
335 << setw(width) << "dndb::z" << " "
336 << setw(width) << "dxdbDirect::x" << " "
337 << setw(width) << "dxdbDirect::y" << " "
338 << setw(width) << "dxdbDirect::z" << " "
339 << setw(width) << "dVdb::x" << " "
340 << setw(width) << "dVdb::y" << " "
341 << setw(width) << "dVdb::z" << " "
342 << setw(width) << "distance::x" << " "
343 << setw(width) << "distance::y" << " "
344 << setw(width) << "distance::z" << " "
345 << setw(width) << "options::x" << " "
346 << setw(width) << "options::y" << " "
347 << setw(width) << "options::z" << " "
348 << setw(width) << "dvdb::x" << " "
349 << setw(width) << "dvdb::y" << " "
350 << setw(width) << "dvdb::z" << endl;
351
352 label passedCPs(0);
353 label lastActive(-1);
355 forAll(boxes, iNURB)
356 {
357 label nb = boxes[iNURB].getControlPoints().size();
358 const boolList& activeCPs = boxes[iNURB].getActiveCPs();
359 for (label iCP = 0; iCP < nb; iCP++)
360 {
361 if (activeCPs[iCP])
362 {
363 label globalCP = passedCPs + iCP;
364 if (globalCP!=lastActive + 1) derivFile << "\n";
365 lastActive = globalCP;
366
367 derivFile
368 << setw(widthDV) << globalCP << " "
369 << setw(width) << derivatives_[3*globalCP] << " "
370 << setw(width) << derivatives_[3*globalCP + 1] << " "
371 << setw(width) << derivatives_[3*globalCP + 2] << " "
372 << setw(width) << flowSens_[globalCP].x() << " "
373 << setw(width) << flowSens_[globalCP].y() << " "
374 << setw(width) << flowSens_[globalCP].z() << " "
375 << setw(width) << dSdbSens_[globalCP].x() << " "
376 << setw(width) << dSdbSens_[globalCP].y() << " "
377 << setw(width) << dSdbSens_[globalCP].z() << " "
378 << setw(width) << dndbSens_[globalCP].x() << " "
379 << setw(width) << dndbSens_[globalCP].y() << " "
380 << setw(width) << dndbSens_[globalCP].z() << " "
381 << setw(width) << dxdbDirectSens_[globalCP].x() << " "
382 << setw(width) << dxdbDirectSens_[globalCP].y() << " "
383 << setw(width) << dxdbDirectSens_[globalCP].z() << " "
384 << setw(width) << dVdbSens_[globalCP].x() << " "
385 << setw(width) << dVdbSens_[globalCP].y() << " "
386 << setw(width) << dVdbSens_[globalCP].z() << " "
387 << setw(width) << distanceSens_[globalCP].x() << " "
388 << setw(width) << distanceSens_[globalCP].y() << " "
389 << setw(width) << distanceSens_[globalCP].z() << " "
390 << setw(width) << optionsSens_[globalCP].x() << " "
391 << setw(width) << optionsSens_[globalCP].y() << " "
392 << setw(width) << optionsSens_[globalCP].z() << " "
393 << setw(width) << bcSens_[globalCP].x() << " "
394 << setw(width) << bcSens_[globalCP].y() << " "
395 << setw(width) << bcSens_[globalCP].z() << endl;
396 }
397 }
398 passedCPs += nb;
399 }
400 }
401}
402
403
404// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
405
406} // End namespace incompressible
407} // End namespace Foam
408
409// ************************************************************************* //
Istream and Ostream manipulators taking arguments.
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
tmp< Field< cmptType > > component(const direction) const
Return a component field of the field.
Definition: Field.C:545
const Internal::FieldType & primitiveField() const
Return a const-reference to the internal field.
Internal & ref(const bool updateAccessTime=true)
Return a reference to the dimensioned internal field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
static unsigned int defaultPrecision() noexcept
Return the default precision.
Definition: IOstream.H:342
Output to file stream, using an OSstream.
Definition: OFstream.H:57
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
label size() const noexcept
The number of elements in the list.
Definition: UPtrListI.H:106
components
Component labeling enumeration.
Definition: Vector.H:81
Base class for adjoint solvers.
Definition: adjointSolver.H:60
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
void reset(autoPtr< T > &&other) noexcept
Delete managed object and set to new given pointer.
Definition: autoPtrI.H:117
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:109
virtual bool write()
Write the output fields.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
Base class for incompressibleAdjoint solvers.
Base class for Field Integral-based sensitivity derivatives.
virtual void clearSensitivities()
Zero sensitivity fields and their constituents.
bool includeDistance_
Include distance variation in sens computation.
volTensorField gradDxDbMult_
grad(dx/db) multiplier
autoPtr< adjointEikonalSolver > eikonalSolver_
Adjoint eikonal equation solver.
void read()
Read options and update solver pointers if necessary.
scalarField divDxDbMult_
div(dx/db) multiplier
vectorField optionsDxDbMult_
dx/db multiplier coming from fvOptions
Abstract base class for adjoint-based sensitivities in incompressible flows.
Calculation of adjoint based sensitivities at vol B-Splines control points using the FI approach.
virtual void clearSensitivities()
Zero sensitivity fields and their constituents.
vectorField dVdbSens_
Term depending on delta(V)/delta b.
vectorField dSdbSens_
Term depending on delta(n dS)/delta b.
vectorField bcSens_
Term depending on the differentiation of boundary conditions.
vectorField dndbSens_
Term depending on delta(n)/delta b.
vectorField distanceSens_
Term depending on distance differentiation.
volBSplinesBase & volBSplinesBase_
Reference to underlaying volumetric B-Splines morpher.
autoPtr< boundaryVectorField > dSfdbMult_
Fields related to direct sensitivities.
void interpolate(const GeometricField< Type, pointPatchField, pointMesh > &, GeometricField< Type, fvPatchField, volMesh > &) const
Interpolate from pointField to volField.
const fvMesh & mesh_
Definition: sensitivity.H:69
splitCell * master() const
Definition: splitCell.H:113
A class for managing temporary objects.
Definition: tmp.H:65
const word & solverName() const
Return solver name.
Definition: variablesSet.C:84
Class constructing a number of volumetric B-Splines boxes, read from dynamicMeshDict....
PtrList< NURBS3DVolume > & boxesRef()
Get non-const reference to the vol. B-splines boxes.
label getTotalControlPointsNumber() const
Get cumulative number of control points from all boxes.
void boundControlPointMovement(vectorField &controlPointsMovement) const
Bound control points movement.
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
const volScalarField & T
dynamicFvMesh & mesh
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:49
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Type gSum(const FieldField< Field, Type > &f)
const dimensionSet dimless
Dimensionless.
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:515
messageStream Info
Information stream (stdout output on master, null elsewhere)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
void unzipCol(const FieldField< Field, SymmTensor< Cmpt > > &input, const direction idx, FieldField< Field, Vector< Cmpt > > &result)
Extract a symmTensor field field column (x,y,z) == (0,1,2)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Omanip< int > setw(const int i)
Definition: IOmanip.H:199
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
Field< vector > vectorField
Specialisation of Field<T> for vector.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
A non-counting (dummy) refCount.
Definition: refCount.H:59