Loading...
Searching...
No Matches
lineSearch.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-2022 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\*---------------------------------------------------------------------------*/
29
30#include "lineSearch.H"
31#include "Time.H"
32
33// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34
35namespace Foam
36{
37 defineTypeNameAndDebug(lineSearch, 0);
39}
40
41
42// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
43
46 return dict_.optionalSubDict(type() + "Coeffs");
47}
48
49
50// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
51
52Foam::lineSearch::lineSearch
53(
54 const dictionary& dict,
55 const Time& time,
56 updateMethod& UpdateMethod
57)
58:
59 dict_(dict),
60 lineSearchDict_
61 (
63 (
64 "lineSearch",
65 time.timeName(),
66 "uniform",
67 time,
68 IOobject::READ_IF_PRESENT,
69 IOobject::NO_WRITE,
70 IOobject::NO_REGISTER
71 )
72 ),
73 directionalDeriv_(Zero),
74 direction_(0),
75 oldMeritValue_(Zero),
76 newMeritValue_(Zero),
77 prevMeritDeriv_
78 (
79 lineSearchDict_.getOrDefault<scalar>("prevMeritDeriv", Zero)
80 ),
81 initialStep_(dict.getOrDefault<scalar>("initialStep", 1.)),
82 minStep_(dict.getOrDefault<scalar>("minStep", 0.3)),
83 step_(Zero),
84 iter_(lineSearchDict_.getOrDefault<label>("iter", 0)),
85 innerIter_(0),
86 maxIters_(dict.getOrDefault<label>("maxIters", 4)),
87 extrapolateInitialStep_
88 (
89 dict.getOrDefault<bool>
90 (
91 "extrapolateInitialStep",
92 false
93 )
94 ),
96 updateMethod_(UpdateMethod)
97{}
98
99
100// * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
101
103(
104 const dictionary& dict,
105 const Time& time,
106 updateMethod& UpdateMethod
107)
108{
109 autoPtr<lineSearch> lineSrch(nullptr);
110
111 const word modelType(dict.getOrDefault<word>("type", "none"));
112
113 Info<< "lineSearch type : " << modelType << endl;
114
115 if (modelType != "none")
116 {
117 auto* ctorPtr = dictionaryConstructorTable(modelType);
118
119 if (!ctorPtr)
120 {
122 (
123 dict,
124 "lineSearch",
125 modelType,
126 *dictionaryConstructorTablePtr_
127 ) << exit(FatalIOError);
128 }
129
130 lineSrch.reset((ctorPtr(dict, time, UpdateMethod)).ptr());
131 }
132 else
133 {
134 Info<< "No line search method specified. "
135 << "Proceeding with constant step" << endl;
136 }
138 return lineSrch;
139}
140
141
142// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * //
143
144void Foam::lineSearch::setDeriv(const scalar deriv)
145{
146 directionalDeriv_ = deriv;
147 stepUpdate_->setDeriv(deriv);
148}
149
151void Foam::lineSearch::setNewDeriv(const scalar deriv)
152{
153 // Does nothing in base
154}
155
156
158{
159 newMeritValue_ = value;
160 stepUpdate_->setNewMeritValue(value);
161}
162
163
165{
166 oldMeritValue_ = value;
167 stepUpdate_->setOldMeritValue(value);
168}
169
170
172{
173 innerIter_ = 0;
174 if (extrapolateInitialStep_ && iter_ != 0)
175 {
176 // step_ = 2*(oldMeritValue_-prevMeritValue_)/directionalDeriv_;
177 // Interpolate in order to get same improvement with the previous
178 // optimisation cycle
179 step_ =
180 clamp(step_*prevMeritDeriv_/directionalDeriv_, minStep_, scalar(1));
181 Info<< "\n------- Computing initial step-------" << endl;
182 Info<< "old dphi(0) " << prevMeritDeriv_ << endl;
183 Info<< "dphi(0) " << directionalDeriv_ << endl;
184 Info<< "Setting initial step value " << step_ << endl << endl;
185 }
186 else
187 {
189 }
190}
191
193void Foam::lineSearch::updateStep(const scalar newStep)
194{
195 step_ = newStep;
196}
197
200{
201 correction *= step_;
202}
203
204
206{
207 const bool isRunning = innerIter_ < maxIters_;
208
209 if (isRunning)
210 {
212 }
213
214 return isRunning;
215}
216
219{
220 return false;
221}
222
225{
226 this->operator++();
227}
228
229
231{
232 ++iter_;
233 prevMeritDeriv_ = directionalDeriv_;
234 lineSearchDict_.add<scalar>("prevMeritDeriv", prevMeritDeriv_, true);
235 lineSearchDict_.add<label>("iter", iter_, true);
236 if (lineSearchDict_.time().writeTime())
237 {
238 lineSearchDict_.regIOobject::writeObject
239 (
241 true
242 );
243 }
244
245 return *this;
246}
247
248
250{
251 return operator++();
252}
253
254
255// ************************************************************************* //
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
A simple container for options an IOstream can normally have.
@ ASCII
"ascii" (normal default)
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
void reset(T *p=nullptr) noexcept
Delete managed object and set to new given pointer.
Definition autoPtrI.H:37
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Abstract base class for line search methods.
Definition lineSearch.H:54
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
static autoPtr< lineSearch > New(const dictionary &dict, const Time &time, updateMethod &UpdateMethod)
Return a reference to the selected turbulence model.
Definition lineSearch.C:96
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
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
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_
Inner line search iteration.
Definition lineSearch.H:117
Abstract base class for step update methods used in line search.
Definition stepUpdate.H:51
Abstract base class for optimisation methods.
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition error.H:637
word timeName
Definition getTimeIndex.H:3
Namespace for OpenFOAM.
messageStream Info
Information stream (stdout output on master, null elsewhere).
dimensionSet clamp(const dimensionSet &a, const dimensionSet &range)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition POSIX.C:801
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix by subtracting the matrix multiplied by the current fi...
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.
dictionary dict