KASKADE 7 development version
Public Types | Public Member Functions | Protected Member Functions | List of all members
Kaskade::QPSlackSolver< sparsity, R > Class Template Reference

An iterative solver for particular instances of bound constrained quadratic programs. More...

#include <qp.hh>

Detailed Description

template<QPStructure sparsity = QPStructure::DENSE, class R = double>
class Kaskade::QPSlackSolver< sparsity, R >

An iterative solver for particular instances of bound constrained quadratic programs.

The solver addresses problems of the form

\[ \min_x \frac{1}{2} x^T Qx + c^T x \quad \text{s.t. } x \le b, \]

where

\[ Q = \begin{bmatrix} A+\gamma B^T B & -\gamma B^T \\ -\gamma B & \gamma I \end{bmatrix}\]

is a positive definite block structured matrix. Problems of this type rarely arise on their own, mostly from a penalty or augmented Lagrangian relaxation of general QPs (see QPPenalizedSolver).

The implementation is intended to solve problems where \( A \) is small and dense , with just a couple of variables, typically less than 30.

Template Parameters
structurespecifies which type of linear algebra to use:
  • DENSE uses dense linear algebra and is intended to solve problems where \( A \) is small, with usually not more than 30 variables
  • SPARSE uses sparse solvers, and is intended for problems of moderate size
Realthe real field type of matrix/vector entries. Explicit instantiations are provided for float and double.

Definition at line 356 of file qp.hh.

Inheritance diagram for Kaskade::QPSlackSolver< sparsity, R >:
Kaskade::QPBoundSolverBase< double, QPSlackSolver< QPStructure::DENSE, double > >

Public Types

using Real = R
 
using Vector = typename Base::Vector
 
using Matrix = std::conditional_t< sparsity==QPStructure::DENSE, DynamicMatrix< Entry >, NumaBCRSMatrix< Entry > >
 

Public Member Functions

 QPSlackSolver (Matrix const &A, Matrix const &B, Real gamma, QPConvexificationStrategy convex=QPConvexificationStrategy::DONOTHING)
 Constructs the solver, providing the matrix \( A \) and \( B \). More...
 
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...
 
QPSlackSolver< QPStructure::DENSE, double > & setMinSteps (int count)
 defines the minimal number of steps to perform in solving. More...
 
QPSlackSolver< QPStructure::DENSE, double > & setMaxSteps (int count)
 defines the maximum allowed number of steps to perform in solving. More...
 

Protected Member Functions

Vector solveSubmatrix (std::vector< int > const &idx, Vector const &b) const
 solves a subsystem \( \tilde A x = b \). More...
 
virtual void updateWithColumn (double a, int i, Vector &s) const
 Performs the update \( s \gets s + aAe_i \), where \( e_i \) is the i-th unit vector. More...
 

Member Typedef Documentation

◆ Matrix

template<QPStructure sparsity = QPStructure::DENSE, class R = double>
using Kaskade::QPSlackSolver< sparsity, R >::Matrix = std::conditional_t<sparsity==QPStructure::DENSE, DynamicMatrix<Entry>, NumaBCRSMatrix<Entry> >

Definition at line 364 of file qp.hh.

◆ Real

template<QPStructure sparsity = QPStructure::DENSE, class R = double>
using Kaskade::QPSlackSolver< sparsity, R >::Real = R

Definition at line 362 of file qp.hh.

◆ Vector

template<QPStructure sparsity = QPStructure::DENSE, class R = double>
using Kaskade::QPSlackSolver< sparsity, R >::Vector = typename Base::Vector

Definition at line 363 of file qp.hh.

Constructor & Destructor Documentation

◆ QPSlackSolver()

template<QPStructure sparsity = QPStructure::DENSE, class R = double>
Kaskade::QPSlackSolver< sparsity, R >::QPSlackSolver ( Matrix const &  A,
Matrix const &  B,
Real  gamma,
QPConvexificationStrategy  convex = QPConvexificationStrategy::DONOTHING 
)

Constructs the solver, providing the matrix \( A \) and \( B \).

The argument matrices are copied internally, and are not referenced later.

Throws a NonpositiveMatrixException if \( A \) is not positive definite and

Todo:
write a move constructor

Member Function Documentation

◆ gradientStep()

void Kaskade::QPBoundSolverBase< double , QPSlackSolver< QPStructure::DENSE, double > >::gradientStep ( Vector const &  c,
Vector const &  b,
Vector x 
) const
inherited

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.

◆ setMaxSteps()

QPSlackSolver< QPStructure::DENSE, double > & Kaskade::QPBoundSolverBase< double , QPSlackSolver< QPStructure::DENSE, double > >::setMaxSteps ( int  count)
inherited

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.

◆ setMinSteps()

QPSlackSolver< QPStructure::DENSE, double > & Kaskade::QPBoundSolverBase< double , QPSlackSolver< QPStructure::DENSE, double > >::setMinSteps ( int  count)
inherited

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.

◆ solve() [1/2]

std::tuple< Vector, int > Kaskade::QPBoundSolverBase< double , QPSlackSolver< QPStructure::DENSE, double > >::solve ( Vector const &  c,
Vector const &  b,
double  tol 
) const
inherited

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)\).

Parameters
cthe linear term of the objective
bthe constraints right hand side
tolthe tolerance for termination criterion

The iteration is terminated as soon as the projected gradient \( s \) is small, i.e. \( \|s\|_2 \le \mathsf{tol} \).

Returns
a tuple \( (x,iter) \) of solution vector, and iteration count

If the QP is infeasible, an InfeasibleProblemException is thrown.

◆ solve() [2/2]

std::tuple< Vector, int > Kaskade::QPBoundSolverBase< double , QPSlackSolver< QPStructure::DENSE, double > >::solve ( Vector const &  x,
Vector const &  c,
Vector const &  b,
double  tol 
) const
inherited

Solves the QP for the given right hand side vectors.

Parameters
xan initial guess for the solution. If an unfeasible starting point is provided, it will be projected to the feasible set.
cthe linear term of the objective
bthe constraints right hand side
tolthe tolerance for termination criterion

The iteration is terminated as soon as the projected gradient \( s \) is small, i.e. \( \|s\|_2 \le \mathsf{tol} \).

Returns
a tuple \( (x,iter) \) of solution vector, and iteration count. If the iteration count equals the maximum iteration count (
See also
setMaxSteps), the returned Solution Tools vector need not satisfy the tolerance criterion.

If the QP is infeasible, an InfeasibleProblemException is thrown.

◆ solveSubmatrix()

Vector Kaskade::QPBoundSolverBase< double , QPSlackSolver< QPStructure::DENSE, double > >::solveSubmatrix ( std::vector< int > const &  idx,
Vector const &  b 
) const
protectedinherited

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).

Parameters
idxthe indices defining the submatrix
bthe right hand side (of the same length as the index vector)
Returns
the solution \( x \)

◆ updateWithColumn()

virtual void Kaskade::QPBoundSolverBase< double , QPSlackSolver< QPStructure::DENSE, double > >::updateWithColumn ( double  a,
int  i,
Vector s 
) const
protectedvirtualinherited

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.


The documentation for this class was generated from the following file: