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