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 -------------------------------------------------------------------------------
11 License
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 
28 Class
29  Foam::ISQP
30 
31 Description
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 
50 SourceFiles
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 
63 namespace Foam
64 {
65 
66 /*---------------------------------------------------------------------------*\
67  Class ISQP Declaration
68 \*---------------------------------------------------------------------------*/
69 
70 class ISQP
71 :
72  public LBFGS,
73  public SQPBase
74 {
75 public:
76 
77  //- Enumeration of preconditioners for the subproblem
78  enum preconditioners
79  {
81  };
82 
83  //- Names of preconditioners for the subproblem
85 
86 
87 protected:
88 
89  // Protected data
90 
91  //- Should Lagrange multipliers be allocated
93 
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
100  bool includeExtraVars_;
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'
134  scalar c_;
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
163  label maxLineSearchIters_;
165  //- Maxmimum number of iterations for the deltaX equation
166  label maxDxIters_;
168  //- Relative tolerance of the CG solver, solving for deltaP_
169  scalar relTol_;
171  //- Which preconditioner to use for the solution of the subproblem
172  label preconType_;
173 
174  //- Percentage reduction for the constraints.
175  // Can be used to relax QP problems with no feasible points
176  scalar cRed_;
177 
178  //- File including the l1 merit function
181 
182  // Protected Member Functions
183 
184  //- Update sizes of fields related to the active design variables
185  void updateSizes();
187  //- Allocate multipliers for the bound constraints
189 
190  //- Allocate Lagrange multipliers for the inequality constraints
192 
193  //- Update the vectors accosiated with the Hessian matrix
194  void updateYS();
195 
196  //- Allocate fields related to constraints
197  void initialize();
198 
199 
200  // Functions related to the solution of the primal-dual subproblem
201 
202  //- Solve subproblem using a Newton optimiser
203  void solveSubproblem();
204 
205  //- Compute direction for the Newton problem
206  void computeNewtonDirection();
208  //- Zero the updates computed in the previous optimisation cycle
209  void zeroUpdates();
210 
211  //- Solve the equation for deltaX, which is the expensive part of
212  //- the Newtopn step.
213  // All other corrections can be computed based on this
214  void solveDeltaPEqn();
215 
216  //- Compute the RHS for the deltaX equation
217  // Currently not in use
220  //- CG algorithm for the solution of the deltaP eqn
221  void CGforDeltaP(const scalarField& FDx);
222 
223  //- Procudt of the LHS of the deltaP eqn with a vector
225  (
226  const scalarField& vector
227  );
228 
229  //- Preconditioner-vector product for CG
230  // In case a diagonal preconditioner is used, it is stored in
231  // precon for all CG iterations
233  (
234  const scalarField& vector,
235  autoPtr<scalarField>& precon
236  );
237 
238  //- Diagonal preconditioner of the deltaP eqn
240 
241  //- Use the Sherman-Morrison formula to compute the matrix-free
242  //- product of the approximate inverse of the LHS with a vector
244 
245  //- Recursive (naive) implementation of the rank-1 update
246  // WIP, efficient for up to 2 flow-related constraints
248  (
249  const PtrList<scalarField>& r1Updates,
250  const scalarField& p,
251  const tmp<scalarField>& diag,
252  const scalarField& mult,
253  label n
254  );
255 
256  //- Perform line search and return max residual corresponding to
257  //- the updated solution
258  scalar lineSearch();
259 
260  //- Adjust step to satisfy cireteria
261  void adjustStep
262  (
263  scalar& step,
264  const scalar value,
265  const scalar update
266  );
267 
268  //- Update the current solution using the known direction and the
269  //- given step length
270  void updateSolution(const scalar step);
271 
272 
273  // Residuals of the various KKT conditions
274 
275  //- Compute and return residuals based on the current solution
277 
278  //- Residual of the gradient of the Lagrangian
279  // Size is that of the active design variables
281 
282  //- Product of the inverse Hessian with the residual of the
283  //- gradient of the Lagrangian.
284  // Avoid the formation of the Hessian matrix.
285  // Size is that of the active design variables.
287 
288  //- Residual of the inequality constraint slackness
290 
291  //- Residual of the complementarity slackness for the
292  //- inequality constraints
294 
295  //- Residual of the lower bounds slackness
297 
298  //- Residual of the upper bounds slackness
300 
301  //- Residual of the complementarity slackness for the
302  //- lower bound constraints
304 
305  //- Residual of the complementarity slackness for the
306  //- upper bound constraints
308 
309  //- Residual of the Lagrangian derivative wrt the extra variables
311 
312  //- Residual of the complementarity slackness for the
313  //- extra variables
315 
316 
317  //- Store old fields needed for the next iter
318  void storeOldFields();
319 
320  //- Get the part the merit function that depends on the constraints
321  virtual scalar meritFunctionConstraintPart() const;
322 
323 
324 private:
325 
326  // Private Member Functions
327 
328  //- No copy construct
329  ISQP(const ISQP&) = delete;
330 
331  //- No copy assignment
332  void operator=(const ISQP&) = delete;
333 
334 
335 public:
336 
337  //- Runtime type information
338  TypeName("ISQP");
339 
340 
341  // Constructors
342 
343  //- Construct from components
344  ISQP
345  (
346  const fvMesh& mesh,
347  const dictionary& dict,
348  autoPtr<designVariables>& designVars,
349  const label nConstraints,
350  const word& type
351  );
352 
353 
354  //- Destructor
355  virtual ~ISQP() = default;
356 
357 
358  // Member Functions
359 
360  //- Compute design variables correction
361  void computeCorrection();
362 
363  //- Compute merit function. Could be different than the objective
364  //- in the presence of constraints
365  virtual scalar computeMeritFunction();
366 
367  //- Derivative of the merit function. Could be different than the
368  //- objective derivative in the presence of constraints
369  virtual scalar meritFunctionDirectionalDerivative();
370 
371  //- Update old correction. Useful for quasi-Newton methods coupled with
372  //- line search
373  virtual void updateOldCorrection(const scalarField& oldCorrection);
374 
375  //- Write useful quantities to files
376  virtual bool writeData(Ostream& os) const;
377 
378  //- Write merit function information
379  virtual bool writeAuxiliaryData();
380 };
381 
382 
383 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
384 
385 } // End namespace Foam
386 
387 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
388 
389 #endif
390 
391 // ************************************************************************* //
tmp< scalarField > preconVectorProduct(const scalarField &vector, autoPtr< scalarField > &precon)
Preconditioner-vector product for CG.
Definition: ISQP.C:359
tmp< scalarField > resFuTilda()
Residual of the complementarity slackness for the upper bound constraints.
Definition: ISQP.C:850
tmp< scalarField > diagPreconditioner()
Diagonal preconditioner of the deltaP eqn.
Definition: ISQP.C:385
dictionary dict
virtual void updateOldCorrection(const scalarField &oldCorrection)
Update old correction. Useful for quasi-Newton methods coupled with line search.
Definition: ISQP.C:1122
scalarField p_
The set of design variables being updated during the subproblem.
Definition: ISQP.H:115
label maxLineSearchIters_
Maxmimum number of line search iterations for each iteration of the subproblem.
Definition: ISQP.H:197
scalarField deltaP_
Definition: ISQP.H:164
virtual void update()
Update design variables.
Definition: LBFGS.C:590
label nConstraints() const
Get the number of constraints.
Definition: updateMethod.C:393
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
void updateYS()
Update the vectors accosiated with the Hessian matrix.
Definition: ISQP.C:134
void initialize()
Allocate fields related to constraints.
Definition: ISQP.C:160
tmp< scalarField > ShermanMorrisonRank1Update(const PtrList< scalarField > &r1Updates, const scalarField &p, const tmp< scalarField > &diag, const scalarField &mult, label n)
Recursive (naive) implementation of the rank-1 update.
Definition: ISQP.C:455
virtual scalar computeMeritFunction()
Compute merit function. Could be different than the objective in the presence of constraints.
Definition: ISQP.C:1103
autoPtr< scalarField > deltaUTilda_
Definition: ISQP.H:169
bool doAllocateLagrangeMultipliers_
Should Lagrange multipliers be allocated.
Definition: ISQP.H:93
tmp< scalarField > resFlTilda()
Residual of the complementarity slackness for the lower bound constraints.
Definition: ISQP.C:840
autoPtr< scalarField > uTilda_
Lagrange multipliers for the upper bound constraints.
Definition: ISQP.H:137
bool includeBoundConstraints_
Are bound constraints included?
Definition: ISQP.H:98
autoPtr< scalarField > lTilda_
Lagrange multipliers for the lower bound constraints.
Definition: ISQP.H:127
scalar c_
Multiplier of the additional variables y in the Lagrangian, to make them &#39;expensive&#39;.
Definition: ISQP.H:159
virtual bool writeData(Ostream &os) const
Write useful quantities to files.
Definition: ISQP.C:1129
virtual scalar meritFunctionDirectionalDerivative()
Derivative of the merit function. Could be different than the objective derivative in the presence of...
Definition: ISQP.C:1112
label preconType_
Which preconditioner to use for the solution of the subproblem.
Definition: ISQP.H:212
tmp< scalarField > computeRHSForDeltaX(const scalarField &FDx)
Compute the RHS for the deltaX equation.
Definition: ISQP.C:243
autoPtr< scalarField > deltaZ_
Definition: ISQP.H:172
void computeCorrection()
Compute design variables correction.
Definition: ISQP.C:1077
autoPtr< scalarField > extraVars_
Extra variables for finding solutions even in infeasible problems.
Definition: ISQP.H:148
autoPtr< scalarField > deltaLs_
Definition: ISQP.H:168
void adjustStep(scalar &step, const scalar value, const scalar update)
Adjust step to satisfy cireteria.
Definition: ISQP.C:630
virtual scalar meritFunctionConstraintPart() const
Get the part the merit function that depends on the constraints.
Definition: ISQP.C:981
void CGforDeltaP(const scalarField &FDx)
CG algorithm for the solution of the deltaP eqn.
Definition: ISQP.C:283
virtual ~ISQP()=default
Destructor.
tmp< scalarField > resFExtraVars()
Residual of the Lagrangian derivative wrt the extra variables.
Definition: ISQP.C:860
tmp< scalarField > resFL()
Residual of the gradient of the Lagrangian.
Definition: ISQP.C:723
void updateSizes()
Update sizes of fields related to the active design variables.
Definition: ISQP.C:57
tmp< scalarField > resFlamda()
Residual of the complementarity slackness for the inequality constraints.
Definition: ISQP.C:806
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
autoPtr< scalarField > deltaExtraVars_
Definition: ISQP.H:171
autoPtr< scalarField > deltaLTilda_
Definition: ISQP.H:167
scalar cRed_
Percentage reduction for the constraints.
Definition: ISQP.H:219
label maxDxIters_
Maxmimum number of iterations for the deltaX equation.
Definition: ISQP.H:202
dynamicFvMesh & mesh
void solveDeltaPEqn()
Solve the equation for deltaX, which is the expensive part of the Newtopn step.
Definition: ISQP.C:211
void computeNewtonDirection()
Compute direction for the Newton problem.
Definition: ISQP.C:486
A class for handling words, derived from Foam::string.
Definition: word.H:63
scalarField deltaLamda_
Definition: ISQP.H:165
scalarField deltaGs_
Definition: ISQP.H:166
tmp< scalarField > resFls()
Residual of the lower bounds slackness.
Definition: ISQP.C:812
scalar lineSearch()
Perform line search and return max residual corresponding to the updated solution.
Definition: ISQP.C:538
tmp< scalarField > computeResiduals()
Compute and return residuals based on the current solution.
Definition: ISQP.C:663
Base class for Sequantial Quadratic Programming (SQP) methods.
Definition: SQPBase.H:50
void allocateBoundMultipliers()
Allocate multipliers for the bound constraints.
Definition: ISQP.C:90
autoPtr< scalarField > z_
Lagrange multipliers for positive extra variables.
Definition: ISQP.H:153
tmp< scalarField > resFz()
Residual of the complementarity slackness for the extra variables.
Definition: ISQP.C:870
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
tmp< scalarField > resFGs()
Residual of the inequality constraint slackness.
Definition: ISQP.C:780
OBJstream os(runTime.globalPath()/outputName)
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:415
void solveSubproblem()
Solve subproblem using a Newton optimiser.
Definition: ISQP.C:880
An L-BFGS-based SQP algorithm for computing the update of the design variables in the presence of ine...
Definition: ISQP.H:65
void updateSolution(const scalar step)
Update the current solution using the known direction and the given step length.
Definition: ISQP.C:643
scalarField gs_
Slack variables for the inequality constraints.
Definition: ISQP.H:122
static const Enum< preconditioners > preconditionerNames
Names of preconditioners for the subproblem.
Definition: ISQP.H:83
virtual bool writeAuxiliaryData()
Write merit function information.
Definition: ISQP.C:1141
TypeName("ISQP")
Runtime type information.
autoPtr< scalarField > ls_
Slack variables the lower bound constraints.
Definition: ISQP.H:132
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
bool includeExtraVars_
Are additional design variables included?
Definition: ISQP.H:106
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
scalar eps_
Infinitesimal quantity.
Definition: ISQP.H:180
scalar epsMin_
Final eps quantity to be reached during the solution of the subproblem.
Definition: ISQP.H:186
label maxNewtonIters_
Maxmimum number of Newton iterations for the subproblem.
Definition: ISQP.H:191
scalar relTol_
Relative tolerance of the CG solver, solving for deltaP_.
Definition: ISQP.H:207
tmp< scalarField > matrixVectorProduct(const scalarField &vector)
Procudt of the LHS of the deltaP eqn with a vector.
Definition: ISQP.C:316
The quasi-Newton Limited-memory BFGS formula. Keeps nPrevSteps_ of the y and s vectors and approximat...
Definition: LBFGS.H:50
label n
autoPtr< scalarField > us_
Slack variables the upper bound constraints.
Definition: ISQP.H:142
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
tmp< scalarField > resFus()
Residual of the upper bounds slackness.
Definition: ISQP.C:826
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
tmp< scalarField > invHFL()
Product of the inverse Hessian with the residual of the gradient of the Lagrangian.
Definition: ISQP.C:755
void storeOldFields()
Store old fields needed for the next iter.
Definition: ISQP.C:966
autoPtr< scalarField > deltaUs_
Definition: ISQP.H:170
void allocateLagrangeMultipliers()
Allocate Lagrange multipliers for the inequality constraints.
Definition: ISQP.C:116
void zeroUpdates()
Zero the updates computed in the previous optimisation cycle.
Definition: ISQP.C:189
preconditioners
Enumeration of preconditioners for the subproblem.
Definition: ISQP.H:75
autoPtr< OFstream > meritFunctionFile_
File including the l1 merit function.
Definition: ISQP.H:224
Namespace for OpenFOAM.