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 -------------------------------------------------------------------------------
12 License
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 
30 #include "adjointSolverManager.H"
31 #include "designVariablesUpdate.H"
33 #include "adjointNull.H"
34 #include "IOmanip.H"
35 #include "runTimeSelectionTables.H"
36 
37 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
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 
101 void 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 
140 (
141  const scalarField& oldCorrection
142 )
143 {
144  bool converged(false);
145  // Design variables convergence check
146  if (designVarsThreshold_)
147  {
148  const labelList& activeVarIDs =
149  designVars_->activeDesignVariables();
150  const scalarField correction(oldCorrection, activeVarIDs);
151  const scalarField activeVars(designVars_->getVars(), activeVarIDs);
152  const scalar scaledCorrection =
153  gMax(mag(correction)/(mag(activeVars) + SMALL));
154  DebugInfo
155  << "Current correction " << correction << nl
156  << "Active vars " << activeVars << endl;
157  Info<< "Max. scaled correction of the design variables = "
158  << scaledCorrection
159  << endl;
160  if (scaledCorrection < designVarsThreshold_())
161  {
162  Info<< tab << "Design variables have converged " << endl;
163  converged = true;
164  }
165  }
166  // Objective convergence check
167  if (objectiveThreshold_)
168  {
169  const scalar newObjective = computeObjectiveValue();
170  const scalar oldObjective = updateMethod_->getObjectiveValueOld();
171  const scalar relativeUpdate =
172  mag(newObjective - oldObjective)/(mag(oldObjective) + SMALL);
173  Info<< "Relative change of the objective value = "
174  << relativeUpdate
175  << endl;
176  if (relativeUpdate < objectiveThreshold_())
177  {
178  Info<< tab << "Objective function has converged " << endl;
179  converged = true;
180  }
181  }
182  // Feasibility check
183  const scalarField& constraints = updateMethod_->getConstraintValues();
184  const scalar feasibility = sum(pos(constraints)*constraints);
185  Info<< "Feasibility = " << feasibility << endl;
186  if (converged && feasibility < feasibilityThreshold_)
187  {
188  Info<< "Stopping criteria met and all constraints satisfied." << nl
189  << "Optimisation has converged, stopping ..." << nl << nl
190  << "End" << nl
191  << endl;
192  // Force writing of all objective and constraint functions, to get
193  // the converged results to files
194  for (adjointSolverManager& am : adjointSolvManagers_)
195  {
196  for (adjointSolver& as : am.adjointSolvers())
197  {
198  // Use dummy weighted objective
199  as.getObjectiveManager().writeObjectives();
200  }
201  }
202  writeToCostFile(true);
203  if (UPstream::parRun())
204  {
205  UPstream::exit(0);
206  }
207  else
208  {
209  std::exit(0);
210  }
211  }
212 }
213 
214 
215 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
216 
217 Foam::designVariablesUpdate::designVariablesUpdate
218 (
219  fvMesh& mesh,
220  const dictionary& dict,
221  PtrList<adjointSolverManager>& adjointSolverManagers,
222  autoPtr<designVariables>& designVars
223 )
224 :
225  mesh_(mesh),
226  dict_(dict),
227  adjointSolvManagers_(adjointSolverManagers),
228  designVars_(designVars),
229  updateMethod_
230  (
231  updateMethod::New
232  (
233  mesh_,
234  dict_.subDict("updateMethod"),
235  designVars_,
236  nConstraints(adjointSolverManagers)
237  )
238  ),
239  lineSearch_
240  (
241  lineSearch::New
242  (
243  dict_.subDict("updateMethod").subOrEmptyDict("lineSearch"),
244  mesh.time(),
245  updateMethod_.ref()
246  )
247  ),
248  CPUcostFile_(mesh_.time().globalPath()/"optimisation"/"CPUcost"),
249  nPrimalsPerCycle_(adjointSolverManagers.size()),
250  nAdjointsPerCycle_(nAdjointSolvers()),
251  nPrimalSolutions_(nPrimalsPerCycle_),
252  nAdjointSolutions_(nAdjointsPerCycle_),
253  CPUcost_(0),
254  designVarsThreshold_(nullptr),
255  objectiveThreshold_(nullptr),
256  convergenceCriteriaDefined_(false),
257  feasibilityThreshold_
258  (
259  dict.subOrEmptyDict("convergence").
260  getOrDefault<scalar>("feasibilityThreshold", 1.e-06)
261  )
262 {
263  dictionary convergenceDict = dict.subOrEmptyDict("convergence");
264  if (convergenceDict.found("designVariables"))
265  {
267  (new scalar(convergenceDict.get<scalar>("designVariables")));
268  }
269  if (convergenceDict.found("objective"))
270  {
272  (new scalar(convergenceDict.get<scalar>("objective")));
273  }
275  // Check whether eta of maxInitChange are set
276  if (!designVars_().isMaxInitChangeSet() && !updateMethod_().initialEtaSet())
277  {
279  << "Neither eta (updateMethod) or maxInitChange (designVariables) "
280  << "is set."
281  << exit(FatalError);
282  }
283 
284  label nConstr(nConstraints(adjointSolvManagers_));
285  // Sanity checks for combinations of number of constraints and
286  // optimisation methods
287  if
288  (
289  nConstr
290  && !isA<constrainedOptimisationMethod>(updateMethod_())
291  )
292  {
293  const auto& cnstrTable =
294  *(constrainedOptimisationMethod::dictionaryConstructorTablePtr_);
295 
296  // Has constraints but is not a constraint optimisation method
298  << "Found " << nConstr << " adjoint solvers corresponding to "
299  << "constraints but the optimisation method ("
300  << updateMethod_().type()
301  << ") is not a constrainedOptimisationMethod." << nl
302  << "Available constrainedOptimisationMethods:" << nl
303  << cnstrTable.sortedToc()
304  << exit(FatalError);
305  }
306  else if
307  (
308  !nConstr
309  && isA<constrainedOptimisationMethod>(updateMethod_())
310  )
311  {
312  // Does not have constraints but is a constrained optimisation method
314  << "Did not find any adjoint solvers corresponding to "
315  << "constraints but the optimisation method ("
316  << updateMethod_().type()
317  << ") is a constrainedOptimisationMethod." << nl << nl
318  << "This can cause some constraintOptimisationMethods to misbehave."
319  << nl << nl
320  << "Either the isConstraint bool is not set in one of the adjoint "
321  << "solvers or you should consider using an updateMethod "
322  << "that is not a constrainedOptimisationMethod"
323  << nl << endl;
324  }
325 
327  {
328  Info<< "Optimisation will run until the max. scaled correction "
329  << "of the design variables is < " << designVarsThreshold_()
330  << endl;
331  }
333  {
334  Info<< "Optimisation will run until the relative update of the "
335  << "objective function is < " << objectiveThreshold_()
336  << endl;
337  }
339  {
340  Info<< "No convergence criterion defined for optimsation" << nl
341  << "It will run for " << mesh_.time().endTime().value()
342  << " optimisation cycles " << nl << endl;
343  }
344  Info<< "Feasibility threshold is " << feasibilityThreshold_ << endl;
345 
346  // Write header of the cost file
348 }
349 
350 
351 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * //
352 
354 {
355  // Compute update of the design variables
356  tmp<scalarField> tcorrection(computeDirection());
357  scalarField& correction = tcorrection.ref();
358 
359  // Set the old value of the objective function
360  setOldObjectiveValue();
361 
362  // Update the design variables
363  designVars_->update(correction);
365  // If direction has been scaled (say by setting the initial eta), the
366  // old correction has to be updated
367  postUpdate(correction);
368 }
369 
370 
372 {
373  // Multiply with line search step, if necessary
375  if (lineSearch_.valid())
376  {
377  lineSearch_->updateCorrection(correction);
378  }
379 
380  // Update the design variables
381  designVars_->update(correction);
382 }
383 
384 
386 {
387  updateGradientsAndValues();
388  updateMethod_->computeCorrection();
389  scalarField& correction = updateMethod_->returnCorrection();
390 
391  // Compute eta if needed
392  if (!updateMethod_->initialEtaSet() || designVars_->resetEta())
393  {
394  const scalar eta(designVars_->computeEta(correction));
395  updateMethod_->modifyStep(eta);
396  updateMethod_->initialEtaSet() = true;
397  }
398 
399  // Intentionally copies result to new field
401 }
402 
403 
405 {
406  scalarField objectiveSens;
407  PtrList<scalarField> constraintSens;
408  scalar objectiveValue(Zero);
409  DynamicList<scalar> constraintValues;
410 
411  for (adjointSolverManager& adjSolvManager : adjointSolvManagers_)
412  {
413  const scalar opWeight = adjSolvManager.operatingPointWeight();
414 
415  // Aggregate sensitivities of solvers corresponding to objectives
416  // (i.e. not constraints)
417  tmp<scalarField> tadjointSolverManagerSens =
418  adjSolvManager.aggregateSensitivities();
419 
420  // Aggregate objective values of solvers corresponding to objectives
421  // (i.e. not constraints)
422  objectiveValue += opWeight*adjSolvManager.objectiveValue();
423 
424  // Get constraint sensitivities
425  PtrList<scalarField> adjointSolverManagerConstSens =
426  adjSolvManager.constraintSensitivities();
427 
428  // Post process sensitivities if needed.
429  // Done here since each adjointSolverManager might post-process
430  // its sensitivities in a different way
431  designVars_->postProcessSens
432  (
433  tadjointSolverManagerSens.ref(),
434  adjointSolverManagerConstSens,
435  adjSolvManager.adjointSolversNames(),
436  adjSolvManager.isMaster()
437  );
438 
439  if (objectiveSens.empty())
440  {
441  objectiveSens.setSize(tadjointSolverManagerSens().size(), Zero);
442  }
443 
444  // Accumulate sensitivities
445  objectiveSens += opWeight*tadjointSolverManagerSens();
446  forAll(adjointSolverManagerConstSens, sI)
447  {
448  constraintSens.
449  push_back(adjointSolverManagerConstSens.set(sI, nullptr));
450  }
451  constraintValues.push_back(adjSolvManager.constraintValues());
452  }
453  // Add contraint values and gradients from the design variables
454  tmp<scalarField> designVarsConstValues = designVars_->constraintValues();
455  PtrList<scalarField> designVarsConstDerivs =
456  designVars_->constraintDerivatives();
457  if (designVarsConstValues && designVarsConstDerivs.size())
458  {
459  if (designVarsConstValues().size() != designVarsConstDerivs.size())
460  {
462  << "Size of design variables constraints and derivatives differ"
463  << endl
464  << exit(FatalError);
465  }
466  constraintValues.push_back(designVarsConstValues());
467  constraintSens.push_back(std::move(designVarsConstDerivs));
468  }
469 
470  // Update objective/constraint values/gradients, known by the update method
471  updateMethod_->setObjectiveDeriv(objectiveSens);
472  updateMethod_->setConstraintDeriv(constraintSens);
473  updateMethod_->setObjectiveValue(objectiveValue);
474  updateMethod_->setConstraintValues
475  (scalarField(std::move(constraintValues)));
476 }
477 
478 
480 {
481  scalar objectiveValue(Zero);
482  for (adjointSolverManager& adjSolvManager : adjointSolvManagers_)
483  {
484  const scalar opWeight = adjSolvManager.operatingPointWeight();
485  objectiveValue += opWeight*adjSolvManager.objectiveValue();
486  }
487  return objectiveValue;
488 }
489 
492 {
493  updateMethod_->setObjectiveValueOld(computeObjectiveValue());
494 }
495 
496 
498 {
499  // Compute new objective and constraint values and update the ones
500  // in updateMethod
501  scalar objectiveValue(Zero);
502  DynamicList<scalar> constraintValues;
503 
504  for (adjointSolverManager& adjSolvManager : adjointSolvManagers_)
505  {
506  const scalar opWeight = adjSolvManager.operatingPointWeight();
507 
508  objectiveValue += opWeight*adjSolvManager.objectiveValue();
509  constraintValues.push_back(adjSolvManager.constraintValues());
510  }
511 
512  // Add constraints directly imposed to the design variables
513  tmp<scalarField> designVarsConstValues = designVars_->constraintValues();
514  if (designVarsConstValues)
515  {
516  constraintValues.push_back(designVarsConstValues());
517  }
518  updateMethod_->setObjectiveValue(objectiveValue);
519  updateMethod_->setConstraintValues
520  (scalarField(std::move(constraintValues)));
521 
522  return updateMethod_->computeMeritFunction();
523 }
524 
525 
527 {
528  return updateMethod_->meritFunctionDirectionalDerivative();
529 }
530 
531 
533 (
534  const scalarField& oldCorrection
535 )
536 {
537  updateMethod_->updateOldCorrection(oldCorrection);
538 }
539 
540 
542 {
543  updateMethod_->writeAuxiliaryData();
544  designVars_->writeDesignVars();
545  writeToCostFile();
546 }
547 
548 
549 void Foam::designVariablesUpdate::postUpdate(const scalarField& oldCorrection)
550 {
551  updateOldCorrection(oldCorrection);
552  write();
553  designVars_->evolveNumber();
554  if (lineSearch_.valid())
555  {
556  lineSearch_().postUpdate();
557  // Append additional empty line at the end of the instantaneous
558  // objective file to indicate the end of the block corresponding to
559  // this optimisation cycle
560  for (adjointSolverManager& am : adjointSolvManagers_)
561  {
562  for (adjointSolver& as : am.adjointSolvers())
563  {
564  PtrList<objective>& objectives =
565  as.getObjectiveManager().getObjectiveFunctions();
566  for (objective& obj : objectives)
567  {
568  obj.writeInstantaneousSeparator();
569  }
570  }
571  }
572  }
573  if (convergenceCriteriaDefined_)
574  {
575  checkConvergence(oldCorrection);
576  }
577 }
578 
579 
580 // ************************************************************************* //
PtrList< adjointSolverManager > & adjointSolvManagers_
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix by subtracting the matrix multiplied by the current fi...
A class encapsulating functionality neccessary to perform an optimisation loop, such as updating the ...
void writeCPUcostHeader()
Write CPU cost header.
const Type & value() const noexcept
Return const reference to value.
dictionary dict
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1)
uint8_t direction
Definition: direction.H:46
void postUpdate(const scalarField &oldCorrection)
Steps to be executed after the susccessfull update of the design varibles, i.e. the last step of line...
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
label nAdjointSolvers() const
Get total number of adjoint solvers.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
bool convergenceCriteriaDefined_
Is at least a single convergence criterion defined.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:1049
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.
static unsigned int defaultPrecision() noexcept
Return the default precision.
Definition: IOstream.H:423
Base class for adjoint solvers.
Definition: adjointSolver.H:53
constexpr char tab
The tab &#39;\t&#39; character(0x09)
Definition: Ostream.H:49
virtual dimensionedScalar endTime() const
Return end time.
Definition: Time.C:795
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
autoPtr< designVariables > & designVars_
label nConstraints(PtrList< adjointSolverManager > &adjointSolverManagers) const
Get the number of adjoint solvers that correspond to constraints.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
void reset(T *p=nullptr) noexcept
Delete managed object and set to new given pointer.
Definition: autoPtrI.H:37
dimensionedScalar pos(const dimensionedScalar &ds)
autoPtr< scalar > objectiveThreshold_
The relative update of the objective value w.r.t. to its last value should be smaller than this value...
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
tmp< scalarField > computeDirection()
Compute update direction.
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
dynamicFvMesh & mesh
static void exit(int errNo=1)
Shutdown (finalize) MPI as required and exit program with errNo.
Definition: UPstream.C:55
scalar feasibilityThreshold_
In case of a constrained optimisation, the sum of positive constraints should be lower than this valu...
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
void write()
Write useful quantities to files.
scalar computeObjectiveValue()
Sum objective values from all adjointSolverManagers.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
#define DebugInfo
Report an information message using Foam::Info.
Istream and Ostream manipulators taking arguments.
Abstract base class for objective functions. No point in making this runTime selectable since its chi...
Definition: objective.H:55
rDeltaT ref()
Type gMax(const FieldField< Field, Type > &f)
defineTypeNameAndDebug(combustionModel, 0)
scalar meritFunctionDirectionalDerivative()
Derivative of the merit function.
scalar computeMeritFunction()
Compute the merit function of the optimisation problem.
void push_back(const T &val)
Copy append an element to the end of this list.
Definition: DynamicListI.H:555
void writeToCostFile(bool zeroAdjointSolns=false)
Write to cost file.
Class for managing adjoint solvers, which may be more than one per operating point.
autoPtr< scalar > designVarsThreshold_
The maximum of the correction/designVariables values must be lower that this threshold to consider th...
#define WarningInFunction
Report a warning using Foam::Warning.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
void updateOldCorrection(const scalarField &)
Update old correction. Needed for quasi-Newton Methods.
void update()
Update design variables.
messageStream Info
Information stream (stdout output on master, null elsewhere)
label n
Omanip< int > setw(const int i)
Definition: IOmanip.H:199
Macros to ease declaration of run-time selection tables.
autoPtr< updateMethod > updateMethod_
Method to update the design variables based on the provided sensitivity derivatives.
A class for managing temporary objects.
Definition: HashPtrTable.H:50
void updateGradientsAndValues()
Compute cumulative objective and constraint gradients.
PtrList< adjointSolverManager > & adjointSolverManagers
Definition: createFields.H:8
void setOldObjectiveValue()
Set the old objective value known by the updateMethod.
Namespace for OpenFOAM.
void checkConvergence(const scalarField &oldCorrection)
Check if the optimisation loop has converged based on the provided criteria.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127