Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
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
29Class
30 Foam::designVariablesUpdate
31
32Description
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
43SourceFiles
44 designVariablesUpdate.C
45
46\*---------------------------------------------------------------------------*/
47
48#ifndef designVariablesUpdate_H
49#define designVariablesUpdate_H
50
52#include "designVariables.H"
53#include "updateMethod.H"
54#include "lineSearch.H"
55
56// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
57
58namespace Foam
59{
61/*---------------------------------------------------------------------------*\
62 Class designVariablesUpdate Declaration
63\*---------------------------------------------------------------------------*/
64
65class designVariablesUpdate
67protected:
69 // Protected data
70
72 const dictionary dict_;
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
88 //- Primal solutions per optimisation cycle
90
91 //- Adjoint solutions per optimisation cycle
93
94 //- Primal evaluations performed so far
96
97 //- Adjoint evaluations performed so far
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)
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
142private:
143
144 // Private Member Functions
145
146 //- No copy construct
147 designVariablesUpdate(const designVariablesUpdate&) = delete;
148
149 //- No copy assignment
150 void operator=(const designVariablesUpdate&) = delete;
151
152
153public:
154
155 //- Runtime type information
156 TypeName("designVariablesUpdate");
157
158
159 // Constructors
160
161 //- Construct from components
162 designVariablesUpdate
163 (
164 fvMesh& mesh,
165 const dictionary& dict,
167 autoPtr<designVariables>& designVars
168 );
169
170 //- Destructor
171 virtual ~designVariablesUpdate() = default;
172
173
174 // Member Functions
175
176 //- Update design variables
177 void update();
178
179 //- Update design variables based on a given correction
180 void update(const scalarField& correction);
181
182 //- Compute update direction
184
185 //- Compute cumulative objective and constraint gradients
187
188 //- Sum objective values from all adjointSolverManagers
190
191 //- Set the old objective value known by the updateMethod
192 // Used to check convergence of the optimisation loop
194
195 //- Compute the merit function of the optimisation problem.
196 // Could be different than the objective function in case of
197 // constraint optimisation
198 scalar computeMeritFunction();
199
200 //- Derivative of the merit function
202
203 //- Update old correction. Needed for quasi-Newton Methods
204 void updateOldCorrection(const scalarField&);
205
206 //- Write useful quantities to files
207 void write();
209 //- Steps to be executed after the susccessfull update of the design
210 //- varibles, i.e. the last step of line search or the simple update
211 //- in the fixedStep approach
212 void postUpdate(const scalarField& oldCorrection);
213
214 //- Check if the optimisation loop has converged based on the provided
215 //- criteria
216 // May terminate the program
217 void checkConvergence();
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// ************************************************************************* //
Output to file stream as an OSstream, normally using std::ofstream for the actual output.
Definition OFstream.H:75
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
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.
autoPtr< lineSearch > & getLineSearch()
Get a reference to the line search object.
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.
virtual ~designVariablesUpdate()=default
Destructor.
label nPrimalsPerCycle_
Primal solutions per optimisation cycle.
autoPtr< designVariables > & getDesignVariables()
Get access to design variables.
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_
TypeName("designVariablesUpdate")
Runtime type information.
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
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
A class for managing temporary objects.
Definition tmp.H:75
dynamicFvMesh & mesh
Namespace for OpenFOAM.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix by subtracting the matrix multiplied by the current fi...
dictionary dict
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition typeInfo.H:68
PtrList< adjointSolverManager > & adjointSolverManagers
Definition createFields.H:8