Loading...
Searching...
No Matches
shapeDesignVariables.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-2023 PCOpt/NTUA
9 Copyright (C) 2013-2023 FOSS GP
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
30#include "cellQuality.H"
31#include "createZeroField.H"
33#include "volFieldsFwd.H"
35#include "IOmanip.H"
37// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38
39namespace Foam
40{
44 (
47 designVariables
48 );
49}
50
51
52// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
55{
56 return size();
57}
58
59
61{
63}
64
65
68(
69 const label patchI,
70 const label varID
71) const
72{
73 const dictionary dxdbDict = dict_.subOrEmptyDict("dxdbSolver");
74 const label iters = dxdbDict.getOrDefault<label>("iters", 1000);
75 const scalar tolerance =
76 dxdbDict.getOrDefault<scalar>("tolerance", 1.e-07);
78 (
80 (
82 (
83 mesh_,
84 "m",
86 )
87 )
88 );
89 volVectorField& m = tm.ref();
90
91 // Solve for dxdb
92 //~~~~~~~~~~~~~~~~
93 m.boundaryFieldRef()[patchI] == dxdbFace(patchI, varID);
94
95 // Iterate the direct differentiation of the grid displacement equation
96 for (label iter = 0; iter < iters; ++iter)
97 {
98 Info<< "Mesh Movement Propagation for varID" << varID
99 << ", Iteration : "<< iter << endl;
100
101 fvVectorMatrix mEqn
102 (
104 );
105
106 scalar residual = mag(mEqn.solve().initialResidual());
107
109 << "Max dxdb " << gMax(mag(m)()) << endl;
110
111 mesh_.time().printExecutionTime(Info);
112
113 // Check convergence
114 if (residual < tolerance)
115 {
116 Info<< "\n***Reached dxdb convergence limit, iteration " << iter
117 << "***\n\n";
118 break;
120 }
121
122 return tm;
123}
124
125
127{
128 if (dxdbVolSens_.empty())
129 {
130 dxdbVolSens_.setSize(sensSize(), Zero);
131 dxdbSurfSens_.setSize(sensSize(), Zero);
132 dSdbSens_.setSize(sensSize(), Zero);
133 dndbSens_.setSize(sensSize(), Zero);
134 dxdbDirectSens_.setSize(sensSize(), Zero);
135 dVdbSens_.setSize(sensSize(), Zero);
137 optionsSens_.setSize(sensSize(), Zero);
138 bcSens_.setSize(sensSize(), Zero);
139 }
140}
141
142
144{
145 dxdbVolSens_ = Zero;
146 dxdbSurfSens_ = Zero;
147 dSdbSens_ = Zero;
148 dndbSens_ = Zero;
149 dxdbDirectSens_ = Zero;
150 dVdbSens_ = Zero;
151 distanceSens_ = Zero;
153 bcSens_ = Zero;
154}
155
156
157// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
158
159Foam::shapeDesignVariables::shapeDesignVariables
160(
161 fvMesh& mesh,
162 const dictionary& dict
163)
164:
166 parametertisedPatches_
167 (
168 mesh_.boundaryMesh().patchSet(dict.get<wordRes>("patches"))
169 ),
170 displMethodPtr_
171 (
172 displacementMethod::New(mesh_, parametertisedPatches_.toc())
173 ),
174 pointsInit_(nullptr),
175 writeEachMesh_(dict.getOrDefault<bool>("writeEachMesh", true)),
176 dxdbVolSens_(),
177 dxdbSurfSens_(),
178 dSdbSens_(),
179 dndbSens_(),
180 dxdbDirectSens_(),
181 dVdbSens_(),
182 distanceSens_(),
183 optionsSens_(),
184 bcSens_(),
185 derivativesFolder_
186 (
187 word("optimisation")/word("derivatives")
188 /word(mesh.name() == polyMesh::defaultRegion ? word() : mesh.name())
189 )
190{
192 {
193 FatalErrorInFunction
194 << "None of the provided parameterised patches is valid"
195 << endl
196 << exit(FatalError);
199}
200
201
202// * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
203
205(
206 fvMesh& mesh,
207 const dictionary& dict
208)
209{
210 const word modelType(dict.get<word>("shapeType"));
211
212 Info<< "shapeDesignVariables type : " << modelType << endl;
213
214 auto* ctorPtr = dictionaryConstructorTable(modelType);
215
216 if (!ctorPtr)
217 {
219 (
220 "shapeType",
221 modelType,
222 *dictionaryConstructorTablePtr_
223 ) << exit(FatalError);
224 }
226 return autoPtr<shapeDesignVariables>(ctorPtr(mesh, dict));
227}
228
229
230// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * //
231
233{
235 {
236 parametertisedPatches_ =
237 mesh_.boundaryMesh().patchSet(dict.get<wordRes>("patches"));
238 displMethodPtr_->setPatchIDs(parametertisedPatches_.toc());
239
240 writeEachMesh_ =
241 dict.getOrDefault<bool>("writeEachMesh", true);
242
243 return true;
244 }
245
246 return false;
247}
248
249
251{
253
254 if (!pointsInit_)
256 pointsInit_.reset(new pointField(mesh_.nPoints(), Zero));
257 }
258 pointsInit_() = mesh_.points();
259}
260
261
263{
265 mesh_.movePoints(pointsInit_());
266}
267
268
270{
271 // Move mesh
272 displMethodPtr_->update();
273
274 if (writeEachMesh_)
275 {
276 Info<< " Writing new mesh points for mesh region "
277 << mesh_.name() << endl;
279 (
281 (
282 "points",
283 mesh_.pointsInstance(),
285 mesh_,
289 ),
290 mesh_.points()
291 );
292 points.write();
294
295 // Check mesh quality
296 mesh_.checkMesh(true);
297}
298
299
301(
302 adjointSensitivity& adjointSens
303)
304{
305 // Return field
306 tmp<scalarField> tsens(tmp<scalarField>::New(sensSize(), Zero));
307 scalarField& sens = tsens.ref();
308
309 // Reset sensitivity components to zero
310 allocateSensFields();
311 zeroSensFields();
312
313 // Grab multipliers from the adjoint sensitivities
314 const autoPtr<volTensorField>& gradDxDbMult = adjointSens.gradDxDbMult();
315 const autoPtr<scalarField>& divDxDbMult = adjointSens.divDxDbMult();
316 const autoPtr<boundaryVectorField>& dxdbMult = adjointSens.dxdbMult();
317 const autoPtr<boundaryVectorField>& dSdbMult = adjointSens.dSfdbMult();
318 const autoPtr<boundaryVectorField>& dndbMult = adjointSens.dnfdbMult();
319 const autoPtr<boundaryVectorField>& dxdbDirectMult =
320 adjointSens.dxdbDirectMult();
321 const autoPtr<boundaryVectorField>& bcDxDbmult = adjointSens.bcDxDbMult();
322 const autoPtr<vectorField>& optionsDxDbMult = adjointSens.optionsDxDbMult();
323 const volScalarField::Internal& V = mesh_.V();
324 autoPtr<adjointEikonalSolver>& eikonalSolver =
325 adjointSens.getAdjointEikonalSolver();
326
327 autoPtr<volTensorField> distanceSens(nullptr);
328 if (adjointSens.includeDistance())
329 {
330 distanceSens.reset
331 (
332 new volTensorField(eikonalSolver->getFISensitivityTerm())
333 );
334 }
335
336 // Loop over active design variables only
337 for (const label varI : activeSensitivities())
338 {
339 // FI approach, integrate terms including variations of the grid
340 // sensitivities
341 if (adjointSens.computeDxDbInternalField())
342 {
343 // Parameterization info
344 tmp<volVectorField> tvolDxDbI = dCdb(varI);
345 const volVectorField& volDxDbI = tvolDxDbI();
346 tmp<volTensorField> gradDxDb = fvc::grad(volDxDbI);
347
348 // Contributions from the adjoint-related part
349 dxdbVolSens_[varI] = gSum((gradDxDbMult() && gradDxDb())*V);
350
351 // Contributions from the distance related part
352 if (adjointSens.includeDistance())
353 {
354 distanceSens_[varI] = gSum((distanceSens() && gradDxDb)*V);
355 }
356 // Contributions from the multiplier of divDxDb
357 if (divDxDbMult)
358 {
359 dVdbSens_[varI] +=
360 gSum(divDxDbMult()*fvc::div(volDxDbI)().primitiveField()*V);
361 }
362
363 // Contributions from fvOptions
364 optionsSens_[varI] +=
365 gSum((optionsDxDbMult() & volDxDbI.primitiveField())*V);
366 }
367
368 // Contribution from boundary terms
369 // Most of them (with the expection of dxdbMult) exist in both the
370 // FI and E-SI approaches
371 for (const label patchI : parametertisedPatches_)
372 {
373 if (dSdbMult)
374 {
375 tmp<vectorField> pdSdb = dSdb(patchI, varI);
376 dSdbSens_[varI] += gSum(dSdbMult()[patchI] & pdSdb);
377 }
378
379 if (dndbMult)
380 {
381 tmp<vectorField> pdndb = dndb(patchI, varI);
382 dndbSens_[varI] += gSum((dndbMult()[patchI] & pdndb));
383 }
384
385 tmp<vectorField> pdxdb = dxdbFace(patchI, varI);
386 // Main contribution in the E-SI approach
387 if (dxdbMult)
388 {
389 dxdbSurfSens_[varI] += gSum(dxdbMult()[patchI] & pdxdb());
390 }
391 if (dxdbDirectMult)
392 {
393 dxdbDirectSens_[varI] +=
394 gSum((dxdbDirectMult()[patchI] & pdxdb()));
395 }
396 if (bcDxDbmult)
397 {
398 bcSens_[varI] += gSum((bcDxDbmult()[patchI] & pdxdb()));
399 }
400 }
401 }
402
403 sens =
404 dxdbVolSens_ + dxdbSurfSens_ + dSdbSens_ + dndbSens_ + dxdbDirectSens_
405 + dVdbSens_ + distanceSens_ + optionsSens_ + bcSens_;
407 writeSensitivities(sens, adjointSens);
408
409 return tsens;
410}
411
412
414(
415 const scalarField& sens,
416 const adjointSensitivity& adjointSens
417)
418{
419 OFstream derivFile
420 (
421 derivativesFolder_/
422 type() + adjointSens.getAdjointSolver().solverName()
423 + adjointSens.getSuffix() + mesh_.time().timeName()
424 );
425 unsigned int widthDV = max(int(name(dxdbVolSens_.size()).size()), int(6));
426 unsigned int width = IOstream::defaultPrecision() + 7;
427 derivFile
428 << setw(widthDV) << "#varID" << " "
429 << setw(width) << "total"<< " "
430 << setw(width) << "dxdbVol" << " "
431 << setw(width) << "dxdbSurf" << " "
432 << setw(width) << "dSdb" << " "
433 << setw(width) << "dndb" << " "
434 << setw(width) << "dxdbDirect" << " "
435 << setw(width) << "dVdb" << " "
436 << setw(width) << "distance" << " "
437 << setw(width) << "options" << " "
438 << setw(width) << "dvdb" << endl;
439
440 for (const label varI : activeSensitivities())
441 {
442 derivFile
443 << setw(widthDV) << varI << " "
444 << setw(width) << sens[varI] << " "
445 << setw(width) << dxdbVolSens_[varI] << " "
446 << setw(width) << dxdbSurfSens_[varI] << " "
447 << setw(width) << dSdbSens_[varI] << " "
448 << setw(width) << dndbSens_[varI] << " "
449 << setw(width) << dxdbDirectSens_[varI] << " "
450 << setw(width) << dVdbSens_[varI] << " "
451 << setw(width) << distanceSens_[varI] << " "
452 << setw(width) << optionsSens_[varI] << " "
453 << setw(width) << bcSens_[varI] << " "
454 << endl;
455 }
456}
457
458
460(
461 const label varID
462) const
463{
464 // Deliberately returning a zero-sized field
465 return tmp<vectorField>::New(0);
466}
467
468
470(
471 const label patchI,
472 const label varID
473) const
474{
476 return nullptr;
477}
478
479
481(
482 const label patchI,
483 const label varID
484) const
485{
487 return nullptr;
488}
489
490
492(
493 const label patchI,
494 const label varID
495) const
498 return nullptr;
499}
500
501
503Foam::shapeDesignVariables::dCdb(const label varID) const
504{
506 return nullptr;
507}
508
509
510// ************************************************************************* //
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.
DimensionedField< scalar, volMesh > Internal
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
const Internal::FieldType & primitiveField() const noexcept
Return a const-reference to the internal field values.
label size() const noexcept
The number of elements in table.
Definition HashTable.H:358
@ NO_REGISTER
Do not request registration (bool: false).
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
static unsigned int defaultPrecision() noexcept
Return the default precision.
Definition IOstream.H:437
Output to file stream as an OSstream, normally using std::ofstream for the actual output.
Definition OFstream.H:75
const Type & initialResidual() const noexcept
Return initial residual.
bool get(const label i) const
Definition UList.H:868
label size() const noexcept
Definition UList.H:706
Abstract base class for adjoint-based sensitivities.
autoPtr< adjointEikonalSolver > & getAdjointEikonalSolver()
Return the adjoint eikonal solver.
const autoPtr< volTensorField > & gradDxDbMult() const
const autoPtr< vectorField > & optionsDxDbMult() const
bool includeDistance() const
Should the adjoint eikonal PDE should be solved.
const autoPtr< boundaryVectorField > & dxdbDirectMult() const
const adjointSolver & getAdjointSolver() const
Const access to adjoint solver.
const autoPtr< boundaryVectorField > & dnfdbMult() const
virtual bool computeDxDbInternalField() const
Should the parameterization compute the internalField of dxdb.
const autoPtr< boundaryVectorField > & dxdbMult() const
const autoPtr< scalarField > & divDxDbMult() const
const autoPtr< boundaryVectorField > & dSfdbMult() const
const word & getSuffix() const
Get suffix.
const autoPtr< boundaryVectorField > & bcDxDbMult() const
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
void reset(T *p=nullptr) noexcept
Delete managed object and set to new given pointer.
Definition autoPtrI.H:37
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
Abstract base class for defining design variables.
virtual void resetDesignVariables()
Reset to the starting point of line search.
virtual bool readDict(const dictionary &dict)
Read dictionary if changed.
labelList activeDesignVariables_
Which of the design variables will be updated.
virtual void storeDesignVariables()
Store design variables, as the starting point for line search.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
dictionary subOrEmptyDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX, const bool mandatory=false) const
Find and return a sub-dictionary as a copy, otherwise return an empty dictionary.
Definition dictionary.C:521
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Abstract base class for displacement methods, which are a set or wrapper classes allowing to change t...
SolverPerformance< Type > solve(const dictionary &)
Solve returning the solution statistics.
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh").
Definition polyMesh.H:411
Abstract base class for defining design variables for shape optimisation.
bool writeEachMesh_
Write the mesh points irrespective of whether this is a write time.
scalarField dSdbSens_
Term depending on delta(n dS)/delta b.
virtual tmp< vectorField > dSdb(const label patchI, const label varID) const
Get dSdb for a given design variable and patch.
virtual void resetDesignVariables()
Reset to starting point of line search.
scalarField dVdbSens_
Term depending on delta(V)/delta b.
void zeroSensFields()
Zero the fields assosiated with the computation of sensitivities.
scalarField optionsSens_
Term depending on fvOptions.
fileName derivativesFolder_
Name of the sensitivity derivatives folder.
virtual void moveMesh()
Move mesh based on displacementMethod.
virtual tmp< volVectorField > solveMeshMovementEqn(const label patchI, const label varID) const
Compute dxdb at the mesh cell centers by solving a Laplace PDE.
virtual bool readDict(const dictionary &dict)
Read dictionary if changed.
autoPtr< displacementMethod > displMethodPtr_
Mesh movement mechanism.
void allocateSensFields()
Allocate the fields assosiated with the computation of sensitivities.
virtual void writeSensitivities(const scalarField &sens, const adjointSensitivity &adjointSens)
Write final sensitivity derivatives to files.
scalarField dxdbDirectSens_
Term depending on delta(x)/delta b for objectives that directly depend on x.
autoPtr< pointField > pointsInit_
Store old points. Useful for line search.
static autoPtr< shapeDesignVariables > New(fvMesh &mesh, const dictionary &dict)
Construct and return the selected shapeDesignVariables.
scalarField dxdbSurfSens_
Flow related term.
scalarField distanceSens_
Term depending on distance differentiation.
labelHashSet parametertisedPatches_
Patches to be moved by the design variables.
scalarField dxdbVolSens_
Flow related term.
scalarField bcSens_
Term depending on the differenation of boundary conditions.
virtual void storeDesignVariables()
Store design variables, as the starting point for line search.
virtual tmp< vectorField > dxdbVol(const label varID) const
Get dxdb for all mesh points.
virtual const labelList & activeSensitivities() const
Active variables for which to compute sensitivities.
virtual label sensSize() const
Size of the sensitivity derivatives.
virtual tmp< vectorField > dxdbFace(const label patchI, const label varID) const
Get dxdb for a given design variable and patch.
virtual tmp< vectorField > dndb(const label patchI, const label varID) const
Get dndb for a given design variable and patch.
virtual tmp< volVectorField > dCdb(const label varID) const
Get dCdb for a given design variable.
virtual tmp< scalarField > assembleSensitivities(adjointSensitivity &adjointSens)
Add part of sensitivity derivatives related to geometry variations.
scalarField dndbSens_
Term depending on delta(n)/delta b.
const word & solverName() const
Return the solver name.
Definition solverI.H:30
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition tmpI.H:235
static tmp< volVectorField > autoCreateMeshMovementField(const fvMesh &mesh, const word &name, const dimensionSet &dims)
Auto create variable for mesh movement.
A List of wordRe with additional matching capabilities.
Definition wordRes.H:56
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
dynamicFvMesh & mesh
#define FatalErrorInLookup(lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalError.
Definition error.H:607
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition error.H:688
auto & name
const pointField & points
#define DebugInfo
Report an information message using Foam::Info.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcGrad.C:47
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcDiv.C:42
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
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:40
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
Type gSum(const FieldField< Field, Type > &f)
GeometricField< vector, fvPatchField, volMesh > volVectorField
List< label > labelList
A List of labels.
Definition List.H:62
bool mkDir(const fileName &pathName, mode_t mode=0777)
Make a directory and return an error if it could not be created.
Definition POSIX.C:616
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
messageStream Info
Information stream (stdout output on master, null elsewhere).
vectorIOField pointIOField
pointIOField is a vectorIOField.
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition POSIX.C:801
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:519
GeometricField< tensor, fvPatchField, volMesh > volTensorField
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
fvMatrix< vector > fvVectorMatrix
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
vectorField pointField
pointField is a vectorField.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
Type gMax(const FieldField< Field, Type > &f)
loopControl iters(runTime, aMesh.solutionDict(), "solution")
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.
dictionary dict
Forwards and collection of common volume field types.