KASKADE 7 development version
|
An iterative solver for small instances of bound constrained quadratic programs. More...
#include <qp.hh>
An iterative solver for small instances of bound constrained quadratic programs.
The solver addresses problems of the form
\[ \min_x \frac{1}{2} x^T Ax + c^T x \quad \text{s.t. } x \le b, \]
where \( A \) is positive definite. The implementation is intended to solve small problems of just a couple of variables, typically less than 30.
The actual matrix operations have to be provided by a derived class.
Real | the real field type of matrix/vector entries. |
Implementation | the type of the derived class (Barton-Nackman trick) The implementation class needs to overwrite the solveSubmatrix method. |
A simple projected gradient method is used.
Public Types | |
using | Vector = Dune::BlockVector< Dune::FieldVector< Real, 1 > > |
using | Matrix = DynamicMatrix< Dune::FieldMatrix< Real, 1, 1 > > |
Public Member Functions | |
virtual | ~QPBoundSolverBase () |
std::tuple< Vector, int > | solve (Vector const &c, Vector const &b, double tol) const |
Solves the QP for the given right hand side vectors. More... | |
std::tuple< Vector, int > | solve (Vector const &x, Vector const &c, Vector const &b, double tol) const |
Solves the QP for the given right hand side vectors. More... | |
void | gradientStep (Vector const &c, Vector const &b, Vector &x) const |
Performs a single projected gradient step. More... | |
Implementation & | setMinSteps (int count) |
defines the minimal number of steps to perform in solving. More... | |
Implementation & | setMaxSteps (int count) |
defines the maximum allowed number of steps to perform in solving. More... | |
Protected Member Functions | |
QPBoundSolverBase (int dim) | |
Constructor (to be used by derived classes only) More... | |
Vector | solveSubmatrix (std::vector< int > const &idx, Vector const &b) const |
solves a subsystem \( \tilde A x = b \). More... | |
virtual void | updateWithColumn (Real a, int i, Vector &s) const |
Performs the update \( s \gets s + aAe_i \), where \( e_i \) is the i-th unit vector. More... | |
using Kaskade::QPBoundSolverBase< Real, Implementation >::Matrix = DynamicMatrix<Dune::FieldMatrix<Real,1,1> > |
using Kaskade::QPBoundSolverBase< Real, Implementation >::Vector = Dune::BlockVector<Dune::FieldVector<Real,1> > |
|
inlinevirtual |
|
inlineprotected |
void Kaskade::QPBoundSolverBase< Real, Implementation >::gradientStep | ( | Vector const & | c, |
Vector const & | b, | ||
Vector & | x | ||
) | const |
Performs a single projected gradient step.
This is a generalization of a Richardson smoother with linesearch, and can be used when the QP is just a subproblem, e.g. in multigrid methods, and rather well conditioned.
Implementation & Kaskade::QPBoundSolverBase< Real, Implementation >::setMaxSteps | ( | int | count | ) |
defines the maximum allowed number of steps to perform in solving.
The default value is 200. If the maximum iteration count is below the minimum iteration count, the number of steps performed in solving is undefined.
Implementation & Kaskade::QPBoundSolverBase< Real, Implementation >::setMinSteps | ( | int | count | ) |
defines the minimal number of steps to perform in solving.
The default value is 0. If the minimum iteration count is set to less than the maximum iteration count, the number of steps performed is undefined.
std::tuple< Vector, int > Kaskade::QPBoundSolverBase< Real, Implementation >::solve | ( | Vector const & | c, |
Vector const & | b, | ||
double | tol | ||
) | const |
Solves the QP for the given right hand side vectors.
Solves the QP with a projected gradient method starting at the feasible point \( x = \min(0,b)\).
c | the linear term of the objective |
b | the constraints right hand side |
tol | the tolerance for termination criterion |
The iteration is terminated as soon as the projected gradient \( s \) is small, i.e. \( \|s\|_2 \le \mathsf{tol} \).
If the QP is infeasible, an InfeasibleProblemException is thrown.
std::tuple< Vector, int > Kaskade::QPBoundSolverBase< Real, Implementation >::solve | ( | Vector const & | x, |
Vector const & | c, | ||
Vector const & | b, | ||
double | tol | ||
) | const |
Solves the QP for the given right hand side vectors.
x | an initial guess for the solution. If an unfeasible starting point is provided, it will be projected to the feasible set. |
c | the linear term of the objective |
b | the constraints right hand side |
tol | the tolerance for termination criterion |
The iteration is terminated as soon as the projected gradient \( s \) is small, i.e. \( \|s\|_2 \le \mathsf{tol} \).
If the QP is infeasible, an InfeasibleProblemException is thrown.
|
protected |
solves a subsystem \( \tilde A x = b \).
This needs to be overwritten by derived classes.
For a set idx of indices, the submatrix \( \tilde A \) in Matlab notation is A(idx,idx).
idx | the indices defining the submatrix |
b | the right hand side (of the same length as the index vector) |
|
protectedvirtual |
Performs the update \( s \gets s + aAe_i \), where \( e_i \) is the i-th unit vector.
The default implementation performs a complete matrix-vector product, which is performance-wise suboptimal. Override this for better performance.