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 -------------------------------------------------------------------------------
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 
30 #include "lineSearch.H"
31 #include "Time.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
39 }
40 
41 
42 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
43 
45 {
46  return dict_.optionalSubDict(type() + "Coeffs");
47 }
48 
49 
50 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
51 
52 Foam::lineSearch::lineSearch
53 (
54  const dictionary& dict,
55  const Time& time,
56  updateMethod& UpdateMethod
57 )
58 :
59  dict_(dict),
60  lineSearchDict_
61  (
62  IOobject
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  ),
95  stepUpdate_(stepUpdate::New(dict)),
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 
144 void Foam::lineSearch::setDeriv(const scalar deriv)
145 {
146  directionalDeriv_ = deriv;
147  stepUpdate_->setDeriv(deriv);
148 }
149 
151 void Foam::lineSearch::setNewDeriv(const scalar deriv)
152 {
153  // Does nothing in base
154 }
155 
156 
157 void Foam::lineSearch::setNewMeritValue(const scalar value)
158 {
159  newMeritValue_ = value;
160  stepUpdate_->setNewMeritValue(value);
161 }
162 
163 
164 void Foam::lineSearch::setOldMeritValue(const scalar value)
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  {
188  step_ = initialStep_;
189  }
190 }
191 
193 void 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  {
211  ++innerIter_;
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 // ************************************************************************* //
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
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
virtual void updateStep()=0
Update the line search step based on the specific line search strategy, e.g. bisection, quadratic fit, etc.
const dictionary dict_
Subdict within updateMethod.
Definition: lineSearch.H:62
virtual void updateCorrection(scalarField &correction)
Update the correction.
Definition: lineSearch.C:192
const dictionary & coeffsDict()
Optional coeffs dict.
Definition: lineSearch.C:37
Abstract base class for optimisation methods.
Definition: updateMethod.H:50
"ascii" (normal default)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
virtual void postUpdate()
Execute steps at the end of the line search iterations.
Definition: lineSearch.C:217
A simple container for options an IOstream can normally have.
void setNewMeritValue(const scalar value)
Set new objective value.
Definition: lineSearch.C:150
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
virtual void setNewDeriv(const scalar deriv)
Set new directional derivative.
Definition: lineSearch.C:144
virtual void setDeriv(const scalar deriv)
Set directional derivative.
Definition: lineSearch.C:137
static autoPtr< lineSearch > New(const dictionary &dict, const Time &time, updateMethod &UpdateMethod)
Return a reference to the selected turbulence model.
Definition: lineSearch.C:96
word timeName
Definition: getTimeIndex.H:3
const dictionary & optionalSubDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary, otherwise return this dictionary.
Definition: dictionary.C:560
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
A class for handling words, derived from Foam::string.
Definition: word.H:63
void setOldMeritValue(const scalar value)
Set old objective value.
Definition: lineSearch.C:157
virtual lineSearch & operator++()
Increment iteration number and store merit value corresponding to the previous optimisation cycle...
Definition: lineSearch.C:223
virtual void reset()
Reset step to initial value.
Definition: lineSearch.C:164
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
defineTypeNameAndDebug(combustionModel, 0)
Abstract base class for line search methods.
Definition: lineSearch.H:53
Abstract base class for step update methods used in line search.
Definition: stepUpdate.H:50
messageStream Info
Information stream (stdout output on master, null elsewhere)
virtual bool loop()
Return true if lineSearch should continue and if so increment inner.
Definition: lineSearch.C:198
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition: error.H:635
Namespace for OpenFOAM.
dimensionSet clamp(const dimensionSet &a, const dimensionSet &range)
Definition: dimensionSet.C:271
IOerror FatalIOError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL IO ERROR&#39; header text and ...
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127
virtual bool computeGradient() const
Does line search need to update the gradient?
Definition: lineSearch.C:211