LBFGS.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-2019 PCOpt/NTUA
9 Copyright (C) 2013-2019 FOSS GP
10 Copyright (C) 2019-2020 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 "LBFGS.H"
32
33// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34
35namespace Foam
36{
39 (
41 LBFGS,
43 );
44}
45
46
47// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
48
50{
51 // Set active design variables, if necessary
53 {
55 }
56
57 // Allocate vectors
58 label nVars(activeDesignVars_.size());
59 for (label i = 0; i < nPrevSteps_; i++)
60 {
61 y_.set(i, new scalarField(nVars, Zero));
62 s_.set(i, new scalarField(nVars, Zero));
63 }
64}
65
66
68{
69 if (counter_ > nPrevSteps_)
70 {
71 // Reorder list by moving pointers down the line
72 labelList newOrder(nPrevSteps_, -1);
73 newOrder[0] = nPrevSteps_ - 1;
74 for (label i = 1; i < nPrevSteps_; ++i)
75 {
76 newOrder[i] = i - 1;
77 }
78 list.reorder(newOrder);
79
80 // Fill in last element with the provided field
81 list[nPrevSteps_ - 1] = f;
82 }
83 else
84 {
85 list[counter_ - 1] = f;
86 }
87}
88
89
91{
92 // Update list of y. Can only be done here since objectiveDerivatives_
93 // was not known at the end of the previous loop
94 scalarField yRecent
95 (objectiveDerivatives_ - derivativesOld_, activeDesignVars_);
96 pivotFields(y_, yRecent);
97 // Update list of s.
98 // correction_ holds the previous correction
99 scalarField sActive(correctionOld_, activeDesignVars_);
100 pivotFields(s_, sActive);
101
103 << "y fields" << nl << y_ << endl;
105 << "s fields" << nl << s_ << endl;
106}
107
108
110{
111 Info<< "Using steepest descent to update design variables" << endl;
112 correction_ = -eta_*objectiveDerivatives_;
113}
114
115
117{
118 // L-BFGS two loop recursion
119 //~~~~~~~~~~~~~~~~~~~~~~~~~~
120 label nSteps(min(counter_, nPrevSteps_));
121 label nLast(nSteps - 1);
122 scalarField q(objectiveDerivatives_, activeDesignVars_);
123 scalarField a(nSteps, Zero);
124 scalarField r(nSteps, Zero);
125 for (label i = nLast; i > -1; --i)
126 {
127 r[i] = 1./globalSum(y_[i]*s_[i]);
128 a[i] = r[i]*globalSum(s_[i]*q);
129 q -= a[i]*y_[i];
130 }
131
132 scalar gamma =
133 globalSum(y_[nLast]*s_[nLast])/globalSum(y_[nLast]*y_[nLast]);
134 q *= gamma;
135
136 scalarField b(activeDesignVars_.size(), Zero);
137 for (label i = 0; i < nSteps; ++i)
138 {
139 b = r[i]*globalSum(y_[i]*q);
140 q += s_[i]*(a[i] -b);
141 }
142
143 // Update correction
144 forAll(activeDesignVars_, varI)
145 {
146 correction_[activeDesignVars_[varI]] = -etaHessian_*q[varI];
147 }
148}
149
150
152{
153 // In the first few iterations, use steepest descent but update the Hessian
154 // matrix
155 if (counter_ < nSteepestDescent_)
156 {
157 steepestDescentUpdate();
158 }
159 // else use LBFGS formula to update the design variables
160 else
161 {
162 LBFGSUpdate();
163 }
164
165 // Store fields for the next iteration
166 derivativesOld_ = objectiveDerivatives_;
167 correctionOld_ = correction_;
168}
169
170
172{
173 if (optMethodIODict_.headerOk())
174 {
175 optMethodIODict_.readEntry("y", y_);
176 optMethodIODict_.readEntry("s", s_);
177 optMethodIODict_.readEntry("derivativesOld", derivativesOld_);
178 optMethodIODict_.readEntry("counter", counter_);
179 optMethodIODict_.readEntry("eta", eta_);
180 optMethodIODict_.readEntry("correctionOld", correctionOld_);
181
182 correction_ = scalarField(correctionOld_.size(), Zero);
183
184 if (activeDesignVars_.empty())
185 {
186 activeDesignVars_ = identity(derivativesOld_.size());
187 }
188 }
189}
190
191
192// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
193
195(
196 const fvMesh& mesh,
197 const dictionary& dict
198)
199:
201
202 // Construct null matrix since we dont know the dimension yet
203 etaHessian_
204 (
205 coeffsDict().getOrDefault<scalar>("etaHessian", 1)
206 ),
207 nSteepestDescent_
208 (
209 coeffsDict().getOrDefault<label>("nSteepestDescent", 1)
210 ),
211 activeDesignVars_(0),
212 nPrevSteps_
213 (
214 coeffsDict().getOrDefault<label>("nPrevSteps", 10)
215 ),
216 y_(nPrevSteps_),
217 s_(nPrevSteps_),
218 derivativesOld_(0),
219 counter_(Zero)
220{
221 if
222 (
223 !coeffsDict().readIfPresent("activeDesignVariables", activeDesignVars_)
224 )
225 {
226 // If not, all available design variables will be used. Number is not
227 // know at the moment
228 Info<< "\t Did not find explicit definition of active design variables. "
229 << "Treating all available ones as active " << endl;
230 }
231
232 // Read old Hessian, correction and derivatives, if present
233 readFromDict();
234}
235
236
237// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
238
240{
241 if (counter_ == 0)
242 {
243 allocateMatrices();
244 }
245 else
246 {
247 updateVectors();
248 }
249
250 update();
251 ++counter_;
252}
253
254
256{
258 correctionOld_ = oldCorrection;
259}
260
261
263{
264 optMethodIODict_.add<PtrList<scalarField>>("y", y_, true);
265 optMethodIODict_.add<PtrList<scalarField>>("s", s_, true);
266 optMethodIODict_.add<scalarField>("derivativesOld", derivativesOld_, true);
267 optMethodIODict_.add<scalarField>("correctionOld", correctionOld_, true);
268 optMethodIODict_.add<label>("counter", counter_, true);
269
271}
272
273
274// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
The quasi-Newton Limited-memory BFGS formula. Keeps nPrevSteps_ of the y and s vectors and approximat...
Definition: LBFGS.H:58
PtrList< scalarField > y_
The previous differences of derivatives. Holds nPrevSteps_ fields.
Definition: LBFGS.H:87
void updateVectors()
Update y and s vectors.
Definition: LBFGS.C:90
void pivotFields(PtrList< scalarField > &list, const scalarField &f)
Definition: LBFGS.C:67
void computeCorrection()
Compute design variables correction.
Definition: LBFGS.C:239
void LBFGSUpdate()
Update based on LBFGS.
Definition: LBFGS.C:116
void readFromDict()
Read old info from dict.
Definition: LBFGS.C:171
PtrList< scalarField > s_
The previous correction. Holds nPrevSteps_ fields.
Definition: LBFGS.H:90
labelList activeDesignVars_
Map to active design variables.
Definition: LBFGS.H:81
label nPrevSteps_
Number of old corrections and grad differences kept.
Definition: LBFGS.H:84
void allocateMatrices()
Allocate matrices in the first optimisation cycle.
Definition: LBFGS.C:49
virtual void write()
Write old info to dict.
Definition: LBFGS.C:262
virtual void updateOldCorrection(const scalarField &oldCorrection)
Definition: LBFGS.C:255
void update()
Update design variables.
Definition: LBFGS.C:151
void steepestDescentUpdate()
Update based on steepest descent.
Definition: LBFGS.C:109
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
bool empty() const noexcept
True if the UList is empty (ie, size() is zero)
Definition: UListI.H:427
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
void reorder(const labelUList &oldToNew, const bool check=false)
Definition: UPtrList.C:69
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
Abstract base class for optimisation methods.
Definition: updateMethod.H:55
scalarField objectiveDerivatives_
Derivatives of the objective functions.
Definition: updateMethod.H:68
dictionary coeffsDict()
Return optional dictionary with parameters specific to each method.
Definition: updateMethod.C:254
virtual void write()
Write useful quantities to files.
Definition: updateMethod.C:398
virtual void updateOldCorrection(const scalarField &oldCorrection)
Definition: updateMethod.C:390
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
const scalar gamma
Definition: EEqn.H:9
mesh update()
dynamicFvMesh & mesh
#define DebugInfo
Report an information message using Foam::Info.
Namespace for OpenFOAM.
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i)
Definition: labelList.C:38
messageStream Info
Information stream (stdout output on master, null elsewhere)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
labelList f(nPoints)
dictionary dict
volScalarField & b
Definition: createFields.H:27
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
propsDict readIfPresent("fields", acceptFields)