DBFGS.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 "DBFGS.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37  defineTypeNameAndDebug(DBFGS, 0);
39  (
40  updateMethod,
41  DBFGS,
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 Hessian to be a diagonal matrix
63 
64  // Allocate correct size and content to Hessian matrices
65  // has a max. capability of approximately 34000 variables.
66  HessianOld_ = temp;
67  Hessian_ = temp;
68 }
69 
70 
72 {
73  // Vectors needed to construct the inverse Hessian 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  HessianOld_[varI][varI] /= scaleFactor;
88  }
89  }
90 
91  scalar sBs = globalSum(leftMult(s, HessianOld_)*s);
92 
93  // Check curvature condition and apply dampening is necessary
94  scalar theta(1);
95  if (ys < gamma_*sBs)
96  {
98  << " y*s is below threshold. Using damped form" << endl;
99  theta = (scalar(1)-gamma_)*sBs/(sBs - ys);
100  }
101 
102  DebugInfo
103  << "Hessian curvature index " << ys << endl;
104 
105  scalarField r(theta*y + (scalar(1)-theta)*rightMult(HessianOld_, s));
106 
107  // Construct the inverse Hessian
108  Hessian_ =
109  HessianOld_
110  - outerProd(rightMult(HessianOld_, s), leftMult(s/sBs, HessianOld_))
111  + outerProd(r, r/globalSum(s*r));
112 }
113 
114 
116 {
117  SquareMatrix<scalar> HessianInv = inv(Hessian_);
118 
119  // In the first few iterations, use steepest descent but update the Hessian
120  // matrix
121  if (counter_ < nSteepestDescent_)
122  {
123  Info<< "Using steepest descent to update design variables" << endl;
124  correction_ = -eta_*objectiveDerivatives_;
125  }
126  // Else use DBFGS formula to update design variables
127  else
128  {
129  // Compute correction for active design variables
130  scalarField activeDerivs(activeDesignVars_.size(), Zero);
131  activeDerivs.map(objectiveDerivatives_, activeDesignVars_);
132  scalarField activeCorrection
133  (
134  -etaHessian_*rightMult(HessianInv, activeDerivs)
135  );
136 
137  // Transfer correction to the global list
138  correction_ = Zero;
139  forAll(activeDesignVars_, varI)
140  {
141  correction_[activeDesignVars_[varI]] = activeCorrection[varI];
142  }
143  }
144 
145  // Store fields for the next iteration
146  derivativesOld_ = objectiveDerivatives_;
147  correctionOld_ = correction_;
148  HessianOld_ = Hessian_;
149 }
150 
151 
153 {
154  if (optMethodIODict_.headerOk())
155  {
156  optMethodIODict_.readEntry("HessianOld", HessianOld_);
157  optMethodIODict_.readEntry("derivativesOld", derivativesOld_);
158  optMethodIODict_.readEntry("correctionOld", correctionOld_);
159  optMethodIODict_.readEntry("counter", counter_);
160  optMethodIODict_.readEntry("eta", eta_);
161 
162  label n = HessianOld_.n();
163  Hessian_ = SquareMatrix<scalar>(n, Zero);
164  correction_ = scalarField(correctionOld_.size(), Zero);
165  }
166 }
167 
168 
169 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
170 
171 Foam::DBFGS::DBFGS
172 (
173  const fvMesh& mesh,
174  const dictionary& dict
175 )
176 :
178 
179  // Construct null matrix since we dont know the dimension yet
180  etaHessian_
181  (
182  coeffsDict().lookupOrDefault<scalar>("etaHessian", 1)
183  ),
184  nSteepestDescent_
185  (
186  coeffsDict().lookupOrDefault<label>("nSteepestDescent", 1)
187  ),
188  activeDesignVars_(0),
189  scaleFirstHessian_
190  (
191  coeffsDict().lookupOrDefault<bool>("scaleFirstHessian", false)
192  ),
193  curvatureThreshold_
194  (
195  coeffsDict().lookupOrDefault<scalar>("curvatureThreshold", 1e-10)
196  ),
197  Hessian_(),
198  HessianOld_(),
199  derivativesOld_(0),
200  correctionOld_(0),
201  counter_(0),
202  gamma_(coeffsDict().lookupOrDefault<scalar>("gamma", 0.2))
203 
204 {
205  if
206  (
207  !coeffsDict().readIfPresent("activeDesignVariables", activeDesignVars_)
208  )
209  {
210  // If not, all available design variables will be used. Number is not
211  // know at the moment
212  Info<< "\t Did not find explicit definition of active design variables."
213  " Treating all available ones as active " << endl;
214  }
215 
216  // read old hessian, correction and derivatives, if present
217  readFromDict();
218 }
219 
220 
221 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
222 
224 {
225  if (counter_ == 0)
226  {
227  allocateMatrices();
228  }
229  else
230  {
231  updateHessian();
232  }
233 
234  update();
235  ++counter_;
236 }
237 
238 
240 {
241  updateMethod::updateOldCorrection(oldCorrection);
242  correctionOld_ = oldCorrection;
243 }
244 
245 
247 {
248  optMethodIODict_.add<SquareMatrix<scalar>>("HessianOld", HessianOld_, true);
249  optMethodIODict_.add<scalarField>("derivativesOld", derivativesOld_, true);
250  optMethodIODict_.add<scalarField>("correctionOld", correctionOld_, true);
251  optMethodIODict_.add<label>("counter", counter_, true);
252 
254 }
255 
256 
257 // ************************************************************************* //
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
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::DBFGS::allocateMatrices
void allocateMatrices()
Allocate matrices in the first optimisation cycle.
Definition: DBFGS.C:49
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
DBFGS.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
Foam::DBFGS::update
void update()
Update design variables.
Definition: DBFGS.C:115
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::inv
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
Definition: dimensionedSphericalTensor.C:73
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::DBFGS::updateOldCorrection
virtual void updateOldCorrection(const scalarField &oldCorrection)
Definition: DBFGS.C:239
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::DBFGS::HessianOld_
SquareMatrix< scalar > HessianOld_
The previous Hessian.
Definition: DBFGS.H:81
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::DBFGS::write
virtual void write()
Write old info to dict.
Definition: DBFGS.C:246
Foam::DBFGS::updateHessian
void updateHessian()
Update approximation of the inverse Hessian.
Definition: DBFGS.C:71
Foam::SquareMatrix< scalar >
Foam::updateMethod::write
virtual void write()
Write usefull quantities to files.
Definition: updateMethod.C:398
DebugInfo
#define DebugInfo
Report an information message using Foam::Info.
Definition: messageStream.H:350
Foam::DBFGS::Hessian_
SquareMatrix< scalar > Hessian_
The Hessian. Should have the size of the active design variables.
Definition: DBFGS.H:78
Foam::DBFGS::activeDesignVars_
labelList activeDesignVars_
Map to active design variables.
Definition: DBFGS.H:69
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::List::setSize
void setSize(const label newSize)
Alias for resize(const label)
Definition: ListI.H:146
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:294
Foam::DBFGS::computeCorrection
void computeCorrection()
Compute design variables correction.
Definition: DBFGS.C:223
Foam::I
static const Identity< scalar > I
Definition: Identity.H:95
y
scalar y
Definition: LISASMDCalcMethod1.H:14
Foam::DBFGS::readFromDict
void readFromDict()
Read old info from dict.
Definition: DBFGS.C:152