KASKADE 7 development version
|
An augmented Lagrangian solver for small instances of a particular class of quadratic programs. More...
#include <qp.hh>
An augmented Lagrangian solver for small instances of a particular class of quadratic programs.
The solver addresses problems of the form
\[ \min_x \frac{1}{2} x^T Ax + c^T x \quad \text{s.t. } Bx \le b, \]
where \( A \) is positive semidefinite, and approximates the solution by minimizing the Powell-Hestenes-Rockafellar augmented Lagrangian
\[ \min_{x} \frac{1}{2} x^T Ax + c^T x + \frac{\gamma}{2}\left\|Bx-\Big(b-\frac{\lambda}{\gamma}\Big)\right\|_+^2 \]
using a semismooth Newton method with exact linesearch and first-order multiplier updates.
The penalty parameter \( \gamma \) affects both, solution quality and solution effort. For \( \gamma\to\infty \), the augmented Lagrangian iteration converges faster, but the penalized problem's condition gets worse. A reasonable default value is chosen based on theoretical estimates of the convergence rate, aiming at a contraction factor of 0.01 (a very rough estimate, but should work as a default).
A good multiplier estimate \(\lambda\) improves the solution accuracy significantly. Compared to the standard augmented Lagrangian where a linear term \( \lambda^T (Bx-b) \) is added, the PHR version with shiftet constraints has the advantage that in the inner penalized QPs that are to be solved the strongly active constraints remain strongly active. In the standard formulation, they become weakly active in the course of the AL iteration, which can induce frequent switching of constraints between active and inactive.
d | the dimension of (vector valued) optimization variables |
sparsity | the type of linear algebra data structures to use, either DENSE or SPARSE |
Real | the real field type of matrix/vector entries. Explicit instantiations are provided for float and double. |
Public Types | |
using | VectorX = Dune::BlockVector< Dune::FieldVector< Real, d > > |
using | VectorB = Dune::BlockVector< Dune::FieldVector< Real, 1 > > |
using | MatrixA = std::conditional_t< sparsity==QPStructure::DENSE, DynamicMatrix< AEntry >, NumaBCRSMatrix< AEntry > > |
using | MatrixB = std::conditional_t< sparsity==QPStructure::DENSE, DynamicMatrix< BEntry >, NumaBCRSMatrix< BEntry > > |
Public Member Functions | |
QPALSolver (MatrixA const &A, MatrixB const &B) | |
Constructs the solver, providing the matrices \( A \) and \( B \). More... | |
QPALSolver< d, sparsity, Real > & | setMinSteps (int i) |
Sets the minimum number of augmented Lagrangian steps to perform. More... | |
QPALSolver< d, sparsity, Real > & | setMaxSteps (int i) |
Sets the maximum number of projected gradient steps to perform. More... | |
QPALSolver< d, sparsity, Real > & | setPenalty (Real gamma) |
Sets the penalty parameter \( \gamma \). More... | |
Real | getPenalty () const |
std::tuple< VectorX, VectorB, int > | solve (VectorX const &c, VectorB const &b, double tol) const |
Solves the QP approximately for the given right hand side vectors. More... | |
std::tuple< VectorX, VectorB, int > | solve (VectorX const &c, VectorB const &b, double tol, VectorX x) const |
Solves the QP approximately for the given right hand side vectors. More... | |
Static Public Attributes | |
static constexpr QPStructure | sparse = sparsity |
Given info on whether dense or sparse arithmetic is used (i.e. the choice of template parameter used in this class). More... | |
using Kaskade::QPALSolver< d, sparsity, Real >::MatrixA = std::conditional_t<sparsity==QPStructure::DENSE, DynamicMatrix<AEntry>, NumaBCRSMatrix<AEntry> > |
using Kaskade::QPALSolver< d, sparsity, Real >::MatrixB = std::conditional_t<sparsity==QPStructure::DENSE, DynamicMatrix<BEntry>, NumaBCRSMatrix<BEntry> > |
using Kaskade::QPALSolver< d, sparsity, Real >::VectorB = Dune::BlockVector<Dune::FieldVector<Real,1> > |
using Kaskade::QPALSolver< d, sparsity, Real >::VectorX = Dune::BlockVector<Dune::FieldVector<Real,d> > |
Kaskade::QPALSolver< d, sparsity, Real >::QPALSolver | ( | MatrixA const & | A, |
MatrixB const & | B | ||
) |
Constructs the solver, providing the matrices \( A \) and \( B \).
The argument matrices are copied internally, and are not referenced later.
|
inline |
QPALSolver< d, sparsity, Real > & Kaskade::QPALSolver< d, sparsity, Real >::setMaxSteps | ( | int | i | ) |
Sets the maximum number of projected gradient steps to perform.
i | the maximum number of steps, needs to be nonnegative |
If i is less than the minimum number of steps, the latter ist decreased to the same value.
QPALSolver< d, sparsity, Real > & Kaskade::QPALSolver< d, sparsity, Real >::setMinSteps | ( | int | i | ) |
Sets the minimum number of augmented Lagrangian steps to perform.
i | the minimum number of steps, needs to be nonnegative |
If i is larger than the maximum number of steps, the latter ist increased to the same value.
QPALSolver< d, sparsity, Real > & Kaskade::QPALSolver< d, sparsity, Real >::setPenalty | ( | Real | gamma | ) |
Sets the penalty parameter \( \gamma \).
The default on construction is \( \gamma = \).
std::tuple< VectorX, VectorB, int > Kaskade::QPALSolver< d, sparsity, Real >::solve | ( | VectorX const & | c, |
VectorB const & | b, | ||
double | tol | ||
) | const |
Solves the QP approximately for the given right hand side vectors.
The starting point is \( x = 0 \).
c | the linear term of the objective |
b | the constraints right hand side |
tol | the tolerance for termination |
The iteration is terminated as soon as the projected gradient \( s \) is small, i.e. \( \|s\|_2 \le \mathsf{tol}\). Together with the solution vector, a multiplier estimate based on the first order update \( \lambda \gets \lambda - (Bx-s) \) is computed.
If the QP is infeasible, an InfeasibleProblemException is thrown.
std::tuple< VectorX, VectorB, int > Kaskade::QPALSolver< d, sparsity, Real >::solve | ( | VectorX const & | c, |
VectorB const & | b, | ||
double | tol, | ||
VectorX | x | ||
) | const |
Solves the QP approximately for the given right hand side vectors.
x | the start iterate |
c | the linear term of the objective |
b | the constraints right hand side |
tol | the tolerance for termination |
The iteration is terminated as soon as the projected gradient \( s \) is small, i.e. \( \|s\|_2 \le \mathsf{tol}\). Together with the solution vector, a multiplier estimate based on the first order update \( \lambda \gets \lambda - (Bx-s) \) is computed.
If the QP is infeasible, an InfeasibleProblemException is thrown.
|
staticconstexpr |