KASKADE 7 development version
|
An iterative solver for small instances of a particular class of quadratic programs. More...
#include <qp.hh>
An iterative 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 definite. The implementation works with dense matrices and is intended to solve small problems of just a couple of variables, typically less than 30, and moreover well conditioned problems.
Real | the real field type of matrix/vector entries. Explicit instantiations are provided for float and double. |
d | the dimension of (vector valued) optimization variables |
The solver forms the Lagrange dual problem
\[ \max_\lambda - \frac{1}{2} \lambda^TBA^{-1}B^T\lambda + (BA^{-1}c+b)^T \lambda \quad\text{s.t. } \lambda \le 0, \]
which is bound constrained and can be solved by a projected gradient method. As the primal problem is convex, strong duality holds, and the unique primal minimizer \( x^* \) can be computed from a dual maximizer \( \lambda^* \) by
\[ x^* = A^{-1}(B^T\lambda^*-c). \]
Public Types | |
using | VectorX = Dune::BlockVector< Dune::FieldVector< Real, d > > |
using | VectorB = Dune::BlockVector< Dune::FieldVector< Real, 1 > > |
using | MatrixA = DynamicMatrix< Dune::FieldMatrix< Real, d, d > > |
using | MatrixB = DynamicMatrix< Dune::FieldMatrix< Real, 1, d > > |
using | Self = QPSolver< d, Real > |
Public Member Functions | |
QPSolver (MatrixA const &A, MatrixB const &B, QPConvexificationStrategy convex=QPConvexificationStrategy::DONOTHING) | |
Constructs the solver, providing the matrices \( A \) and \( B \). More... | |
QPSolver (Self const &qp) | |
Self & | operator= (Self const &qp) |
std::tuple< VectorX, VectorB, int > | solve (VectorX const &c, VectorB const &b, double tol) const |
Solves the QP for the given right hand side vectors. More... | |
MatrixB const & | matrixB () const |
Provides access to the internal copy of \( B\). More... | |
Static Public Attributes | |
static constexpr QPStructure | sparse = QPStructure::DENSE |
Given info on whether dense or sparse arithmetic is used. More... | |
using Kaskade::QPSolver< d, Real >::MatrixA = DynamicMatrix<Dune::FieldMatrix<Real,d,d> > |
using Kaskade::QPSolver< d, Real >::MatrixB = DynamicMatrix<Dune::FieldMatrix<Real,1,d> > |
using Kaskade::QPSolver< d, Real >::Self = QPSolver<d,Real> |
using Kaskade::QPSolver< d, Real >::VectorB = Dune::BlockVector<Dune::FieldVector<Real,1> > |
using Kaskade::QPSolver< d, Real >::VectorX = Dune::BlockVector<Dune::FieldVector<Real,d> > |
Kaskade::QPSolver< d, Real >::QPSolver | ( | MatrixA const & | A, |
MatrixB const & | B, | ||
QPConvexificationStrategy | convex = QPConvexificationStrategy::DONOTHING |
||
) |
Constructs the solver, providing the matrices \( A \) and \( B \).
The argument matrices are copied internally, and are not referenced later.
Throws a NonpositiveMatrixException if \( A \) is not positive definite and convex
is QPConvexificationStrategy::DONOTHING.
Kaskade::QPSolver< d, Real >::QPSolver | ( | Self const & | qp | ) |
|
inline |
Self & Kaskade::QPSolver< d, Real >::operator= | ( | Self const & | qp | ) |
std::tuple< VectorX, VectorB, int > Kaskade::QPSolver< d, Real >::solve | ( | VectorX const & | c, |
VectorB const & | b, | ||
double | tol | ||
) | const |
Solves the QP for the given right hand side vectors.
Solves the QP by solving the dual QP with a projected gradient method starting at the dually feasible point \( \lambda = 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 \) of the dual problem is small, i.e. \( \|s\|_2 \le \mathsf{tol}\).
If the QP is infeasible, an InfeasibleProblemException is thrown.
|
staticconstexpr |