Loading...
Searching...
No Matches
lineSearch.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::lineSearch
31
32Description
33 Abstract base class for line search methods
34
35SourceFiles
36 lineSearch.C
37
38\*---------------------------------------------------------------------------*/
39
40#ifndef lineSearch_H
41#define lineSearch_H
42
44#include "IOdictionary.H"
45#include "scalarField.H"
46#include "stepUpdate.H"
47#include "updateMethod.H"
48
49// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
50
51namespace Foam
52{
54/*---------------------------------------------------------------------------*\
55 Class lineSearch Declaration
56\*---------------------------------------------------------------------------*/
57
58class lineSearch
59{
60protected:
61
62 // Protected data
63
64 //- Subdict within updateMethod
65 const dictionary dict_;
66
67 //- IOdictionary under time/uniform for continuation
69
70 //- Directional derivative of the merit function
71 scalar directionalDeriv_;
73 //- Update direction
75
76 //- Old merit value from this opt cycle
78
79 //- New merit value from this opt cycle
80 scalar newMeritValue_;
81
82 //- Merit directional deriv from the previous opt cycle
83 scalar prevMeritDeriv_;
84
85 //- Correction multiplier at the first step of line search
86 scalar initialStep_;
88 //- Minimum allowed correction multiplier
89 scalar minStep_;
90
91 //- Correction multiplier
92 scalar step_;
93
94 //- Optimisation cycle
95 label iter_;
96
97 //- Inner line search iteration
98 label innerIter_;
99
100 //- Maximum line search iterations
101 label maxIters_;
103 //- Whether to extrapolate the correction multiplier for
104 //- this optimisation cycle based on the previous ones.
105 // Useful for non-quasi Newton methods
108 //- Mechanism to update method if line search conditions are not set
110
111 //- Reference to the update method related to the line search
113
114
115 // Protected Member Functions
116
117 //- Optional coeffs dict
118 const dictionary& coeffsDict();
119
120
121private:
123 // Private Member Functions
124
125 //- No copy construct
126 lineSearch(const lineSearch&) = delete;
127
128 //- No copy assignment
129 void operator=(const lineSearch&) = delete;
131
132public:
133
134 //- Runtime type information
135 TypeName("lineSearch");
136
137
138 // Declare run-time constructor selection table
139
141 (
142 autoPtr,
143 lineSearch,
145 (
146 const dictionary& dict,
147 const Time& time,
148 updateMethod& UpdateMethod
149 ),
150 (dict, time, UpdateMethod)
151 );
152
153
154 // Constructors
155
156 //- Construct from components
157 lineSearch
158 (
159 const dictionary& dict,
160 const Time& time,
161 updateMethod& UpdateMethod
162 );
163
164 // Selectors
165
166 //- Return a reference to the selected turbulence model
168 (
169 const dictionary& dict,
170 const Time& time,
171 updateMethod& UpdateMethod
172 );
173
174
175 //- Destructor
176 virtual ~lineSearch() = default;
177
178
179 // Member Functions
180
181 //- Set directional derivative
182 virtual void setDeriv(const scalar deriv);
183
184 //- Set new directional derivative.
185 // Does nothing in base. Only used by methods that require the
186 // gradient information to be computed at the new point
187 virtual void setNewDeriv(const scalar deriv);
188
189 //- Set direction
190 inline void setDirection(const scalarField& direction)
191 {
193 }
194
195 //- Set new objective value
196 void setNewMeritValue(const scalar value);
197
198 //- Set old objective value
199 void setOldMeritValue(const scalar value);
200
201 //- Reset step to initial value
202 virtual void reset();
203
204 //- Get inner line search iteration
205 inline label innerIter() const
206 {
207 return innerIter_;
208 }
209
210 //- Get max number of iterations
211 inline label maxIters() const
212 {
213 return maxIters_;
214 }
215
216 //- Get current step
217 inline scalar step() const
219 return step_;
220 }
221
222 //- Return the correction of the design variables
223 virtual bool converged() = 0;
224
225 //- Update the line search step based on the specific line search
226 //- strategy, e.g. bisection, quadratic fit, etc.
227 virtual void updateStep() = 0;
228
229 //- Update the step using the supplied value
230 virtual void updateStep(const scalar newStep);
231
232 //- Update the correction.
233 // Multiplies with step in base
235
236 //- Return true if lineSearch should continue and if so increment inner
237 // iteration
238 virtual bool loop();
240 //- Does line search need to update the gradient?
241 virtual bool computeGradient() const;
242
243 //- Execute steps at the end of the line search iterations
244 virtual void postUpdate();
245
246
247 // Member operators
248
249 //- Increment iteration number and store merit value corresponding to
250 //- the previous optimisation cycle
251 virtual lineSearch& operator++();
252
253 //- Postfix increment. Necessary for compilation
254 virtual lineSearch& operator++(int);
255};
256
257
258// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
259
260} // End namespace Foam
261
262// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
263
264#endif
265
266// ************************************************************************* //
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
const dictionary & coeffsDict()
Optional coeffs dict.
Definition lineSearch.C:37
virtual void setDeriv(const scalar deriv)
Set directional derivative.
Definition lineSearch.C:137
scalar newMeritValue_
New merit value from this opt cycle.
Definition lineSearch.H:87
scalar oldMeritValue_
Old merit value from this opt cycle.
Definition lineSearch.H:82
bool extrapolateInitialStep_
Whether to extrapolate the correction multiplier for this optimisation cycle based on the previous on...
Definition lineSearch.H:130
virtual void updateCorrection(scalarField &correction)
Update the correction.
Definition lineSearch.C:192
TypeName("lineSearch")
Runtime type information.
static autoPtr< lineSearch > New(const dictionary &dict, const Time &time, updateMethod &UpdateMethod)
Return a reference to the selected turbulence model.
Definition lineSearch.C:96
scalar step() const
Get current step.
Definition lineSearch.H:278
virtual void postUpdate()
Execute steps at the end of the line search iterations.
Definition lineSearch.C:217
updateMethod & updateMethod_
Reference to the update method related to the line search.
Definition lineSearch.H:140
declareRunTimeSelectionTable(autoPtr, lineSearch, dictionary,(const dictionary &dict, const Time &time, updateMethod &UpdateMethod),(dict, time, UpdateMethod))
autoPtr< stepUpdate > stepUpdate_
Mechanism to update method if line search conditions are not set.
Definition lineSearch.H:135
void setOldMeritValue(const scalar value)
Set old objective value.
Definition lineSearch.C:157
void setNewMeritValue(const scalar value)
Set new objective value.
Definition lineSearch.C:150
virtual bool converged()=0
Return the correction of the design variables.
const dictionary dict_
Subdict within updateMethod.
Definition lineSearch.H:62
virtual void setNewDeriv(const scalar deriv)
Set new directional derivative.
Definition lineSearch.C:144
virtual bool computeGradient() const
Does line search need to update the gradient?
Definition lineSearch.C:211
virtual void updateStep()=0
Update the line search step based on the specific line search strategy, e.g. bisection,...
label iter_
Optimisation cycle.
Definition lineSearch.H:112
scalarField direction_
Update direction.
Definition lineSearch.H:77
scalar step_
Correction multiplier.
Definition lineSearch.H:107
scalar directionalDeriv_
Directional derivative of the merit function.
Definition lineSearch.H:72
virtual lineSearch & operator++()
Increment iteration number and store merit value corresponding to the previous optimisation cycle.
Definition lineSearch.C:223
scalar minStep_
Minimum allowed correction multiplier.
Definition lineSearch.H:102
label maxIters_
Maximum line search iterations.
Definition lineSearch.H:122
virtual void reset()
Reset step to initial value.
Definition lineSearch.C:164
IOdictionary lineSearchDict_
IOdictionary under time/uniform for continuation.
Definition lineSearch.H:67
scalar initialStep_
Correction multiplier at the first step of line search.
Definition lineSearch.H:97
virtual bool loop()
Return true if lineSearch should continue and if so increment inner.
Definition lineSearch.C:198
scalar prevMeritDeriv_
Merit directional deriv from the previous opt cycle.
Definition lineSearch.H:92
label innerIter() const
Get inner line search iteration.
Definition lineSearch.H:262
label innerIter_
Inner line search iteration.
Definition lineSearch.H:117
void setDirection(const scalarField &direction)
Set direction.
Definition lineSearch.H:239
label maxIters() const
Get max number of iterations.
Definition lineSearch.H:270
virtual ~lineSearch()=default
Destructor.
Abstract base class for optimisation methods.
Namespace for OpenFOAM.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
uint8_t direction
Definition direction.H:49
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix by subtracting the matrix multiplied by the current fi...
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