Loading...
Searching...
No Matches
optimisationManager.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 Copyright (C) 2019-2021 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
31#include "dictionary.H"
33
34// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35
36namespace Foam
37{
38 defineTypeNameAndDebug(optimisationManager, 0);
39 defineRunTimeSelectionTable(optimisationManager, dictionary);
40}
41
42
43// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
44
46(
48 const label startTimeIndex,
49 const scalar endTime
50)
51{
52 // Does nothing in base
53}
54
55
57{
58 // Compute direction of update
59 tmp<scalarField> tdirection = dvUpdate_->computeDirection();
60 const scalarField& direction = tdirection.ref();
61
62 // Grab reference to line search
63 autoPtr<lineSearch>& lineSrch = dvUpdate_->getLineSearch();
64
65 // Store starting point
66 dvUpdate_->getDesignVariables()->storeDesignVariables();
67
68 // Compute merit function before update
69 scalar meritFunction = dvUpdate_->computeMeritFunction();
70 lineSrch->setOldMeritValue(meritFunction);
71 dvUpdate_->setOldObjectiveValue();
72
73 // Get merit function derivative
74 scalar dirDerivative =
75 dvUpdate_->meritFunctionDirectionalDerivative();
76 lineSrch->setDeriv(dirDerivative);
77 lineSrch->setDirection(direction);
78
79 // Reset initial step.
80 // Might be interpolated from previous optimisation cycles
81 lineSrch->reset();
82
83 // Perform line search
84 while (lineSrch->loop())
85 {
86 Info<< "\n- - - - - - - - - - - - - - -" << endl;
87 Info<< "Line search iteration " << lineSrch->innerIter() << endl;
88 Info<< "- - - - - - - - - - - - - - -\n" << endl;
89
90 // Update design variables. Multiplication with line search step
91 // happens inside the update(direction) function
92 moveDesignVariables(direction);
93
94 const dimensionedScalar startTime = mesh_.time();
95 const label startTimeIndex = mesh_.time().timeIndex();
96 const scalar primalEndTime = mesh_.time().endTime().value();
97 // Solve all primal equations
98 solvePrimalEquations();
99
100 // Compute and set new merit function
101 meritFunction = dvUpdate_->computeMeritFunction();
102 lineSrch->setNewMeritValue(meritFunction);
103
104 if (lineSrch->computeGradient())
105 {
106 // Reset adjoint sensitivities in all adjoint solver managers
107 clearSensitivities();
108
109 // Solve all adjoint equations
110 solveAdjointEquations();
111
112 // Update objective and gradient information known by updateMethod
113 dvUpdate_->updateGradientsAndValues();
114
115 // Update the directional derivative
116 dirDerivative =
117 dvUpdate_->meritFunctionDirectionalDerivative();
118 lineSrch->setNewDeriv(dirDerivative);
119 }
120
121 if (lineSrch->converged())
122 {
123 // If line search criteria have been met, proceed
124 Info<< "Line search converged in " << lineSrch->innerIter()
125 << " iterations." << endl;
126 scalarField scaledCorrection(lineSrch->step()*direction);
127 dvUpdate_->postUpdate(scaledCorrection);
128 break;
129 }
130 else
131 {
132 // If maximum number of iteration has been reached, continue
133 if (lineSrch->innerIter() == lineSrch->maxIters())
134 {
135 Info<< "Line search reached max. number of iterations.\n"
136 << "Proceeding to the next optimisation cycle" << endl;
137 scalarField scaledCorrection(lineSrch->step()*direction);
138 dvUpdate_->postUpdate(scaledCorrection);
139 }
140 // Reset to initial design variables and update step
141 else
142 {
143 //- Reset time if necessary
144 this->resetTime(startTime, startTimeIndex, primalEndTime);
145 dvUpdate_->getDesignVariables()->resetDesignVariables();
146 lineSrch->updateStep();
147 }
148 }
149 }
150
151 // If line search did not need to get the new gradient, do it now in order
152 // to have it for the next optimisation cycle
153 if (!lineSrch->computeGradient())
154 {
155 // Reset adjoint sensitivities in all adjoint solver managers
156 clearSensitivities();
157
158 // Solve all adjoint equations
160 }
161}
162
163
165{
166 // Update design variables
167 if (shouldUpdateDesignVariables_)
168 {
169 moveDesignVariables();
170 }
171
172 // Solve primal equations
173 solvePrimalEquations();
174
175 // Check the convergence criteria of the optimisation loop
176 dvUpdate_->checkConvergence();
177
178 // Reset adjoint sensitivities in all adjoint solver managers
179 clearSensitivities();
180
181 // Solve all adjoint equations
183}
184
185
187{
188 // Update design variables
189 dvUpdate_->update();
190}
191
192
194(
196)
197{
198 // Update design variables
199 dvUpdate_->update(direction);
200}
201
202
204{
205 dictionary& primalSolversDict = subDict("primalSolvers");
206 const wordList& primalSolverNames = primalSolversDict.toc();
207
208 // Construct primal solvers
209 primalSolvers_.setSize(primalSolverNames.size());
210 forAll(primalSolvers_, solveri)
211 {
212 dictionary& solverDict =
213 primalSolversDict.subDict(primalSolverNames[solveri]);
214 if (primalSolvers_.size() > 1)
215 {
216 solverDict.add<bool>("useSolverNameForFields", true);
217 }
218 primalSolvers_.set
219 (
220 solveri,
222 (
223 mesh_,
224 managerType_,
225 solverDict,
226 primalSolverNames[solveri]
227 )
228 );
229 }
230
231 // Construct adjointSolverManagers
232 const dictionary& adjointManagersDict = subDict("adjointManagers");
233 const wordList adjointManagerNames = adjointManagersDict.toc();
234 adjointSolverManagers_.setSize(adjointManagerNames.size());
235
236 // Determine the number of adjoint solvers which are not null (i.e. need to
237 // allocate adjoint fields). Used to determine whether the adjoint field
238 // names should be appended by the solver name
239 label nNotNullAdjointSolvers(0);
240 for (const word& adjManager : adjointManagerNames)
241 {
242 const dictionary& adjSolversDict =
243 adjointManagersDict.subDict(adjManager).subDict("adjointSolvers");
244 const wordList adjointSolverNames = adjSolversDict.toc();
245 for (const word& adjSolver : adjointSolverNames)
246 {
247 if (adjSolversDict.subDict(adjSolver).get<word>("type") != "null")
248 {
249 ++nNotNullAdjointSolvers;
250 }
251 }
252 }
253 Info<< "Found "
254 << nNotNullAdjointSolvers
255 << " adjoint solvers that allocate fields"
256 << endl;
257 bool overrideUseSolverName(nNotNullAdjointSolvers > 1);
258
259 forAll(adjointSolverManagers_, manageri)
260 {
261 adjointSolverManagers_.set
262 (
263 manageri,
265 (
266 mesh_,
267 designVars_,
268 managerType_,
269 adjointManagersDict.subDict(adjointManagerNames[manageri]),
270 overrideUseSolverName
271 )
272 );
273 }
274
275 // Sanity checks on the naming convention
276 if (primalSolvers_.size() > 1)
277 {
278 for (const primalSolver& solveri : primalSolvers_)
279 {
280 if (!solveri.useSolverNameForFields())
281 {
283 << "Multiple primal solvers are present but "
284 << "useSolverNameForFields is set to false in "
285 << "primal solver " << solveri.solverName() << nl
286 << "This is considered fatal."
287 << exit(FatalError);
288 }
289 }
290 }
291
292 if (nNotNullAdjointSolvers > 1)
293 {
294 for (const adjointSolverManager& amI : adjointSolverManagers_)
295 {
296 const PtrList<adjointSolver>& adjointSolvers = amI.adjointSolvers();
297 for (const adjointSolver& asI : adjointSolvers)
298 {
299 if (!asI.useSolverNameForFields())
300 {
302 << "Multiple adjoint solvers are present but "
303 << "useSolverNameForFields is set to false in "
304 << "adjoint solver " << asI.solverName() << nl
305 << "This is considered fatal."
306 << exit(FatalError);
307 }
308 }
309 }
310 }
311 if (designVars_)
312 {
313 designVars_().addFvOptions(primalSolvers_, adjointSolverManagers_);
314 }
316
317
318// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
319
320Foam::optimisationManager::optimisationManager(fvMesh& mesh)
321:
323 (
325 (
326 "optimisationDict",
327 mesh.time().system(),
328 mesh,
329 IOobject::READ_MODIFIED,
330 IOobject::NO_WRITE,
331 IOobject::REGISTER
332 )
333 ),
334 mesh_(mesh),
335 time_(const_cast<Time&>(mesh.time())),
336 designVars_(nullptr),
337 primalSolvers_(),
338 adjointSolverManagers_(),
339 managerType_(get<word>("optimisationManager")),
340 dvUpdate_(nullptr),
341 shouldUpdateDesignVariables_(true)
342{
343 // The "designVariables" sub-dictionary is optional
344 const dictionary odict(this->subOrEmptyDict("optimisation"));
345 const dictionary* designVarsDictPtr = odict.findDict("designVariables");
346
347 if (designVarsDictPtr)
348 {
349 designVars_ = designVariables::New(mesh_, *designVarsDictPtr);
350 }
352
353
354// * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * //
355
357(
358 fvMesh& mesh
359)
360{
361 const IOdictionary dict
362 (
364 (
365 "optimisationDict",
366 mesh.time().system(),
367 mesh,
371 )
372 );
373
374 const word modelType(dict.get<word>("optimisationManager"));
375
376 Info<< "optimisationManager type : " << modelType << endl;
377
378 auto* ctorPtr = dictionaryConstructorTable(modelType);
379
380 if (!ctorPtr)
381 {
383 (
384 dict,
385 "optimisationManager",
386 modelType,
387 *dictionaryConstructorTablePtr_
388 ) << exit(FatalIOError);
389 }
390
391 return autoPtr<optimisationManager>(ctorPtr(mesh));
393
394
395// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * //
396
398{
399 if (regIOobject::read())
400 {
401 // Note: Only changing existing solvers - not adding any new
402 const dictionary& primalSolversDict = subDict("primalSolvers");
403 for (primalSolver& sol : primalSolvers_)
404 {
405 sol.readDict(primalSolversDict.subDict(sol.solverName()));
406 }
407
408 const dictionary& adjointManagersDict = subDict("adjointManagers");
409 for (adjointSolverManager& man : adjointSolverManagers_)
410 {
411 man.readDict(adjointManagersDict.subDict(man.managerName()));
412 }
413
414 if (designVars_)
415 {
416 designVars_->readDict
417 (subDict("optimisation").subDict("designVariables"));
418 }
419
420 return true;
421 }
423 return false;
424}
425
426
428{
429 // Update design variables using either a line-search scheme or
430 // a fixed-step update
431 if (dvUpdate_->getLineSearch())
432 {
433 lineSearchUpdate();
434 }
435 else
436 {
438 }
439}
440
441
443{
444 // Solve all primal equations
445 forAll(primalSolvers_, psI)
446 {
447 primalSolvers_[psI].solve();
448 }
449}
450
451
453{
454 // Solve all adjoint solver equations
455 forAll(adjointSolverManagers_, amI)
456 {
457 adjointSolverManagers_[amI].solveAdjointEquations();
458 }
459}
460
461
463{
464 // Compute senstivities from all adjoint solvers
465 forAll(adjointSolverManagers_, amI)
466 {
467 adjointSolverManagers_[amI].computeAllSensitivities();
468 }
469}
470
471
473{
474 for (adjointSolverManager& adjSolvManager : adjointSolverManagers_)
475 {
476 adjSolvManager.clearSensitivities();
477 }
478}
479
480
482{
483 forAll(adjointSolverManagers_, amI)
484 {
485 PtrList<adjointSolver>& adjointSolvers =
486 adjointSolverManagers_[amI].adjointSolvers();
487
488 forAll(adjointSolvers, asI)
489 {
490 adjointSolvers[asI].updatePrimalBasedQuantities();
491 }
492 }
493}
494
495
496// ************************************************************************* //
Foam::label startTime
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
IOdictionary(const IOobject &io, const dictionary *fallback=nullptr)
Construct given an IOobject and optional fallback dictionary content.
@ NO_REGISTER
Do not request registration (bool: false).
@ REGISTER
Request registration (bool: true).
@ MUST_READ
Reading required.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const Time & time() const noexcept
Return Time associated with the objectRegistry.
Definition IOobject.C:456
void setSize(label n)
Alias for resize().
Definition List.H:536
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
Class for managing adjoint solvers, which may be more than one per operating point.
Base class for adjoint solvers.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
static autoPtr< designVariables > New(fvMesh &mesh, const dictionary &dict)
Return a reference to the selected design variables.
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 get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
const dictionary * findDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary pointer if present (and it is a dictionary) otherwise return nullptr...
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition dictionary.C:441
entry * add(entry *entryPtr, bool mergeEntry=false)
Add a new entry.
Definition dictionary.C:625
wordList toc() const
Return the table of contents.
Definition dictionary.C:587
entry * set(entry *entryPtr)
Assign a new entry, overwriting any existing entry.
Definition dictionary.C:765
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Abstract base class for optimisation methods.
virtual void clearSensitivities()
Clear all adjoint sensitivities.
virtual void initialize()
Initialization. Construct primal and adjoint solvers.
virtual void resetTime(const dimensionedScalar startTime, const label startTimeIndex, const scalar endTime)
Reset time.
const word managerType_
Type of the optimisation manager (singleRun, (un)steadyOptimisation).
autoPtr< designVariablesUpdate > dvUpdate_
Helper class managing parts of the optimisation.
virtual void moveDesignVariables()
Update design variables.
Switch shouldUpdateDesignVariables_
Switch defining if the design variables should be updated or not.
virtual void solvePrimalEquations()
Solve all primal equations.
virtual void updateDesignVariables()
Update design variables.
virtual void updatePrimalBasedQuantities()
Solve all primal equations.
autoPtr< designVariables > designVars_
Design variables of the optimisation problem.
Time & time_
Reference to the time.
PtrList< adjointSolverManager > adjointSolverManagers_
List of adjoint solver managers to be included in the optimisation.
virtual void solveAdjointEquations()
Solve all adjoint equations.
void lineSearchUpdate()
Update design variables using a line-search.
PtrList< primalSolver > primalSolvers_
List of primal solvers to be included in the optimisation.
static autoPtr< optimisationManager > New(fvMesh &mesh)
Return a reference to the selected turbulence model.
fvMesh & mesh_
Reference to the mesh.
virtual void computeSensitivities()
Compute all adjoint sensitivities.
virtual bool read()
Changes in case of run-time update of optimisationDict.
void fixedStepUpdate()
Update design variables using a fixed step.
Base class for primal solvers.
static autoPtr< primalSolver > New(fvMesh &mesh, const word &managerType, const dictionary &dict, const word &solverName)
Return a reference to the selected primal solver.
virtual bool read()
Read object.
A class for managing temporary objects.
Definition tmp.H:75
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition tmpI.H:235
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 FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition error.H:637
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
Namespace for OpenFOAM.
List< word > wordList
List of word.
Definition fileName.H:60
int system(const std::string &command, const bool bg=false)
Execute the specified command via the shell.
Definition POSIX.C:1704
messageStream Info
Information stream (stdout output on master, null elsewhere).
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
uint8_t direction
Definition direction.H:49
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299