Loading...
Searching...
No Matches
updateMethod.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
28Class
29 Foam::updateMethod
30
31Description
32 Abstract base class for optimisation methods
33
34SourceFiles
35 updateMethod.C
36
37\*---------------------------------------------------------------------------*/
38
39#ifndef updateMethod_H
40#define updateMethod_H
41
42#include "localIOdictionary.H"
43#include "designVariables.H"
45
46// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47
48namespace Foam
49{
51/*---------------------------------------------------------------------------*\
52 Class updateMethod Declaration
53\*---------------------------------------------------------------------------*/
54
55class updateMethod
56:
59protected:
61 // Protected data
62
63 const fvMesh& mesh_;
64
66
67 //- Design variables
69
70 //- Number of constraints
71 label nConstraints_;
72
73 //- Map to active design variables
76 //- Derivatives of the objective functions
78
79 //- Derivatives of the constraints
81
82 //- Objective value
83 scalar objectiveValue_;
84
85 //- Old objective value
86 // Used for convergence checking
88
89 //- Constraint values
91
92 //- Design variables correction
94
95 //- Cumulative design variables correction throughout the optimisation
96 //- loop
98
99 //- Step multiplying the correction
100 scalar eta_;
101
102 //- Optimisation cycle count
103 label counter_;
104
105 //- Is initially set?
106 bool initialEtaSet_;
108 //- Folder storing the corrections to file
109 // For some optimisation methods with a very high number of
110 // design variables (e.g. topology), it doesn't make much sense
111 // to write all updates in the updateMethodDict. Hence, a
112 // separate file is used to write the corrections, in case they are
113 // needed for post-processing
115
116 //- Whether to use gSum or sum in the inner products
117 bool globalSum_;
119 // Scalar -- matrix multiplications
121 (
122 const scalarField&,
124 );
125
127 (
129 const scalarField&
130 );
131
133 (
134 const scalarField&,
135 const scalarField&
136 );
137
140 //- Compute either global or local sum, based on globalSum flag
141 scalar globalSum(const scalarField& field);
142
143 //- Compute either global or local sum, based on globalSum flag
145
146 //- Compute either global or local sum, based on globalSum flag
147 label globalSum(const label);
148
149 //- Helper function to either read a scalarField of certain size
150 //- from a dictionary, or construct a zero field
151 tmp<scalarField> readOrZeroField(const word& name, const label size);
152
153
154private:
155
156 // Private Member Functions
157
158 //- No copy construct
159 updateMethod(const updateMethod&) = delete;
160
161 //- No copy assignment
162 void operator=(const updateMethod&) = delete;
163
164
165public:
166
167 //- Runtime type information
168 TypeName("updateMethod");
169
170
171 // Declare run-time constructor selection table
172
174 (
175 autoPtr,
176 updateMethod,
178 (
179 const fvMesh& mesh,
180 const dictionary& dict,
181 autoPtr<designVariables>& designVars,
182 const label nConstraints,
183 const word& type
184 ),
185 (mesh, dict, designVars, nConstraints, type)
186 );
187
188
189 // Constructors
190
191 //- Construct from components
192 updateMethod
193 (
194 const fvMesh& mesh,
195 const dictionary& dict,
196 autoPtr<designVariables>& designVars,
197 const label nConstraints,
198 const word& type
199 );
200
201
202 // Selectors
203
204 //- Return a reference to the selected turbulence model
206 (
207 const fvMesh& mesh,
208 const dictionary& dict,
210 const label nConstraints
211 );
212
213
214 //- Destructor
215 virtual ~updateMethod() = default;
216
217
218 // Member Functions
219
220 //- Return optional dictionary with parameters specific to each method
221 dictionary coeffsDict(const word& type) const;
222
223 //- Set objective derivative
224 void setObjectiveDeriv(const scalarField& derivs);
225
226 //- Set constraints derivative
227 void setConstraintDeriv(const PtrList<scalarField>& derivs);
228
229 //- Set objective value
230 void setObjectiveValue(const scalar value);
231
232 //- Set old objective value
233 void setObjectiveValueOld(const scalar value);
234
235 //- Set values of constraints
236 void setConstraintValues(const scalarField& values);
237
238 //- Get objective value
239 scalar getObjectiveValue() const;
240
241 //- Get old objective value
243
244 //- Get values of constraints
245 const scalarField& getConstraintValues() const;
246
247 //- Get optimisation cycle
248 label getCycle() const;
249
250 //- Set step for optimisation methods
251 void setStep(const scalar eta);
252
253 //- Multiply step
254 void modifyStep(const scalar multiplier);
255
256 //- Set globalSum variable.
257 // Should be set by the optimisationManager owning the updateMethod
258 void setGlobalSum(const bool useGlobalSum);
259
260 //- Set the number of constraints
261 void setConstaintsNumber(const label nConstraints);
263 //- Get the number of constraints
264 label nConstraints() const;
265
266 //- Return the correction of the design variables
267 virtual void computeCorrection() = 0;
268
269 //- Return the correction of the design variables
270 //const scalarField& returnCorrection() const;
271
272 //- Return the correction of the design variables
274
275 void writeCorrection();
276
277 //- Compute merit function. Could be different than the objective
278 //- in the presence of constraints
279 virtual scalar computeMeritFunction();
280
281 //- Directional derivative of the merit function, in the direction of
282 //- the correction. Could be different than the objective directional
283 //- derivative in the presence of constraints
284 virtual scalar meritFunctionDirectionalDerivative();
285
286 //- Return whether initial eta was set
287 bool& initialEtaSet();
288
289 //- Update old correction. Useful for quasi-newton methods coupled with
290 //- line search
291 virtual void updateOldCorrection(const scalarField& oldCorrection);
292
293 //- Write continuation data under uniform
294 virtual bool writeData(Ostream& os) const;
295
296 //- Write non-continuation data, usually under the optimisation folder
297 virtual bool writeAuxiliaryData();
298};
299
300
301// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
302
303} // End namespace Foam
304
305// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
306
307#endif
308
309// ************************************************************************* //
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
label size() const noexcept
The number of elements in list.
Definition DLListBase.H:194
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
A templated (N x N) square matrix of objects of <Type>, containing N*N elements, derived from Matrix.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
const word & name() const
Name function is needed to disambiguate those inherited from regIOobject and dictionary.
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
localIOdictionary(const IOobject &io, const dictionary *fallback=nullptr)
Construct given an IOobject and optional fallback dictionary content.
A class for managing temporary objects.
Definition tmp.H:75
declareRunTimeSelectionTable(autoPtr, updateMethod, dictionary,(const fvMesh &mesh, const dictionary &dict, autoPtr< designVariables > &designVars, const label nConstraints, const word &type),(mesh, dict, designVars, nConstraints, type))
void setConstraintValues(const scalarField &values)
Set values of constraints.
scalarField objectiveDerivatives_
Derivatives of the objective functions.
label nConstraints_
Number of constraints.
virtual scalar meritFunctionDirectionalDerivative()
Directional derivative of the merit function, in the direction of the correction. Could be different ...
scalarField correction_
Design variables correction.
bool globalSum_
Whether to use gSum or sum in the inner products.
const fvMesh & mesh_
dictionary coeffsDict(const word &type) const
Return optional dictionary with parameters specific to each method.
scalar globalSum(const scalarField &field)
Compute either global or local sum, based on globalSum flag.
scalarField cValues_
Constraint values.
scalar objectiveValue_
Objective value.
virtual bool writeAuxiliaryData()
Write non-continuation data, usually under the optimisation folder.
PtrList< scalarField > constraintDerivatives_
Derivatives of the constraints.
void setStep(const scalar eta)
Set step for optimisation methods.
virtual bool writeData(Ostream &os) const
Write continuation data under uniform.
const scalarField leftMult(const scalarField &, const SquareMatrix< scalar > &)
tmp< scalarField > readOrZeroField(const word &name, const label size)
Helper function to either read a scalarField of certain size from a dictionary, or construct a zero f...
scalarField & returnCorrection()
Return the correction of the design variables.
TypeName("updateMethod")
Runtime type information.
SquareMatrix< scalar > inv(SquareMatrix< scalar > A)
const labelList & activeDesignVars_
Map to active design variables.
scalar getObjectiveValue() const
Get objective value.
label getCycle() const
Get optimisation cycle.
bool initialEtaSet_
Is initially set?
const scalarField rightMult(const SquareMatrix< scalar > &, const scalarField &)
void setObjectiveValue(const scalar value)
Set objective value.
autoPtr< scalar > objectiveValueOld_
Old objective value.
static autoPtr< updateMethod > New(const fvMesh &mesh, const dictionary &dict, autoPtr< designVariables > &designVars, const label nConstraints)
Return a reference to the selected turbulence model.
const autoPtr< scalar > & getObjectiveValueOld() const
Get old objective value.
void setObjectiveValueOld(const scalar value)
Set old objective value.
label counter_
Optimisation cycle count.
void setConstaintsNumber(const label nConstraints)
Set the number of constraints.
void setConstraintDeriv(const PtrList< scalarField > &derivs)
Set constraints derivative.
virtual scalar computeMeritFunction()
Compute merit function. Could be different than the objective in the presence of constraints.
virtual ~updateMethod()=default
Destructor.
const dictionary dict_
const scalarField & getConstraintValues() const
Get values of constraints.
label nConstraints() const
Get the number of constraints.
void modifyStep(const scalar multiplier)
Multiply step.
scalarField cumulativeCorrection_
Cumulative design variables correction throughout the optimisation loop.
virtual void updateOldCorrection(const scalarField &oldCorrection)
Update old correction. Useful for quasi-newton methods coupled with line search.
scalar eta_
Step multiplying the correction.
bool & initialEtaSet()
Return whether initial eta was set.
autoPtr< designVariables > & designVars_
Design variables.
word correctionFolder_
Folder storing the corrections to file.
SquareMatrix< scalar > outerProd(const scalarField &, const scalarField &)
void setObjectiveDeriv(const scalarField &derivs)
Set objective derivative.
virtual void computeCorrection()=0
Return the correction of the design variables.
void setGlobalSum(const bool useGlobalSum)
Set globalSum variable.
A class for handling words, derived from Foam::string.
Definition word.H:66
rDeltaTY field()
dynamicFvMesh & mesh
OBJstream os(runTime.globalPath()/outputName)
Namespace for OpenFOAM.
List< label > labelList
A List of labels.
Definition List.H:62
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Macros to ease declaration of run-time selection tables.
#define declareRunTimeSelectionTable(ptrWrapper, baseType, argNames, argList, parList)
Declare a run-time selection (variables and adder classes).
dictionary dict
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition typeInfo.H:68