SQP.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::SQP
31
32Description
33 The quasi-Newton SQP formula for constrained optimisation
34
35SourceFiles
36 SQP.C
37
38\*---------------------------------------------------------------------------*/
39
40#ifndef SQP_H
41#define SQP_H
42
44#include "OFstream.H"
45
46// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47
48namespace Foam
49{
50
51/*---------------------------------------------------------------------------*\
52 Class SQP Declaration
53\*---------------------------------------------------------------------------*/
55class SQP
56:
58{
59protected:
60
61 // Protected data
62
63 //- Step for the Newton method
64 scalar etaHessian_;
65
66 //- Map to active design variables
68
69 //- Scale the initial unitary Hessian approximation
71
72 //- Curvature threshold
73 scalar dumpingThreshold_;
74
75 //- Derivatives of the Lagrangian function
77
78 //- The Hessian inverse. Should have the size of the active design
79 //- variables
81
82 //- The previous Hessian inverse
84
85 //- The previous objective derivatives
87
88 //- The previous constraint derivatives
90
91 //- The previous correction
93
94 //- Lagrange multipliers
96
97 //- Optimisation cycle count
98 label counter_;
99
100 //- Name of the objective folder
102
103 //- File including the l1 merit function
105
106 //- Penalty value for the merit function
107 scalar mu_;
108
109 //- Safety factor
110 scalar delta_;
111
112
113private:
114
115 // Private Member Functions
116
117 //- No copy construct
118 SQP(const SQP&) = delete;
119
120 //- No copy assignment
121 void operator=(const SQP&) = delete;
122
123 //- Make folder holding the Lagrangian file
124 void makeFolder();
125
126 //- Allocate fields and matrices when size of design variables
127 //- is known
128 void allocateMatrices();
129
130 //- Update the Hessian matrix
131 void updateHessian();
132
133 //- Compute new lamdas and update correction
134 void computeLagrangeMultipliersAndCorrect();
135
136 //- Store old fields needed for the next iter
137 void storeOldFields();
138
139 //- Read old values from dict, if present
140 void readFromDict();
141
142
143public:
144
145 //- Runtime type information
146 TypeName("SQP");
147
148
149 // Constructors
150
151 //- Construct from components
152 SQP(const fvMesh& mesh, const dictionary& dict);
153
154
155 //- Destructor
156 virtual ~SQP() = default;
157
158
159 // Member Functions
160
161 //- Compute design variables correction
162 void computeCorrection();
163
164 //- Compute merit function. Could be different than the objective
165 //- in the presence of constraints
166 virtual scalar computeMeritFunction();
167
168 //- Derivative of the merit function. Could be different than the
169 //- objective derivative in the presence of constraints
170 virtual scalar meritFunctionDirectionalDerivative();
171
172 //- Update old correction. Useful for quasi-Newton methods coupled with
173 //- line search
174 virtual void updateOldCorrection(const scalarField& oldCorrection);
175
176 //- Write useful quantities to files
177 virtual void write();
178};
179
180
181// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
182
183} // End namespace Foam
184
185// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
186
187#endif
188
189// ************************************************************************* //
The quasi-Newton SQP formula for constrained optimisation.
Definition: SQP.H:57
void computeCorrection()
Compute design variables correction.
Definition: SQP.C:340
virtual scalar meritFunctionDirectionalDerivative()
Definition: SQP.C:382
virtual ~SQP()=default
Destructor.
scalarField lamdas_
Lagrange multipliers.
Definition: SQP.H:94
SquareMatrix< scalar > HessianOld_
The previous Hessian inverse.
Definition: SQP.H:82
scalarField correctionOld_
The previous correction.
Definition: SQP.H:91
scalar mu_
Penalty value for the merit function.
Definition: SQP.H:106
labelList activeDesignVars_
Map to active design variables.
Definition: SQP.H:66
scalar etaHessian_
Step for the Newton method.
Definition: SQP.H:63
List< scalarField > constraintDerivativesOld_
The previous constraint derivatives.
Definition: SQP.H:88
label counter_
Optimisation cycle count.
Definition: SQP.H:97
autoPtr< OFstream > meritFunctionFile_
File including the l1 merit function.
Definition: SQP.H:103
scalarField objectiveDerivativesOld_
The previous objective derivatives.
Definition: SQP.H:85
virtual scalar computeMeritFunction()
Definition: SQP.C:365
scalar dumpingThreshold_
Curvature threshold.
Definition: SQP.H:72
scalar delta_
Safety factor.
Definition: SQP.H:109
TypeName("SQP")
Runtime type information.
SquareMatrix< scalar > Hessian_
Definition: SQP.H:79
bool scaleFirstHessian_
Scale the initial unitary Hessian approximation.
Definition: SQP.H:69
virtual void write()
Write useful quantities to files.
Definition: SQP.C:399
virtual void updateOldCorrection(const scalarField &oldCorrection)
Definition: SQP.C:392
fileName objFunctionFolder_
Name of the objective folder.
Definition: SQP.H:100
scalarField LagrangianDerivatives_
Derivatives of the Lagrangian function.
Definition: SQP.H:75
A templated (N x N) square matrix of objects of <Type>, containing N*N elements, derived from Matrix.
Definition: SquareMatrix.H:66
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
Abstract base class for optimisation methods supporting constraints. Does not add functionality to up...
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
A class for handling file names.
Definition: fileName.H:76
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
dynamicFvMesh & mesh
Namespace for OpenFOAM.
dictionary dict
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73