48 const label startTimeIndex,
66 dvUpdate_->getDesignVariables()->storeDesignVariables();
69 scalar meritFunction = dvUpdate_->computeMeritFunction();
70 lineSrch->setOldMeritValue(meritFunction);
71 dvUpdate_->setOldObjectiveValue();
74 scalar dirDerivative =
75 dvUpdate_->meritFunctionDirectionalDerivative();
76 lineSrch->setDeriv(dirDerivative);
84 while (lineSrch->loop())
86 Info<<
"\n- - - - - - - - - - - - - - -" <<
endl;
87 Info<<
"Line search iteration " << lineSrch->innerIter() <<
endl;
88 Info<<
"- - - - - - - - - - - - - - -\n" <<
endl;
95 const label startTimeIndex = mesh_.time().timeIndex();
96 const scalar primalEndTime = mesh_.time().endTime().value();
98 solvePrimalEquations();
101 meritFunction = dvUpdate_->computeMeritFunction();
102 lineSrch->setNewMeritValue(meritFunction);
104 if (lineSrch->computeGradient())
107 clearSensitivities();
110 solveAdjointEquations();
113 dvUpdate_->updateGradientsAndValues();
117 dvUpdate_->meritFunctionDirectionalDerivative();
118 lineSrch->setNewDeriv(dirDerivative);
121 if (lineSrch->converged())
124 Info<<
"Line search converged in " << lineSrch->innerIter()
125 <<
" iterations." <<
endl;
127 dvUpdate_->postUpdate(scaledCorrection);
133 if (lineSrch->innerIter() == lineSrch->maxIters())
135 Info<<
"Line search reached max. number of iterations.\n"
136 <<
"Proceeding to the next optimisation cycle" <<
endl;
138 dvUpdate_->postUpdate(scaledCorrection);
144 this->resetTime(
startTime, startTimeIndex, primalEndTime);
145 dvUpdate_->getDesignVariables()->resetDesignVariables();
146 lineSrch->updateStep();
153 if (!lineSrch->computeGradient())
156 clearSensitivities();
167 if (shouldUpdateDesignVariables_)
169 moveDesignVariables();
173 solvePrimalEquations();
176 dvUpdate_->checkConvergence();
179 clearSensitivities();
205 dictionary& primalSolversDict = subDict(
"primalSolvers");
206 const wordList& primalSolverNames = primalSolversDict.
toc();
210 forAll(primalSolvers_, solveri)
213 primalSolversDict.
subDict(primalSolverNames[solveri]);
214 if (primalSolvers_.size() > 1)
216 solverDict.
add<
bool>(
"useSolverNameForFields",
true);
226 primalSolverNames[solveri]
232 const dictionary& adjointManagersDict = subDict(
"adjointManagers");
233 const wordList adjointManagerNames = adjointManagersDict.toc();
234 adjointSolverManagers_.
setSize(adjointManagerNames.size());
239 label nNotNullAdjointSolvers(0);
240 for (
const word& adjManager : adjointManagerNames)
243 adjointManagersDict.subDict(adjManager).subDict(
"adjointSolvers");
244 const wordList adjointSolverNames = adjSolversDict.toc();
245 for (
const word& adjSolver : adjointSolverNames)
247 if (adjSolversDict.subDict(adjSolver).get<word>(
"type") !=
"null")
249 ++nNotNullAdjointSolvers;
254 << nNotNullAdjointSolvers
255 <<
" adjoint solvers that allocate fields"
257 bool overrideUseSolverName(nNotNullAdjointSolvers > 1);
259 forAll(adjointSolverManagers_, manageri)
261 adjointSolverManagers_.set
269 adjointManagersDict.subDict(adjointManagerNames[manageri]),
270 overrideUseSolverName
276 if (primalSolvers_.size() > 1)
280 if (!solveri.useSolverNameForFields())
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."
292 if (nNotNullAdjointSolvers > 1)
299 if (!asI.useSolverNameForFields())
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."
313 designVars_().addFvOptions(primalSolvers_, adjointSolverManagers_);
320Foam::optimisationManager::optimisationManager(
fvMesh&
mesh)
335 time_(const_cast<
Time&>(
mesh.time())),
336 designVars_(nullptr),
338 adjointSolverManagers_(),
339 managerType_(get<
word>(
"optimisationManager")),
341 shouldUpdateDesignVariables_(true)
347 if (designVarsDictPtr)
366 mesh.time().system(),
374 const word modelType(
dict.get<
word>(
"optimisationManager"));
376 Info<<
"optimisationManager type : " << modelType <<
endl;
378 auto* ctorPtr = dictionaryConstructorTable(modelType);
385 "optimisationManager",
387 *dictionaryConstructorTablePtr_
391 return autoPtr<optimisationManager>(ctorPtr(
mesh));
402 const dictionary& primalSolversDict = subDict(
"primalSolvers");
405 sol.readDict(primalSolversDict.
subDict(sol.solverName()));
408 const dictionary& adjointManagersDict = subDict(
"adjointManagers");
411 man.readDict(adjointManagersDict.
subDict(man.managerName()));
416 designVars_->readDict
417 (subDict(
"optimisation").subDict(
"designVariables"));
431 if (dvUpdate_->getLineSearch())
445 forAll(primalSolvers_, psI)
455 forAll(adjointSolverManagers_, amI)
465 forAll(adjointSolverManagers_, amI)
474 for (adjointSolverManager& adjSolvManager : adjointSolverManagers_)
476 adjSolvManager.clearSensitivities();
483 forAll(adjointSolverManagers_, amI)
485 PtrList<adjointSolver>& adjointSolvers =
486 adjointSolverManagers_[amI].adjointSolvers();
488 forAll(adjointSolvers, asI)
490 adjointSolvers[asI].updatePrimalBasedQuantities();
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,...
const Time & time() const noexcept
Return Time associated with the objectRegistry.
void setSize(label n)
Alias for resize().
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
void size(const label n)
Older name for setAddressableSize.
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.
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,...
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.
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.
entry * add(entry *entryPtr, bool mergeEntry=false)
Add a new entry.
wordList toc() const
Return the table of contents.
entry * set(entry *entryPtr)
Assign a new entry, overwriting any existing entry.
Mesh data needed to do the Finite Volume discretisation.
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.
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
A class for handling words, derived from Foam::string.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
List< word > wordList
List of word.
int system(const std::string &command, const bool bg=false)
Execute the specified command via the shell.
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.
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)
constexpr char nl
The newline '\n' character (0x0a).
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.
#define forAll(list, i)
Loop across all elements in list.