Loading...
Searching...
No Matches
designVariablesUpdate.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
33#include "adjointNull.H"
34#include "IOmanip.H"
36
37// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38
39namespace Foam
41 defineTypeNameAndDebug(designVariablesUpdate, 0);
42}
43
44
45// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
46
48(
50) const
51{
52 // Figure out number of adjoint solvers corresponding to constraints.
53 // Looks in all operating points
54 label nConstraints(0);
55 for (const adjointSolverManager& adjManagerI : adjointSolvManagers_)
56 {
57 nConstraints += adjManagerI.nConstraints();
58 }
59 // Add constraints that might emerge from the design variables
60 tmp<scalarField> designVarsConstraints = designVars_().constraintValues();
61 if (designVarsConstraints)
62 {
63 nConstraints += designVarsConstraints().size();
64 }
65 return nConstraints;
66}
67
68
70{
71 label n(0);
72 for (adjointSolverManager& adjSolvManager : adjointSolvManagers_)
73 {
74 for (adjointSolver& solver : adjSolvManager.adjointSolvers())
75 {
76 if (!isA<adjointNull>(solver))
77 {
78 ++n;
79 }
80 }
81 }
82 return n;
83}
84
85
87{
88 unsigned int width(IOstream::defaultPrecision() + 5);
89 CPUcostFile_
90 << setw(width) << "#Cycle" << " "
91 << setw(width) << "LineSearchIters" << " "
92 << setw(width) << "CycleCPUcost" << " "
93 << setw(width) << "CyclePrimalSolutions" << " "
94 << setw(width) << "CycleAdjointSolutions" << " "
95 << setw(width) << "TotalCPUcost" << " "
96 << setw(width) << "TotalPrimalSolutions" << " "
97 << setw(width) << "TotalAdjointSolutions" << endl;
98}
99
100
101void Foam::designVariablesUpdate::writeToCostFile(bool zeroAdjointSolns)
102{
103 unsigned int width(IOstream::defaultPrecision() + 5);
104 label cyclePrimalSolutions(nPrimalsPerCycle_);
105 label cycleAdjointSolutions(nAdjointsPerCycle_);
106 label lineSearchIters(1);
107 if (lineSearch_)
108 {
109 lineSearchIters = lineSearch_().innerIter();
110 cyclePrimalSolutions *= lineSearchIters;
111 if (lineSearch_().computeGradient())
112 {
113 cycleAdjointSolutions *= lineSearchIters;
114 }
115 }
116 if (zeroAdjointSolns)
117 {
118 cycleAdjointSolutions = 0;
119 }
120 nPrimalSolutions_ += cyclePrimalSolutions;
121 nAdjointSolutions_ += cycleAdjointSolutions;
122 const scalar elapsedCpuTime = mesh_.time().elapsedCpuTime();
123 const scalar cycleCost = elapsedCpuTime - CPUcost_;
124 CPUcost_ = elapsedCpuTime;
125
126 CPUcostFile_
127 << setw(width) << mesh_.time().timeName() << " "
128 << setw(width) << lineSearchIters << " "
129 << setw(width) << cycleCost << " "
130 << setw(width) << cyclePrimalSolutions << " "
131 << setw(width) << cycleAdjointSolutions << " "
132 << setw(width) << CPUcost_ << " "
133 << setw(width) << nPrimalSolutions_ << " "
134 << setw(width) << nAdjointSolutions_ << endl;
135
136}
137
138
139// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
140
141Foam::designVariablesUpdate::designVariablesUpdate
142(
143 fvMesh& mesh,
144 const dictionary& dict,
146 autoPtr<designVariables>& designVars
147)
148:
149 mesh_(mesh),
150 dict_(dict),
151 adjointSolvManagers_(adjointSolverManagers),
152 designVars_(designVars),
153 updateMethod_
154 (
156 (
157 mesh_,
158 dict_.subDict("updateMethod"),
159 designVars_,
160 nConstraints(adjointSolverManagers)
161 )
162 ),
163 lineSearch_
164 (
166 (
167 dict_.subDict("updateMethod").subOrEmptyDict("lineSearch"),
168 mesh.time(),
169 updateMethod_.ref()
170 )
171 ),
172 CPUcostFile_(mesh_.time().globalPath()/"optimisation"/"CPUcost"),
173 nPrimalsPerCycle_(adjointSolverManagers.size()),
174 nAdjointsPerCycle_(nAdjointSolvers()),
175 nPrimalSolutions_(nPrimalsPerCycle_),
176 nAdjointSolutions_(nAdjointsPerCycle_),
177 CPUcost_(0),
178 designVarsThreshold_(nullptr),
179 objectiveThreshold_(nullptr),
180 convergenceCriteriaDefined_(false),
181 feasibilityThreshold_
182 (
183 dict.subOrEmptyDict("convergence").
184 getOrDefault<scalar>("feasibilityThreshold", 1.e-06)
185 )
186{
187 dictionary convergenceDict = dict.subOrEmptyDict("convergence");
188 if (convergenceDict.found("designVariables"))
189 {
190 designVarsThreshold_.reset
191 (new scalar(convergenceDict.get<scalar>("designVariables")));
192 }
193 if (convergenceDict.found("objective"))
194 {
195 objectiveThreshold_.reset
196 (new scalar(convergenceDict.get<scalar>("objective")));
197 }
199 // Check whether eta of maxInitChange are set
200 if (!designVars_().isMaxInitChangeSet() && !updateMethod_().initialEtaSet())
201 {
203 << "Neither eta (updateMethod) or maxInitChange (designVariables) "
204 << "is set."
205 << exit(FatalError);
206 }
207
208 label nConstr(nConstraints(adjointSolvManagers_));
209 // Sanity checks for combinations of number of constraints and
210 // optimisation methods
211 if
212 (
213 nConstr
215 )
216 {
217 const auto& cnstrTable =
218 *(constrainedOptimisationMethod::dictionaryConstructorTablePtr_);
219
220 // Has constraints but is not a constraint optimisation method
222 << "Found " << nConstr << " adjoint solvers corresponding to "
223 << "constraints but the optimisation method ("
224 << updateMethod_().type()
225 << ") is not a constrainedOptimisationMethod." << nl
226 << "Available constrainedOptimisationMethods:" << nl
227 << cnstrTable.sortedToc()
228 << exit(FatalError);
229 }
230 else if
231 (
232 !nConstr
234 )
235 {
236 // Does not have constraints but is a constrained optimisation method
238 << "Did not find any adjoint solvers corresponding to "
239 << "constraints but the optimisation method ("
240 << updateMethod_().type()
241 << ") is a constrainedOptimisationMethod." << nl << nl
242 << "This can cause some constraintOptimisationMethods to misbehave."
243 << nl << nl
244 << "Either the isConstraint bool is not set in one of the adjoint "
245 << "solvers or you should consider using an updateMethod "
246 << "that is not a constrainedOptimisationMethod"
247 << nl << endl;
248 }
249
251 {
252 Info<< "Optimisation will run until the max. scaled correction "
253 << "of the design variables is < " << designVarsThreshold_()
254 << endl;
255 }
257 {
258 Info<< "Optimisation will run until the relative update of the "
259 << "objective function is < " << objectiveThreshold_()
260 << endl;
261 }
263 {
264 Info<< "No convergence criterion defined for optimsation" << nl
265 << "It will run for " << mesh_.time().endTime().value()
266 << " optimisation cycles " << nl << endl;
267 }
268 Info<< "Feasibility threshold is " << feasibilityThreshold_ << endl;
269
270 // Write header of the cost file
272}
273
274
275// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * //
276
278{
279 // Compute update of the design variables
280 tmp<scalarField> tcorrection(computeDirection());
281 scalarField& correction = tcorrection.ref();
282
283 // Set the old value of the objective function
284 setOldObjectiveValue();
285
286 // Update the design variables
287 designVars_->update(correction);
289 // If direction has been scaled (say by setting the initial eta), the
290 // old correction has to be updated
292}
293
294
296{
297 // Multiply with line search step, if necessary
299 if (lineSearch_)
300 {
301 lineSearch_->updateCorrection(correction);
303
304 // Update the design variables
305 designVars_->update(correction);
306}
307
308
310{
311 updateGradientsAndValues();
312 updateMethod_->computeCorrection();
313 scalarField& correction = updateMethod_->returnCorrection();
314
315 // Compute eta if needed
316 if (!updateMethod_->initialEtaSet() || designVars_->resetEta())
317 {
318 const scalar eta(designVars_->computeEta(correction));
319 updateMethod_->modifyStep(eta);
320 updateMethod_->initialEtaSet() = true;
322
323 // Intentionally copies result to new field
325}
326
327
329{
330 scalarField objectiveSens;
331 PtrList<scalarField> constraintSens;
332 scalar objectiveValue(Zero);
333 DynamicList<scalar> constraintValues;
334
335 for (adjointSolverManager& adjSolvManager : adjointSolvManagers_)
336 {
337 const scalar opWeight = adjSolvManager.operatingPointWeight();
338
339 // Aggregate sensitivities of solvers corresponding to objectives
340 // (i.e. not constraints)
341 tmp<scalarField> tadjointSolverManagerSens =
342 adjSolvManager.aggregateSensitivities();
343
344 // Aggregate objective values of solvers corresponding to objectives
345 // (i.e. not constraints)
346 objectiveValue += opWeight*adjSolvManager.objectiveValue();
347
348 // Get constraint sensitivities
349 PtrList<scalarField> adjointSolverManagerConstSens =
350 adjSolvManager.constraintSensitivities();
351
352 // Post process sensitivities if needed.
353 // Done here since each adjointSolverManager might post-process
354 // its sensitivities in a different way
355 designVars_->postProcessSens
356 (
357 tadjointSolverManagerSens.ref(),
358 adjointSolverManagerConstSens,
359 adjSolvManager.adjointSolversNames(),
360 adjSolvManager.isMaster()
361 );
362
363 if (objectiveSens.empty())
364 {
365 objectiveSens.setSize(tadjointSolverManagerSens().size(), Zero);
366 }
367
368 // Accumulate sensitivities
369 objectiveSens += opWeight*tadjointSolverManagerSens();
370 forAll(adjointSolverManagerConstSens, sI)
371 {
372 constraintSens.
373 push_back(adjointSolverManagerConstSens.set(sI, nullptr));
374 }
375 constraintValues.push_back(adjSolvManager.constraintValues());
376 }
377 // Add contraint values and gradients from the design variables
378 tmp<scalarField> designVarsConstValues = designVars_->constraintValues();
379 PtrList<scalarField> designVarsConstDerivs =
380 designVars_->constraintDerivatives();
381 if (designVarsConstValues && designVarsConstDerivs.size())
382 {
383 if (designVarsConstValues().size() != designVarsConstDerivs.size())
384 {
386 << "Size of design variables constraints and derivatives differ"
387 << endl
388 << exit(FatalError);
389 }
390 constraintValues.push_back(designVarsConstValues());
391 constraintSens.push_back(std::move(designVarsConstDerivs));
392 }
393
394 // Update objective/constraint values/gradients, known by the update method
395 updateMethod_->setObjectiveDeriv(objectiveSens);
396 updateMethod_->setConstraintDeriv(constraintSens);
397 updateMethod_->setObjectiveValue(objectiveValue);
398 updateMethod_->setConstraintValues
399 (scalarField(std::move(constraintValues)));
400}
401
402
404{
405 scalar objectiveValue(Zero);
406 for (adjointSolverManager& adjSolvManager : adjointSolvManagers_)
407 {
408 const scalar opWeight = adjSolvManager.operatingPointWeight();
409 objectiveValue += opWeight*adjSolvManager.objectiveValue();
410 }
411 return objectiveValue;
412}
413
416{
417 updateMethod_->setObjectiveValueOld(computeObjectiveValue());
418}
419
420
422{
423 // Compute new objective and constraint values and update the ones
424 // in updateMethod
425 scalar objectiveValue(Zero);
426 DynamicList<scalar> constraintValues;
427
428 for (adjointSolverManager& adjSolvManager : adjointSolvManagers_)
429 {
430 const scalar opWeight = adjSolvManager.operatingPointWeight();
431
432 objectiveValue += opWeight*adjSolvManager.objectiveValue();
433 constraintValues.push_back(adjSolvManager.constraintValues());
434 }
435
436 // Add constraints directly imposed to the design variables
437 tmp<scalarField> designVarsConstValues = designVars_->constraintValues();
438 if (designVarsConstValues)
439 {
440 constraintValues.push_back(designVarsConstValues());
441 }
442 updateMethod_->setObjectiveValue(objectiveValue);
443 updateMethod_->setConstraintValues
444 (scalarField(std::move(constraintValues)));
445
446 return updateMethod_->computeMeritFunction();
447}
448
451{
452 return updateMethod_->meritFunctionDirectionalDerivative();
453}
454
455
457(
458 const scalarField& oldCorrection
459)
460{
461 updateMethod_->updateOldCorrection(oldCorrection);
462}
463
464
467 updateMethod_->writeAuxiliaryData();
468 designVars_->writeDesignVars();
470}
471
472
474{
475 updateOldCorrection(oldCorrection);
476 write();
477 designVars_->evolveNumber();
478 if (lineSearch_)
479 {
480 lineSearch_->postUpdate();
481
482 // Append additional empty line at the end of the instantaneous
483 // objective file to indicate the end of the block corresponding to
484 // this optimisation cycle
485 for (adjointSolverManager& am : adjointSolvManagers_)
486 {
487 for (adjointSolver& as : am.adjointSolvers())
488 {
490 as.getObjectiveManager().getObjectiveFunctions();
491 for (objective& obj : objectives)
492 {
493 obj.writeInstantaneousSeparator();
494 }
496 }
497 checkConvergence();
498 }
499}
500
501
503{
504 if (!convergenceCriteriaDefined_)
505 {
506 return;
507 }
508
509 bool converged(false);
510 const scalarField& oldCorrection = updateMethod_->returnCorrection();
511 // Design variables convergence check
512 if (designVarsThreshold_)
513 {
514 const labelList& activeVarIDs =
515 designVars_->activeDesignVariables();
516 const scalarField correction(oldCorrection, activeVarIDs);
517 const scalarField activeVars(designVars_->getVars(), activeVarIDs);
518 const scalar scaledCorrection =
519 gMax(mag(correction)/(mag(activeVars) + SMALL));
521 << "Current correction " << correction << nl
522 << "Active vars " << activeVars << endl;
523 Info<< "Max. scaled correction of the design variables = "
524 << scaledCorrection
525 << endl;
526 if (scaledCorrection < designVarsThreshold_())
527 {
528 Info<< tab << "Design variables have converged " << endl;
529 converged = true;
530 }
531 }
532 // Objective convergence check
533 if (objectiveThreshold_ && updateMethod_->getObjectiveValueOld())
534 {
535 const scalar newObjective = computeObjectiveValue();
536 const scalar oldObjective = updateMethod_->getObjectiveValueOld()();
537 const scalar relativeUpdate =
538 mag(newObjective - oldObjective)/(mag(oldObjective) + SMALL);
539 Info<< "Relative change of the objective value = "
540 << relativeUpdate
541 << endl;
542 if (relativeUpdate < objectiveThreshold_())
543 {
544 Info<< tab << "Objective function has converged " << endl;
545 converged = true;
546 }
547 }
548 // Feasibility check
549 const scalarField& constraints = updateMethod_->getConstraintValues();
550 const scalar feasibility = sum(pos(constraints)*constraints);
551 Info<< "Feasibility = " << feasibility << endl;
552 if (converged && feasibility < feasibilityThreshold_)
553 {
554 Info<< "Stopping criteria met and all constraints satisfied." << nl
555 << "Optimisation has converged, stopping ..." << nl << nl
556 << "End" << nl
557 << endl;
558 // Force writing of all objective and constraint functions, to get
559 // the converged results to files
560 for (adjointSolverManager& am : adjointSolvManagers_)
561 {
562 for (adjointSolver& as : am.adjointSolvers())
563 {
564 // Use dummy weighted objective
565 as.getObjectiveManager().writeObjectives();
566 }
567 }
568 writeToCostFile(true);
569 if (UPstream::parRun())
570 {
572 }
573 else
574 {
575 std::exit(0);
576 }
577 }
578}
579
580
581// ************************************************************************* //
Istream and Ostream manipulators taking arguments.
label n
if(patchID !=-1)
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
void push_back(const T &val)
Copy append an element to the end of this list.
static unsigned int defaultPrecision() noexcept
Return the default precision.
Definition IOstream.H:437
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
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie,...
Definition PtrList.H:171
void push_back(T *ptr)
Append an element to the end of the list.
Definition PtrListI.H:131
bool empty() const noexcept
True if List is empty (ie, size() is zero).
Definition UList.H:701
static bool parRun(const bool on) noexcept
Set as parallel run on/off.
Definition UPstream.H:1669
static void exit(int errNo=1)
Shutdown (finalize) MPI as required and exit program with errNo.
Definition UPstream.C:61
label size() const noexcept
The number of entries in the list.
Definition UPtrListI.H:106
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
A class encapsulating functionality neccessary to perform an optimisation loop, such as updating the ...
void checkConvergence()
Check if the optimisation loop has converged based on the provided criteria.
scalar meritFunctionDirectionalDerivative()
Derivative of the merit function.
PtrList< adjointSolverManager > & adjointSolvManagers_
bool convergenceCriteriaDefined_
Is at least a single convergence criterion defined.
void writeCPUcostHeader()
Write CPU cost header.
scalar computeObjectiveValue()
Sum objective values from all adjointSolverManagers.
autoPtr< updateMethod > updateMethod_
Method to update the design variables based on the provided sensitivity derivatives.
void updateOldCorrection(const scalarField &)
Update old correction. Needed for quasi-Newton Methods.
label nAdjointSolvers() const
Get total number of adjoint solvers.
void updateGradientsAndValues()
Compute cumulative objective and constraint gradients.
label nConstraints(PtrList< adjointSolverManager > &adjointSolverManagers) const
Get the number of adjoint solvers that correspond to constraints.
tmp< scalarField > computeDirection()
Compute update direction.
label nAdjointSolutions_
Adjoint evaluations performed so far.
autoPtr< scalar > objectiveThreshold_
The relative update of the objective value w.r.t. to its last value should be smaller than this value...
autoPtr< lineSearch > lineSearch_
Line search mechanism to approximate the update step length.
scalar feasibilityThreshold_
In case of a constrained optimisation, the sum of positive constraints should be lower than this valu...
scalar computeMeritFunction()
Compute the merit function of the optimisation problem.
label nPrimalsPerCycle_
Primal solutions per optimisation cycle.
void write()
Write useful quantities to files.
label nAdjointsPerCycle_
Adjoint solutions per optimisation cycle.
void postUpdate(const scalarField &oldCorrection)
Steps to be executed after the susccessfull update of the design varibles, i.e. the last step of line...
void update()
Update design variables.
autoPtr< designVariables > & designVars_
autoPtr< scalar > designVarsThreshold_
The maximum of the correction/designVariables values must be lower that this threshold to consider th...
label nPrimalSolutions_
Primal evaluations performed so far.
void writeToCostFile(bool zeroAdjointSolns=false)
Write to cost file.
scalar CPUcost_
CPU cost (in seconds).
void setOldObjectiveValue()
Set the old objective value known by the updateMethod.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry (const access) with the given keyword.
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Abstract base class for line search methods.
Definition lineSearch.H:54
Abstract base class for objective functions. No point in making this runTime selectable since its chi...
Definition objective.H:58
Base solver class.
Definition solver.H:48
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
Abstract base class for optimisation methods.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
rDeltaT ref()
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
#define DebugInfo
Report an information message using Foam::Info.
#define WarningInFunction
Report a warning using Foam::Warning.
Namespace for OpenFOAM.
dimensionedScalar pos(const dimensionedScalar &ds)
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.
List< label > labelList
A List of labels.
Definition List.H:62
messageStream Info
Information stream (stdout output on master, null elsewhere).
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
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
uint8_t direction
Definition direction.H:49
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1, const label comm)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix by subtracting the matrix multiplied by the current fi...
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
Type gMax(const FieldField< Field, Type > &f)
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
constexpr char tab
The tab '\t' character(0x09).
Definition Ostream.H:49
runTime write()
Macros to ease declaration of run-time selection tables.
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
PtrList< adjointSolverManager > & adjointSolverManagers
Definition createFields.H:8