Loading...
Searching...
No Matches
ISQP.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) 2020-2021 PCOpt/NTUA
9 Copyright (C) 2020-2021 FOSS GP
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27
28Class
29 Foam::ISQP
30
31Description
32 An L-BFGS-based SQP algorithm for computing the update of the design
33 variables in the presence of inequality constraints. The QP problem is
34 solved using the interior point method (hence the I in the classe name,
35 which is not standard in the terminology used in the literature). The
36 (potentially dense) linear system formulated by the interior point method
37 is solved using Conjugate Gradient and a choice from three preconditioners
38 (diagonal, inverse Hessian, Sherman-Morrison), using matrix-vector products
39 to avoid storing the LHS matrix.
40
41 Bound constraints on the design variables will also be included, if set by
42 the designVariables known by the updateMethod. If the QP problem is
43 infeasible, the algorithm can still be used by setting includeExtraVars_
44 to true, to allow a computation of an update, despite not being able to
45 satisfy all the constraints of the QP problem. Alternatively,
46 targetConstraintReduction can be set to a number lower than 1 to help with
47 the feasibility of the dual problem.
48
49
50SourceFiles
51 ISQP.C
52
53\*---------------------------------------------------------------------------*/
54
55#ifndef ISQP_H
56#define ISQP_H
57
58#include "LBFGS.H"
59#include "SQPBase.H"
60
61// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
62
63namespace Foam
64{
66/*---------------------------------------------------------------------------*\
67 Class ISQP Declaration
68\*---------------------------------------------------------------------------*/
69
70class ISQP
71:
72 public LBFGS,
73 public SQPBase
74{
75public:
76
77 //- Enumeration of preconditioners for the subproblem
79 {
81 };
82
83 //- Names of preconditioners for the subproblem
85
86
87protected:
88
89 // Protected data
90
91 //- Should Lagrange multipliers be allocated
94 //- Are bound constraints included?
96
97 //- Are additional design variables included?
98 // These are introduced to find relaxed solutions, even if the
99 // original problem does not have any feasible points
101
102 //- The set of design variables being updated during the subproblem.
103 // Size is that of the active design variables.
104 // Correction will end up being the difference of this and the old
105 // design variables.
107
108 // Lagrange multipliers and slack variables
109
110 //- Slack variables for the inequality constraints
112
113 //- Lagrange multipliers for the lower bound constraints
116 //- Slack variables the lower bound constraints
118
119 //- Lagrange multipliers for the upper bound constraints
121
122 //- Slack variables the upper bound constraints
124
125 //- Extra variables for finding solutions even in infeasible
126 //- problems
128
129 //- Lagrange multipliers for positive extra variables
131
132 //- Multiplier of the additional variables y in the Lagrangian, to
133 //- make them 'expensive'
135
136
137 // Fields holding updates of the design, Lagrange and slack variables
138
149
150 //- Infinitesimal quantity
151 // Updated during the inner iterations of the subproblem
152 scalar eps_;
154 //- Final eps quantity to be reached during the solution of the
155 //- subproblem
156 scalar epsMin_;
157
158 //- Maxmimum number of Newton iterations for the subproblem
160
161 //- Maxmimum number of line search iterations for each iteration of the
162 //- subproblem
165 //- Maxmimum number of iterations for the deltaX equation
168 //- Relative tolerance of the CG solver, solving for deltaP_
169 scalar relTol_;
171 //- Which preconditioner to use for the solution of the subproblem
173
174 //- Percentage reduction for the constraints.
175 // Can be used to relax QP problems with no feasible points
176 scalar cRed_;
177
178 //- Disable damping
179 // Generally discouraged but could be useful in cases with no
180 // constraints or only bound ones
181 bool disableDamping_;
182
183 //- File including the l1 merit function
185
187 // Protected Member Functions
188
189 //- Update sizes of fields related to the active design variables
190 void updateSizes();
192 //- Allocate multipliers for the bound constraints
194
195 //- Allocate Lagrange multipliers for the inequality constraints
198 //- Update the vectors accosiated with the Hessian matrix
199 void updateYS();
200
201 //- Allocate fields related to constraints
203
204
205 // Functions related to the solution of the primal-dual subproblem
206
207 //- Solve subproblem using a Newton optimiser
208 void solveSubproblem();
209
210 //- Compute direction for the Newton problem
213 //- Zero the updates computed in the previous optimisation cycle
214 void zeroUpdates();
215
216 //- Solve the equation for deltaX, which is the expensive part of
217 //- the Newtopn step.
218 // All other corrections can be computed based on this
220
221 //- Compute the RHS for the deltaX equation
222 // Currently not in use
224
225 //- CG algorithm for the solution of the deltaP eqn
226 void CGforDeltaP(const scalarField& FDx);
228 //- Procudt of the LHS of the deltaP eqn with a vector
230 (
231 const scalarField& vector
232 );
233
234 //- Preconditioner-vector product for CG
235 // In case a diagonal preconditioner is used, it is stored in
236 // precon for all CG iterations
238 (
239 const scalarField& vector,
241 );
242
243 //- Diagonal preconditioner of the deltaP eqn
245
246 //- Use the Sherman-Morrison formula to compute the matrix-free
247 //- product of the approximate inverse of the LHS with a vector
249
250 //- Recursive (naive) implementation of the rank-1 update
251 // WIP, efficient for up to 2 flow-related constraints
253 (
254 const PtrList<scalarField>& r1Updates,
255 const scalarField& p,
257 const scalarField& mult,
258 label n
259 );
260
261 //- Perform line search and return max residual corresponding to
262 //- the updated solution
263 scalar lineSearch();
264
265 //- Adjust step to satisfy cireteria
266 void adjustStep
267 (
268 scalar& step,
269 const scalar value,
270 const scalar update
271 );
272
273 //- Update the current solution using the known direction and the
274 //- given step length
275 void updateSolution(const scalar step);
276
277
278 // Residuals of the various KKT conditions
279
280 //- Compute and return residuals based on the current solution
282
283 //- Residual of the gradient of the Lagrangian
284 // Size is that of the active design variables
286
287 //- Product of the inverse Hessian with the residual of the
288 //- gradient of the Lagrangian.
289 // Avoid the formation of the Hessian matrix.
290 // Size is that of the active design variables.
292
293 //- Residual of the inequality constraint slackness
295
296 //- Residual of the complementarity slackness for the
297 //- inequality constraints
299
300 //- Residual of the lower bounds slackness
302
303 //- Residual of the upper bounds slackness
305
306 //- Residual of the complementarity slackness for the
307 //- lower bound constraints
309
310 //- Residual of the complementarity slackness for the
311 //- upper bound constraints
313
314 //- Residual of the Lagrangian derivative wrt the extra variables
316
317 //- Residual of the complementarity slackness for the
318 //- extra variables
320
321
322 //- Store old fields needed for the next iter
323 void storeOldFields();
324
325 //- Get the part the merit function that depends on the constraints
326 virtual scalar meritFunctionConstraintPart() const;
327
328
329private:
330
331 // Private Member Functions
332
333 //- No copy construct
334 ISQP(const ISQP&) = delete;
335
336 //- No copy assignment
337 void operator=(const ISQP&) = delete;
338
339
340public:
341
342 //- Runtime type information
343 TypeName("ISQP");
344
345
346 // Constructors
347
348 //- Construct from components
349 ISQP
350 (
351 const fvMesh& mesh,
352 const dictionary& dict,
353 autoPtr<designVariables>& designVars,
354 const label nConstraints,
355 const word& type
356 );
357
358
359 //- Destructor
360 virtual ~ISQP() = default;
361
362
363 // Member Functions
364
365 //- Compute design variables correction
366 void computeCorrection();
367
368 //- Compute merit function. Could be different than the objective
369 //- in the presence of constraints
370 virtual scalar computeMeritFunction();
371
372 //- Derivative of the merit function. Could be different than the
373 //- objective derivative in the presence of constraints
374 virtual scalar meritFunctionDirectionalDerivative();
375
376 //- Update old correction. Useful for quasi-Newton methods coupled with
377 //- line search
378 virtual void updateOldCorrection(const scalarField& oldCorrection);
379
380 //- Write useful quantities to files
381 virtual bool writeData(Ostream& os) const;
382
383 //- Write merit function information
384 virtual bool writeAuxiliaryData();
385};
386
387
388// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
389
390} // End namespace Foam
391
392// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
393
394#endif
395
396// ************************************************************************* //
label n
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition Enum.H:57
bool doAllocateLagrangeMultipliers_
Should Lagrange multipliers be allocated.
Definition ISQP.H:93
scalar lineSearch()
Perform line search and return max residual corresponding to the updated solution.
Definition ISQP.C:542
tmp< scalarField > matrixVectorProduct(const scalarField &vector)
Procudt of the LHS of the deltaP eqn with a vector.
Definition ISQP.C:320
void updateSolution(const scalar step)
Update the current solution using the known direction and the given step length.
Definition ISQP.C:647
scalar cRed_
Percentage reduction for the constraints.
Definition ISQP.H:219
tmp< scalarField > resFls()
Residual of the lower bounds slackness.
Definition ISQP.C:816
void CGforDeltaP(const scalarField &FDx)
CG algorithm for the solution of the deltaP eqn.
Definition ISQP.C:287
void solveSubproblem()
Solve subproblem using a Newton optimiser.
Definition ISQP.C:884
void computeCorrection()
Compute design variables correction.
Definition ISQP.C:1107
virtual scalar meritFunctionDirectionalDerivative()
Derivative of the merit function. Could be different than the objective derivative in the presence of...
Definition ISQP.C:1142
bool includeExtraVars_
Are additional design variables included?
Definition ISQP.H:106
tmp< scalarField > resFL()
Residual of the gradient of the Lagrangian.
Definition ISQP.C:727
void initialize()
Allocate fields related to constraints.
Definition ISQP.C:162
virtual bool writeAuxiliaryData()
Write merit function information.
Definition ISQP.C:1171
autoPtr< Function1< scalar > > c_
Multiplier of the additional variables y in the Lagrangian, to make them 'expensive'.
Definition ISQP.H:159
tmp< scalarField > invHFL()
Product of the inverse Hessian with the residual of the gradient of the Lagrangian.
Definition ISQP.C:759
tmp< scalarField > resFz()
Residual of the complementarity slackness for the extra variables.
Definition ISQP.C:874
void adjustStep(scalar &step, const scalar value, const scalar update)
Adjust step to satisfy cireteria.
Definition ISQP.C:634
autoPtr< scalarField > z_
Lagrange multipliers for positive extra variables.
Definition ISQP.H:153
virtual bool writeData(Ostream &os) const
Write useful quantities to files.
Definition ISQP.C:1159
tmp< scalarField > resFExtraVars()
Residual of the Lagrangian derivative wrt the extra variables.
Definition ISQP.C:864
scalarField deltaLamda_
Definition ISQP.H:165
virtual scalar meritFunctionConstraintPart() const
Get the part the merit function that depends on the constraints.
Definition ISQP.C:991
tmp< scalarField > preconVectorProduct(const scalarField &vector, autoPtr< scalarField > &precon)
Preconditioner-vector product for CG.
Definition ISQP.C:363
static const Enum< preconditioners > preconditionerNames
Names of preconditioners for the subproblem.
Definition ISQP.H:83
scalarField p_
The set of design variables being updated during the subproblem.
Definition ISQP.H:115
void computeNewtonDirection()
Compute direction for the Newton problem.
Definition ISQP.C:490
autoPtr< scalarField > us_
Slack variables the upper bound constraints.
Definition ISQP.H:142
scalarField deltaGs_
Definition ISQP.H:166
scalar epsMin_
Final eps quantity to be reached during the solution of the subproblem.
Definition ISQP.H:186
autoPtr< scalarField > deltaUs_
Definition ISQP.H:170
tmp< scalarField > resFuTilda()
Residual of the complementarity slackness for the upper bound constraints.
Definition ISQP.C:854
autoPtr< scalarField > deltaUTilda_
Definition ISQP.H:169
autoPtr< OFstream > meritFunctionFile_
File including the l1 merit function.
Definition ISQP.H:232
autoPtr< scalarField > deltaLs_
Definition ISQP.H:168
void zeroUpdates()
Zero the updates computed in the previous optimisation cycle.
Definition ISQP.C:193
virtual scalar computeMeritFunction()
Compute merit function. Could be different than the objective in the presence of constraints.
Definition ISQP.C:1133
tmp< scalarField > ShermanMorrisonRank1Update(const PtrList< scalarField > &r1Updates, const scalarField &p, const refPtr< scalarField > &diag, const scalarField &mult, label n)
Recursive (naive) implementation of the rank-1 update.
Definition ISQP.C:459
tmp< scalarField > resFGs()
Residual of the inequality constraint slackness.
Definition ISQP.C:784
scalar eps_
Infinitesimal quantity.
Definition ISQP.H:180
scalarField deltaP_
Definition ISQP.H:164
preconditioners
Enumeration of preconditioners for the subproblem.
Definition ISQP.H:76
@ diag
Definition ISQP.H:77
@ invHessian
Definition ISQP.H:77
@ ShermanMorrison
Definition ISQP.H:77
autoPtr< scalarField > deltaExtraVars_
Definition ISQP.H:171
label maxDxIters_
Maxmimum number of iterations for the deltaX equation.
Definition ISQP.H:202
void allocateBoundMultipliers()
Allocate multipliers for the bound constraints.
Definition ISQP.C:91
tmp< scalarField > resFus()
Residual of the upper bounds slackness.
Definition ISQP.C:830
label maxLineSearchIters_
Maxmimum number of line search iterations for each iteration of the subproblem.
Definition ISQP.H:197
autoPtr< scalarField > deltaZ_
Definition ISQP.H:172
scalar relTol_
Relative tolerance of the CG solver, solving for deltaP_.
Definition ISQP.H:207
autoPtr< scalarField > extraVars_
Extra variables for finding solutions even in infeasible problems.
Definition ISQP.H:148
autoPtr< scalarField > ls_
Slack variables the lower bound constraints.
Definition ISQP.H:132
tmp< scalarField > diagPreconditioner()
Diagonal preconditioner of the deltaP eqn.
Definition ISQP.C:389
tmp< scalarField > resFlamda()
Residual of the complementarity slackness for the inequality constraints.
Definition ISQP.C:810
bool includeBoundConstraints_
Are bound constraints included?
Definition ISQP.H:98
bool disableDamping_
Disable damping.
Definition ISQP.H:227
autoPtr< scalarField > lTilda_
Lagrange multipliers for the lower bound constraints.
Definition ISQP.H:127
virtual void updateOldCorrection(const scalarField &oldCorrection)
Update old correction. Useful for quasi-Newton methods coupled with line search.
Definition ISQP.C:1152
void updateYS()
Update the vectors accosiated with the Hessian matrix.
Definition ISQP.C:136
virtual ~ISQP()=default
Destructor.
tmp< scalarField > resFlTilda()
Residual of the complementarity slackness for the lower bound constraints.
Definition ISQP.C:844
autoPtr< scalarField > uTilda_
Lagrange multipliers for the upper bound constraints.
Definition ISQP.H:137
TypeName("ISQP")
Runtime type information.
void storeOldFields()
Store old fields needed for the next iter.
Definition ISQP.C:976
tmp< scalarField > computeResiduals()
Compute and return residuals based on the current solution.
Definition ISQP.C:667
void updateSizes()
Update sizes of fields related to the active design variables.
Definition ISQP.C:58
void solveDeltaPEqn()
Solve the equation for deltaX, which is the expensive part of the Newtopn step.
Definition ISQP.C:215
autoPtr< scalarField > deltaLTilda_
Definition ISQP.H:167
scalarField gs_
Slack variables for the inequality constraints.
Definition ISQP.H:122
label maxNewtonIters_
Maxmimum number of Newton iterations for the subproblem.
Definition ISQP.H:191
label preconType_
Which preconditioner to use for the solution of the subproblem.
Definition ISQP.H:212
tmp< scalarField > ShermanMorrisonPrecon(const scalarField &vector)
Use the Sherman-Morrison formula to compute the matrix-free product of the approximate inverse of the...
Definition ISQP.C:419
void allocateLagrangeMultipliers()
Allocate Lagrange multipliers for the inequality constraints.
Definition ISQP.C:117
tmp< scalarField > computeRHSForDeltaX(const scalarField &FDx)
Compute the RHS for the deltaX equation.
Definition ISQP.C:247
virtual void update()
Update design variables.
Definition LBFGS.C:590
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
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
A class for managing references or pointers (no reference counting).
Definition refPtr.H:54
A class for managing temporary objects.
Definition tmp.H:75
label nConstraints() const
Get the number of constraints.
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
dynamicFvMesh & mesh
OBJstream os(runTime.globalPath()/outputName)
Namespace for OpenFOAM.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Vector< scalar > vector
Definition vector.H:57
dictionary dict
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition typeInfo.H:68