designVariablesUpdate.H
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 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 Class
30  Foam::designVariablesUpdate
31 
32 Description
33  A class encapsulating functionality neccessary to perform an optimisation
34  loop, such as updating the design variables, checking the
35  sufficient reduction/adhetion of objective and constraint values in each
36  optimisation cycle by performing a line-search and checking the overall
37  convergence of the optimisation loop, if the corresponding criteria are
38  provided.
39 
40  Kept separate from optimisationManager to isolate functionality required
41  only when the update of the design variables is performed.
42 
43 SourceFiles
44  designVariablesUpdate.C
45 
46 \*---------------------------------------------------------------------------*/
47 
48 #ifndef designVariablesUpdate_H
49 #define designVariablesUpdate_H
50 
51 #include "adjointSolverManager.H"
52 #include "designVariables.H"
53 #include "updateMethod.H"
54 #include "lineSearch.H"
55 
56 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
57 
58 namespace Foam
59 {
60 
61 /*---------------------------------------------------------------------------*\
62  Class designVariablesUpdate Declaration
63 \*---------------------------------------------------------------------------*/
64 
66 {
67 protected:
68 
69  // Protected data
70 
71  fvMesh& mesh_;
72  const dictionary dict_;
75 
76  //- Method to update the design variables based on the provided
77  //- sensitivity derivatives
79 
80  //- Line search mechanism to approximate the update step length
82 
83  // Variables related to the computation of the CPU cost
84 
85  //- Output file
87 
88  //- Primal solutions per optimisation cycle
89  label nPrimalsPerCycle_;
90 
91  //- Adjoint solutions per optimisation cycle
93 
94  //- Primal evaluations performed so far
95  label nPrimalSolutions_;
96 
97  //- Adjoint evaluations performed so far
98  label nAdjointSolutions_;
99 
100  //- CPU cost (in seconds)
101  scalar CPUcost_;
103  // Convergence criteria
104 
105  //- The maximum of the correction/designVariables values
106  //- must be lower that this threshold to consider the run converged
108 
109  //- The relative update of the objective value w.r.t. to its last
110  //- value should be smaller than this value to considered the run
111  //- converged
113 
114  //- Is at least a single convergence criterion defined
116 
117  //- In case of a constrained optimisation, the sum of positive
118  //- constraints should be lower than this value to consider the
119  //- run converged (i.e. this tolerates some deviation from
120  //- satisfying all constraints)
121  scalar feasibilityThreshold_;
122 
123 
124  // Protected Member Functions
125 
126  //- Get the number of adjoint solvers that correspond to constraints
128  (
130  ) const;
131 
132  //- Get total number of adjoint solvers
133  label nAdjointSolvers() const;
134 
135  //- Write CPU cost header
136  void writeCPUcostHeader();
137 
138  //- Write to cost file
139  void writeToCostFile(bool zeroAdjointSolns = false);
141  //- Check if the optimisation loop has converged based on the provided
142  //- criteria
143  // May terminate the program
144  void checkConvergence(const scalarField& oldCorrection);
145 
146 
147 private:
148 
149  // Private Member Functions
150 
151  //- No copy construct
153 
154  //- No copy assignment
155  void operator=(const designVariablesUpdate&) = delete;
156 
157 
158 public:
159 
160  //- Runtime type information
161  TypeName("designVariablesUpdate");
162 
163 
164  // Constructors
165 
166  //- Construct from components
168  (
169  fvMesh& mesh,
170  const dictionary& dict,
172  autoPtr<designVariables>& designVars
173  );
174 
175  //- Destructor
176  virtual ~designVariablesUpdate() = default;
177 
178 
179  // Member Functions
180 
181  //- Update design variables
182  void update();
183 
184  //- Update design variables based on a given correction
185  void update(const scalarField& correction);
186 
187  //- Compute update direction
189 
190  //- Compute cumulative objective and constraint gradients
192 
193  //- Sum objective values from all adjointSolverManagers
194  scalar computeObjectiveValue();
195 
196  //- Set the old objective value known by the updateMethod
197  // Used to check convergence of the optimisation loop
198  void setOldObjectiveValue();
199 
200  //- Compute the merit function of the optimisation problem.
201  // Could be different than the objective function in case of
202  // constraint optimisation
203  scalar computeMeritFunction();
204 
205  //- Derivative of the merit function
207 
208  //- Update old correction. Needed for quasi-Newton Methods
209  void updateOldCorrection(const scalarField&);
210 
211  //- Write useful quantities to files
212  void write();
213 
214  //- Steps to be executed after the susccessfull update of the design
215  //- varibles, i.e. the last step of line search or the simple update
216  //- in the fixedStep approach
217  void postUpdate(const scalarField& oldCorrection);
218 
219  //- Get access to design variables
221  {
222  return designVars_;
223  }
224 
225  //- Get a reference to the line search object
227  {
228  return lineSearch_;
229  }
230 };
231 
232 
233 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
234 
235 } // End namespace Foam
236 
237 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
238 
239 #endif
240 
241 // ************************************************************************* //
PtrList< adjointSolverManager > & adjointSolvManagers_
OFstream CPUcostFile_
Output file.
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.
dictionary dict
label nPrimalSolutions_
Primal evaluations performed so far.
void postUpdate(const scalarField &oldCorrection)
Steps to be executed after the susccessfull update of the design varibles, i.e. the last step of line...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
label nAdjointSolvers() const
Get total number of adjoint solvers.
label nAdjointSolutions_
Adjoint evaluations performed so far.
Output to file stream, using an OSstream.
Definition: OFstream.H:49
bool convergenceCriteriaDefined_
Is at least a single convergence criterion defined.
label nAdjointsPerCycle_
Adjoint solutions per optimisation cycle.
autoPtr< designVariables > & designVars_
label nConstraints(PtrList< adjointSolverManager > &adjointSolverManagers) const
Get the number of adjoint solvers that correspond to constraints.
autoPtr< lineSearch > & getLineSearch()
Get a reference to the line search object.
autoPtr< scalar > objectiveThreshold_
The relative update of the objective value w.r.t. to its last value should be smaller than this value...
tmp< scalarField > computeDirection()
Compute update direction.
dynamicFvMesh & mesh
scalar feasibilityThreshold_
In case of a constrained optimisation, the sum of positive constraints should be lower than this valu...
void write()
Write useful quantities to files.
scalar CPUcost_
CPU cost (in seconds)
scalar computeObjectiveValue()
Sum objective values from all adjointSolverManagers.
scalar meritFunctionDirectionalDerivative()
Derivative of the merit function.
scalar computeMeritFunction()
Compute the merit function of the optimisation problem.
void writeToCostFile(bool zeroAdjointSolns=false)
Write to cost file.
autoPtr< scalar > designVarsThreshold_
The maximum of the correction/designVariables values must be lower that this threshold to consider th...
TypeName("designVariablesUpdate")
Runtime type information.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
void updateOldCorrection(const scalarField &)
Update old correction. Needed for quasi-Newton Methods.
void update()
Update design variables.
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
autoPtr< designVariables > & getDesignVariables()
Get access to design variables.
virtual ~designVariablesUpdate()=default
Destructor.
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.
autoPtr< lineSearch > lineSearch_
Line search mechanism to approximate the update step length.
Namespace for OpenFOAM.
label nPrimalsPerCycle_
Primal solutions per optimisation cycle.
void checkConvergence(const scalarField &oldCorrection)
Check if the optimisation loop has converged based on the provided criteria.