Loading...
Searching...
No Matches
SQP.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-2021 PCOpt/NTUA
9 Copyright (C) 2013-2021 FOSS GP
10 Copyright (C) 2019-2020 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\*---------------------------------------------------------------------------*/
29
30#include "SQP.H"
31#include "IOmanip.H"
33
34// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35
36namespace Foam
37{
40 (
42 SQP,
44 );
46 (
48 SQP,
50 );
51}
52
53
54// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
55
57{
58 // Vectors needed to construct the (inverse) Hessian matrix
61 scalarField LagrangianDerivativesOld = derivativesOld_;
63 {
65 LagrangianDerivativesOld -= lamdas_[cI] * constraintDerivativesOld_[cI];
66 }
67 y.map(LagrangianDerivatives_ - LagrangianDerivativesOld, activeDesignVars_);
69
70 scalar ys = globalSum(s*y);
71 if (counter_ == 1 && scaleFirstHessian_)
72 {
73 if (ys > scalar(0))
74 {
75 scalar scaleFactor = ys/globalSum(y*y);
76 Info<< "Scaling Hessian with factor " << scaleFactor << endl;
78 {
79 Hessian_()[varI][varI] /= scaleFactor;
80 }
81 }
82 else
83 {
85 << " y*s is negative. Skipping the scaling of the first Hessian"
86 << endl;
87 }
88 }
89 scalar sBs = globalSum(leftMult(s, Hessian_())*s);
90
91 // Check curvature condition
92 scalar theta(1);
93 if (ys < dumpingThreshold_*sBs)
94 {
96 << " y*s is below threshold. Using damped form" << endl;
97 theta = (1 - dumpingThreshold_)*sBs/(sBs - ys);
98 }
99 scalarField r(theta*y + (scalar(1) - theta)*rightMult(Hessian_(), s));
101 << "Unmodified Hessian curvature index " << ys << endl;
103 << "Modified Hessian curvature index " << globalSum(r*s) << endl;
104
105 // Update the Hessian
106 Hessian_() +=
108 + outerProd(r, r/globalSum(s*r));
109}
110
111
113{
114 // Also denoted below as W
115 SquareMatrix<scalar> HessianInv = inv(Hessian_());
116 if (debug > 1)
117 {
118 Info<< "Hessian " << Hessian_() << endl;
119 Info<< "HessianInv " << HessianInv << endl;
120 label n = Hessian_().n();
122 for (label k = 0; k < n; k++)
123 {
124 for (label l = 0; l < n; l++)
125 {
126 scalar elem(Zero);
127 for (label i = 0; i < n; i++)
128 {
129 elem += Hessian_()[k][i] * HessianInv[i][l];
130 }
131 test[k][l]=elem;
132 }
133 }
134 Info<< "Validation " << test << endl;
135 }
136
137 // Compute new Lagrange multipliers
138 label nc = constraintDerivatives_.size();
139 scalarField activeDerivs(activeDesignVars_.size(), Zero);
140
141 // activeDerivs.map(objectiveDerivatives_, activeDesignVars_);
142 activeDerivs.map(LagrangianDerivatives_, activeDesignVars_);
143 scalarField WgradL = rightMult(HessianInv, activeDerivs);
144
145 scalarField lamdaRHS(nc, Zero);
146 forAll(lamdaRHS, cI)
147 {
148 scalarField activeConsDerivs(activeDesignVars_.size(), Zero);
149 activeConsDerivs.map(constraintDerivatives_[cI], activeDesignVars_);
150 lamdaRHS[cI] = globalSum(activeConsDerivs * WgradL) - cValues_[cI];
151 if (debug > 1)
152 {
153 Info<< "lamdaRHS total|deriv part|constraint part "
154 << lamdaRHS[cI] << " " << globalSum(activeConsDerivs * WgradL)
155 << " " << cValues_[cI] << endl;
156 }
157 }
158
159 // lhs for the lamda system
160 SquareMatrix<scalar> AWA(nc, Zero);
162 for (label j = 0; j < nc; j++)
163 {
164 scalarField gradcJ(activeDesignVars_.size(), Zero);
165 gradcJ.map(constraintDerivatives_[j], activeDesignVars_);
166 WA.set(j, new scalarField(rightMult(HessianInv, gradcJ)));
167 for (label i = 0; i < nc; i++)
168 {
169 scalarField gradcI(activeDesignVars_.size(), Zero);
170 gradcI.map(constraintDerivatives_[i], activeDesignVars_);
171 AWA[i][j] = globalSum(gradcI * WA[j]);
172 }
173 }
174 SquareMatrix<scalar> invAWA = inv(AWA);
175 scalarField deltaLamda = rightMult(invAWA, lamdaRHS);
176 if (debug > 1)
177 {
178 Info<< "AWA " << AWA << endl;
179 Info<< "AWAInv " << invAWA << endl;
180 Info<< "lamda update " << deltaLamda << endl;
181 }
182 lamdas_ += deltaLamda;
183
184 // Compute design variables correction
185 scalarField activeCorrection(-WgradL);
186 forAll(WA, cI)
187 {
188 activeCorrection += WA[cI]*deltaLamda[cI];
189 }
190 activeCorrection *= etaHessian_;
191 // Transfer correction to the global list
192 correction_ = Zero;
193 forAll(activeDesignVars_, varI)
194 {
195 correction_[activeDesignVars_[varI]] = activeCorrection[varI];
196 }
197 if (counter_ == 0)
198 {
199 correction_ *= eta_;
200 }
201}
202
203
204void Foam::SQP::storeOldFields()
205{
206 derivativesOld_ = objectiveDerivatives_;
207 if (constraintDerivativesOld_.empty())
208 {
209 constraintDerivativesOld_.setSize(constraintDerivatives_.size());
210 }
211 forAll(constraintDerivativesOld_, cI)
212 {
213 constraintDerivativesOld_[cI] = constraintDerivatives_[cI];
214 }
215 correctionOld_ = correction_;
216}
217
218
220{
221 // Assumes that all constraints are known by all processors
222 // What about constraints directly imposed on distributed design variables?
223 // These should be met in each iteration of the algorithm, so,
224 // most probably, there is no problem
225 return sum(mag(cValues_));
226}
227
228
229// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
230
231Foam::SQP::SQP
232(
233 const fvMesh& mesh,
234 const dictionary& dict,
235 autoPtr<designVariables>& designVars,
236 const label nConstraints,
237 const word& type
238)
239:
240 quasiNewton(mesh, dict, designVars, nConstraints, type),
241 SQPBase(mesh, dict, designVars, *this, type),
242 dumpingThreshold_
243 (
244 coeffsDict(type).getOrDefault<scalar>("dumpingThreshold", 0.2)
245 )
248}
249
250
251// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
252
254{
255 LagrangianDerivatives_ = objectiveDerivatives_;
257
258 // Store fields for the next iteration and write them to file
259 storeOldFields();
260}
261
262
264{
265 // If condition is not met, update mu value
266 if (mu_ < max(mag(lamdas_)) + delta_)
267 {
268 mu_ = max(mag(lamdas_)) + 2*delta_;
269 if (debug > 1)
270 {
271 Info<< "Updated mu value to " << mu_ << endl;
272 }
274 scalar L = objectiveValue_ + mu_*sum(mag(cValues_));
275
276 return L;
277}
278
279
281{
282 scalar deriv =
287}
288
293}
294
295
297{
298 return SQPBase::writeMeritFunction(*this);
299}
300
301
302// ************************************************************************* //
Istream and Ostream manipulators taking arguments.
scalar y
label k
label n
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
Base class for Sequantial Quadratic Programming (SQP) methods.
Definition SQPBase.H:53
scalarField lamdas_
Lagrange multipliers.
Definition SQPBase.H:71
virtual bool addToFile(Ostream &os) const
Write continuation info.
Definition SQPBase.C:102
scalar mu_
Penalty value for the merit function.
Definition SQPBase.H:86
virtual bool writeMeritFunction(const updateMethod &UpdateMethod)
Write info about the merit function.
Definition SQPBase.C:115
List< scalarField > constraintDerivativesOld_
The previous constraint derivatives.
Definition SQPBase.H:66
scalar delta_
Safety factor.
Definition SQPBase.H:91
virtual scalar meritFunctionConstraintPart() const =0
Get the part the merit function that depends on the constraints.
scalarField LagrangianDerivatives_
Derivatives of the Lagrangian function.
Definition SQPBase.H:61
The quasi-Newton SQP formula for constrained optimisation.
Definition SQP.H:54
void computeCorrection()
Compute design variables correction.
Definition SQP.C:246
virtual scalar meritFunctionDirectionalDerivative()
Derivative of the merit function. Could be different than the objective derivative in the presence of...
Definition SQP.C:273
virtual bool writeAuxiliaryData()
Write merit function information.
Definition SQP.C:289
virtual bool writeData(Ostream &os) const
Write useful quantities to files.
Definition SQP.C:283
virtual scalar computeMeritFunction()
Compute merit function. Could be different than the objective in the presence of constraints.
Definition SQP.C:256
scalar dumpingThreshold_
Curvature threshold.
Definition SQP.H:62
A templated (N x N) square matrix of objects of <Type>, containing N*N elements, derived from Matrix.
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
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:133
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Base class for quasi-Newton methods.
Definition quasiNewton.H:52
void computeCorrection()
Compute design variables correction.
Definition quasiNewton.C:83
scalarField correctionOld_
The previous correction.
Definition quasiNewton.H:88
autoPtr< SquareMatrix< scalar > > Hessian_
The Hessian or its inverse, depending on the deriving class.
Definition quasiNewton.H:78
virtual bool writeData(Ostream &os) const
Write useful quantities to files.
virtual void updateHessian()=0
Update approximation of the inverse Hessian.
virtual void update()=0
Update design variables.
bool scaleFirstHessian_
Scale the initial unitary Hessian approximation.
Definition quasiNewton.H:70
scalarField derivativesOld_
The previous derivatives.
Definition quasiNewton.H:83
void allocateHessian()
Allocate the Hessian matrix.
Definition quasiNewton.C:35
Abstract base class for optimisation methods.
scalarField objectiveDerivatives_
Derivatives of the objective functions.
scalarField correction_
Design variables correction.
dictionary coeffsDict(const word &type) const
Return optional dictionary with parameters specific to each method.
scalar globalSum(const scalarField &field)
Compute either global or local sum, based on globalSum flag.
scalarField cValues_
Constraint values.
scalar objectiveValue_
Objective value.
PtrList< scalarField > constraintDerivatives_
Derivatives of the constraints.
const scalarField leftMult(const scalarField &, const SquareMatrix< scalar > &)
const labelList & activeDesignVars_
Map to active design variables.
const scalarField rightMult(const SquareMatrix< scalar > &, const scalarField &)
label counter_
Optimisation cycle count.
label nConstraints() const
Get the number of constraints.
SquareMatrix< scalar > outerProd(const scalarField &, const scalarField &)
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
dynamicFvMesh & mesh
OBJstream os(runTime.globalPath()/outputName)
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))
#define DebugInfo
Report an information message using Foam::Info.
#define WarningInFunction
Report a warning using Foam::Warning.
Namespace for handling debugging switches.
Definition debug.C:45
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
messageStream Info
Information stream (stdout output on master, null elsewhere).
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1, const label comm)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
const vector L(dict.get< vector >("L"))