LBFGS.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
28Class
29 Foam::LBFGS
30
31Description
32 The quasi-Newton Limited-memory BFGS formula. Keeps nPrevSteps_ of the
33 y and s vectors and approximates the inverse areas through them.
34 Values of 3 < nPrevSteps_ < 20 are suggested.
35
36SourceFiles
37 LBFGS.C
38
39\*---------------------------------------------------------------------------*/
40
41#ifndef LBFGS_H
42#define LBFGS_H
43
44#include "updateMethod.H"
45#include "scalarMatrices.H"
46
47// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48
49namespace Foam
50{
51
52/*---------------------------------------------------------------------------*\
53 Class LBFGS Declaration
54\*---------------------------------------------------------------------------*/
56class LBFGS
57:
58 public updateMethod
59{
60private:
61
62 // Private Member Functions
63
64 //- No copy construct
65 LBFGS(const LBFGS&) = delete;
66
67 //- No copy assignment
68 void operator=(const LBFGS&) = delete;
69
70
71protected:
72
73 // Protected data
74
75 //- Step for the Newton method
76 scalar etaHessian_;
77
78 //- Number of first steepest descent steps
80
81 //- Map to active design variables
83
84 //- Number of old corrections and grad differences kept
85 label nPrevSteps_;
86
87 //- The previous differences of derivatives. Holds nPrevSteps_ fields
89
90 //- The previous correction. Holds nPrevSteps_ fields
92
93 //- The previous derivatives
95
96 //- The previous correction
98
99 //- Optimisation cycle counter
100 label counter_;
101
102 //- Allocate matrices in the first optimisation cycle
103 void allocateMatrices();
104
105 //- Move pointers in PtrList to the left and replace last element with
106 //- given field
107 void pivotFields(PtrList<scalarField>& list, const scalarField& f);
108
109 //- Update y and s vectors
110 void updateVectors();
111
112 //- Update design variables
113 void update();
114
115 //- Update based on steepest descent
117
118 //- Update based on LBFGS
119 void LBFGSUpdate();
120
121 //- Read old info from dict
122 void readFromDict();
123
124
125public:
126
127 //- Runtime type information
128 TypeName("LBFGS");
129
130
131 // Constructors
132
133 //- Construct from components
134 LBFGS(const fvMesh& mesh, const dictionary& dict);
135
136
137 //- Destructor
138 virtual ~LBFGS() = default;
139
140
141 // Member Functions
142
143 //- Compute design variables correction
144 void computeCorrection();
145
146 //- Update old correction. Useful for quasi-Newton methods coupled with
147 //- line search
148 virtual void updateOldCorrection(const scalarField& oldCorrection);
149
150 //- Write old info to dict
151 virtual void write();
152};
153
154
155// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
156
157} // End namespace Foam
158
159// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
160
161#endif
162
163// ************************************************************************* //
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
scalarField correctionOld_
The previous correction.
Definition: LBFGS.H:96
PtrList< scalarField > s_
The previous correction. Holds nPrevSteps_ fields.
Definition: LBFGS.H:90
label nSteepestDescent_
Number of first steepest descent steps.
Definition: LBFGS.H:78
labelList activeDesignVars_
Map to active design variables.
Definition: LBFGS.H:81
scalar etaHessian_
Step for the Newton method.
Definition: LBFGS.H:75
label nPrevSteps_
Number of old corrections and grad differences kept.
Definition: LBFGS.H:84
label counter_
Optimisation cycle counter.
Definition: LBFGS.H:99
void allocateMatrices()
Allocate matrices in the first optimisation cycle.
Definition: LBFGS.C:49
scalarField derivativesOld_
The previous derivatives.
Definition: LBFGS.H:93
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
virtual ~LBFGS()=default
Destructor.
TypeName("LBFGS")
Runtime type information.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
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.
labelList f(nPoints)
dictionary dict
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73