BFGS.H
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-------------------------------------------------------------------------------
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
29Class
30 Foam::BFGS
31
32Description
33 The quasi-Newton BFGS formula
34
35SourceFiles
36 BFGS.C
37
38\*---------------------------------------------------------------------------*/
39
40#ifndef BFGS_H
41#define BFGS_H
42
43#include "updateMethod.H"
44#include "scalarMatrices.H"
45
46// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47
48namespace Foam
49{
50
51/*---------------------------------------------------------------------------*\
52 Class BFGS Declaration
53\*---------------------------------------------------------------------------*/
55class BFGS
56:
57 public updateMethod
58{
59private:
60
61 // Private Member Functions
62
63 //- No copy construct
64 BFGS(const BFGS&) = delete;
65
66 //- No copy assignment
67 void operator=(const BFGS&) = delete;
68
69
70protected:
71
72 // Protected data
73
74 //- Step for the Newton method
75 scalar etaHessian_;
76
77 //- Number of first steepest descent steps
79
80 //- Map to active design variables
82
83 //- Scale the iniitial unitary Hessian approximation
85
86 //- Curvature threshold
88
89 //- The Hessian inverse. Should have the size of the active design
90 //- variables
92
93 //- The previous Hessian inverse
95
96 //- The previous derivatives
98
99 //- The previous correction
101
102 //- Optimisation cycle counter
103 label counter_;
104
105 //- Allocate matrices in the first optimisation cycle
106 void allocateMatrices();
107
108 //- Update approximation of the inverse Hessian
109 void updateHessian();
110
111 //- Update design variables
112 void update();
113
114 //- Read old info from dict
115 void readFromDict();
116
117
118public:
119
120 //- Runtime type information
121 TypeName("BFGS");
122
123
124 // Constructors
125
126 //- Construct from components
127 BFGS(const fvMesh& mesh, const dictionary& dict);
128
129
130 //- Destructor
131 virtual ~BFGS() = default;
132
133
134 // Member Functions
135
136 //- Compute design variables correction
137 void computeCorrection();
138
139 //- Update old correction. Useful for quasi-Newton methods coupled with
140 //- line search
141 virtual void updateOldCorrection(const scalarField& oldCorrection);
142
143 //- Write old info to dict
144 virtual void write();
145};
146
147
148// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
149
150} // End namespace Foam
151
152// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
153
154#endif
155
156// ************************************************************************* //
The quasi-Newton BFGS formula.
Definition: BFGS.H:57
void computeCorrection()
Compute design variables correction.
Definition: BFGS.C:217
scalar curvatureThreshold_
Curvature threshold.
Definition: BFGS.H:86
SquareMatrix< scalar > HessianInvOld_
The previous Hessian inverse.
Definition: BFGS.H:93
void readFromDict()
Read old info from dict.
Definition: BFGS.C:144
scalarField correctionOld_
The previous correction.
Definition: BFGS.H:99
TypeName("BFGS")
Runtime type information.
label nSteepestDescent_
Number of first steepest descent steps.
Definition: BFGS.H:77
labelList activeDesignVars_
Map to active design variables.
Definition: BFGS.H:80
scalar etaHessian_
Step for the Newton method.
Definition: BFGS.H:74
label counter_
Optimisation cycle counter.
Definition: BFGS.H:102
void allocateMatrices()
Allocate matrices in the first optimisation cycle.
Definition: BFGS.C:49
bool scaleFirstHessian_
Scale the iniitial unitary Hessian approximation.
Definition: BFGS.H:83
scalarField derivativesOld_
The previous derivatives.
Definition: BFGS.H:96
virtual void write()
Write old info to dict.
Definition: BFGS.C:240
SquareMatrix< scalar > HessianInv_
Definition: BFGS.H:90
virtual void updateOldCorrection(const scalarField &oldCorrection)
Definition: BFGS.C:233
void update()
Update design variables.
Definition: BFGS.C:109
virtual ~BFGS()=default
Destructor.
void updateHessian()
Update approximation of the inverse Hessian.
Definition: BFGS.C:67
A templated (N x N) square matrix of objects of <Type>, containing N*N elements, derived from Matrix.
Definition: SquareMatrix.H:66
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
dynamicFvMesh & mesh
Namespace for OpenFOAM.
dictionary dict
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73