BFGS.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 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 "BFGS.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37  defineTypeNameAndDebug(BFGS, 0);
39  (
40  updateMethod,
41  BFGS,
42  dictionary
43  );
44 }
45 
46 
47 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
48 
50 {
51  // Set active design variables, if necessary
52  if (activeDesignVars_.empty())
53  {
56  {
57  activeDesignVars_[dvI] = dvI;
58  }
59  }
60 
61  // Set previous HessianInv to be a diagonal matrix
63 
64  // Allocate correct size and content to HessianInv matrices
65  // has a max. capability of approximately 34000 variables.
66  HessianInvOld_ = temp;
67  HessianInv_ = temp;
68 }
69 
70 
72 {
73  // Vectors needed to construct the inverse HessianInv matrix
74  scalarField y(activeDesignVars_.size(), Zero);
75  scalarField s(activeDesignVars_.size(), Zero);
76  y.map(objectiveDerivatives_ - derivativesOld_, activeDesignVars_);
77  s.map(correctionOld_, activeDesignVars_);
78 
79  scalar ys = globalSum(s*y);
80 
81  if (counter_ == 1 && scaleFirstHessian_)
82  {
83  scalar scaleFactor = ys/globalSum(y*y);
84  Info<< "Scaling Hessian with factor " << scaleFactor << endl;
85  forAll(activeDesignVars_, varI)
86  {
87  HessianInvOld_[varI][varI] *= scaleFactor;
88  }
89  }
90 
91  // Check curvature condition
92  if (ys < curvatureThreshold_)
93  {
95  << " y*s is below threshold! y*s=" << ys << endl;
96  }
97 
98  DebugInfo
99  << "Hessian curvature index " << ys << endl;
100 
101  // Construct the inverse HessianInv
102  HessianInv_ =
103  HessianInvOld_
104  + (ys + globalSum(leftMult(y, HessianInvOld_)*y))/sqr(ys)*outerProd(s, s)
105  - (scalar(1)/ys)*
106  (
107  outerProd(rightMult(HessianInvOld_, y), s)
108  + outerProd(s, leftMult(y, HessianInvOld_))
109  );
110 }
111 
112 
114 {
115  // In the first few iterations, use steepest descent but update the Hessian
116  // matrix
117  if (counter_ < nSteepestDescent_)
118  {
119  Info<< "Using steepest descent to update design variables" << endl;
120  correction_ = -eta_*objectiveDerivatives_;
121  }
122  // else use BFGS formula to update the design variables
123  else
124  {
125  // Compute correction for active design variables
126  scalarField activeDerivs(activeDesignVars_.size(), Zero);
127  activeDerivs.map(objectiveDerivatives_, activeDesignVars_);
128  scalarField activeCorrection
129  (
130  -etaHessian_*rightMult(HessianInv_, activeDerivs)
131  );
132 
133  // Transfer correction to the global list
134  correction_ = Zero;
135  forAll(activeDesignVars_, varI)
136  {
137  correction_[activeDesignVars_[varI]] = activeCorrection[varI];
138  }
139  }
140 
141  // Store fields for the next iteration
142  derivativesOld_ = objectiveDerivatives_;
143  correctionOld_ = correction_;
144  HessianInvOld_ = HessianInv_;
145 }
146 
147 
149 {
150  if (optMethodIODict_.headerOk())
151  {
152  optMethodIODict_.readEntry("HessianInvOld", HessianInvOld_);
153  optMethodIODict_.readEntry("derivativesOld", derivativesOld_);
154  optMethodIODict_.readEntry("correctionOld", correctionOld_);
155  optMethodIODict_.readEntry("counter", counter_);
156  optMethodIODict_.readEntry("eta", eta_);
157 
158  label n = HessianInvOld_.n();
159  HessianInv_ = SquareMatrix<scalar>(n, Zero);
160  correction_ = scalarField(correctionOld_.size(), Zero);
161  }
162 }
163 
164 
165 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
166 
167 Foam::BFGS::BFGS
168 (
169  const fvMesh& mesh,
170  const dictionary& dict
171 )
172 :
174  etaHessian_
175  (
176  coeffsDict().lookupOrDefault<scalar>("etaHessian", 1)
177  ),
178  nSteepestDescent_
179  (
180  coeffsDict().lookupOrDefault<label>("nSteepestDescent", 1)
181  ),
182  activeDesignVars_(0),
183  scaleFirstHessian_
184  (
185  coeffsDict().lookupOrDefault<bool>("scaleFirstHessian", false)
186  ),
187  curvatureThreshold_
188  (
189  coeffsDict().lookupOrDefault<scalar>("curvatureThreshold", 1e-10)
190  ),
191  // Construct null matrix since we dont know the dimension yet
192  HessianInv_(),
193  HessianInvOld_(),
194  derivativesOld_(0),
195  correctionOld_(0),
196  counter_(Zero)
197 {
198  if
199  (
200  !coeffsDict().readIfPresent("activeDesignVariables", activeDesignVars_)
201  )
202  {
203  // If not, all available design variables will be used. Number is not
204  // know at the moment
205  Info<< "\t Did not find explicit definition of active design variables."
206  << " Treating all available ones as active" << endl;
207  }
208 
209  // Read old hessian, correction and derivatives, if present
210  readFromDict();
211 }
212 
213 
214 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
215 
217 {
218  if (counter_ == 0)
219  {
220  allocateMatrices();
221  }
222  else
223  {
224  updateHessian();
225  }
226 
227  update();
228  ++counter_;
229 }
230 
231 
233 {
234  updateMethod::updateOldCorrection(oldCorrection);
235  correctionOld_ = oldCorrection;
236 }
237 
238 
240 {
241  optMethodIODict_.add<SquareMatrix<scalar>>
242  (
243  "HessianInvOld",
244  HessianInvOld_,
245  true
246  );
247  optMethodIODict_.add<scalarField>("derivativesOld", derivativesOld_, true);
248  optMethodIODict_.add<scalarField>("correctionOld", correctionOld_, true);
249  optMethodIODict_.add<label>("counter", counter_, true);
250 
252 }
253 
254 
255 // ************************************************************************* //
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::BFGS::readFromDict
void readFromDict()
Read old info from dict.
Definition: BFGS.C:148
s
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Definition: gmvOutputSpray.H:25
update
mesh update()
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
Foam::BFGS::write
virtual void write()
Write old info to dict.
Definition: BFGS.C:239
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::Field::map
void map(const UList< Type > &mapF, const labelUList &mapAddressing)
1 to 1 map from the given field
Definition: Field.C:263
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:62
Foam::updateMethod::updateOldCorrection
virtual void updateOldCorrection(const scalarField &oldCorrection)
Definition: updateMethod.C:390
Foam::Field< scalar >
Foam::updateMethod
Abstract base class for optimisation methods.
Definition: updateMethod.H:54
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::BFGS::updateOldCorrection
virtual void updateOldCorrection(const scalarField &oldCorrection)
Definition: BFGS.C:232
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
Foam::updateMethod::correction_
scalarField correction_
Design variables correction.
Definition: updateMethod.H:80
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:84
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::SquareMatrix< scalar >
Foam::updateMethod::write
virtual void write()
Write usefull quantities to files.
Definition: updateMethod.C:398
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
DebugInfo
#define DebugInfo
Report an information message using Foam::Info.
Definition: messageStream.H:350
Foam::BFGS::allocateMatrices
void allocateMatrices()
Allocate matrices in the first optimisation cycle.
Definition: BFGS.C:49
Foam::BFGS::update
void update()
Update design variables.
Definition: BFGS.C:113
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::BFGS::activeDesignVars_
labelList activeDesignVars_
Map to active design variables.
Definition: BFGS.H:80
Foam::List::setSize
void setSize(const label newSize)
Alias for resize(const label)
Definition: ListI.H:146
BFGS.H
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::BFGS::computeCorrection
void computeCorrection()
Compute design variables correction.
Definition: BFGS.C:216
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:294
Foam::I
static const Identity< scalar > I
Definition: Identity.H:95
Foam::BFGS::HessianInv_
SquareMatrix< scalar > HessianInv_
Definition: BFGS.H:90
y
scalar y
Definition: LISASMDCalcMethod1.H:14
Foam::BFGS::HessianInvOld_
SquareMatrix< scalar > HessianInvOld_
The previous Hessian inverse.
Definition: BFGS.H:93
Foam::BFGS::updateHessian
void updateHessian()
Update approximation of the inverse Hessian.
Definition: BFGS.C:71