Eigen::QLD Class Reference

A wrapper of the ql algorithm by Professor Schittkowski, with some convention changes on the way the constraints are written. More...

#include <eigen-qld/QLD.h>

Inheritance diagram for Eigen::QLD:
Collaboration diagram for Eigen::QLD:

Public Member Functions

EIGEN_QLD_API QLD ()
 
EIGEN_QLD_API QLD (int nrvar, int nreq, int nrineq, int ldq=-1, bool verbose=false)
 
EIGEN_QLD_API void problem (int nrvar, int nreq, int nrineq, int ldq=-1)
 
EIGEN_QLD_API const VectorXd & multipliers () const
 
template<typename MatObj , typename VecObj , typename MatEq , typename VecEq , typename MatIneq , typename VecIneq , typename VecVar >
bool solve (const MatrixBase< MatObj > &Q, const MatrixBase< VecObj > &c, const MatrixBase< MatEq > &Aeq, const MatrixBase< VecEq > &beq, const MatrixBase< MatIneq > &Aineq, const MatrixBase< VecIneq > &bineq, const MatrixBase< VecVar > &xl, const MatrixBase< VecVar > &xu, bool isDecomp=false, double eps=1e-12)
 
- Public Member Functions inherited from Eigen::QLDDirect
EIGEN_QLD_API QLDDirect ()
 
EIGEN_QLD_API QLDDirect (int nrvar, int nreq, int nrineq, int ldq=-1, int lda=-1, bool verbose=false)
 
EIGEN_QLD_API void fdOut (int fd)
 
EIGEN_QLD_API int fdOut () const
 
EIGEN_QLD_API void verbose (bool v)
 
EIGEN_QLD_API bool verbose () const
 
EIGEN_QLD_API int fail () const
 
EIGEN_QLD_API void problem (int nrvar, int nreq, int nrineq, int ldq=-1, int lda=-1)
 
EIGEN_QLD_API const VectorXd & result () const
 
EIGEN_QLD_API const VectorXd & multipliers () const
 
template<typename MatObj , typename VecObj , typename MatConstr , typename VecConstr , typename VecVar >
bool solve (const MatrixBase< MatObj > &Q, const MatrixBase< VecObj > &c, const MatrixBase< MatConstr > &A, const MatrixBase< VecConstr > &b, const MatrixBase< VecVar > &xl, const MatrixBase< VecVar > &xu, int nreq, bool isDecomp=false, double eps=1e-12)
 

Detailed Description

A wrapper of the ql algorithm by Professor Schittkowski, with some convention changes on the way the constraints are written.

Constructor & Destructor Documentation

◆ QLD() [1/2]

EIGEN_QLD_API Eigen::QLD::QLD ( )

Create a solver without allocating memory. A call to QLD::problem will be necessary before calling QLD::solve

◆ QLD() [2/2]

EIGEN_QLD_API Eigen::QLD::QLD ( int  nrvar,
int  nreq,
int  nrineq,
int  ldq = -1,
bool  verbose = false 
)

Create a solver and allocate the memory necessary to solve a problem with the given dimensions. See also QLD::problem

Member Function Documentation

◆ multipliers()

EIGEN_QLD_API const VectorXd& Eigen::QLD::multipliers ( ) const

Return the lagrange multipliers associated with results, with the multipliers corresponding to equality constraints first, then to inequality constraints, then lower bounds, then upper bounds.

Warning
These are the lagrange multipliers as computed by the underlying ql solver, for which the constraints have the form \( A_{eq} + b_{eq} = 0\) and \( A_{ineq} + b_{ineq} >= 0\).
For now, QLD::solve is changing automatically the constraints, but not the multipliers. After a call to QLD::solveNoOverhead, the multipliers have the same meaning as for ql.

◆ problem()

EIGEN_QLD_API void Eigen::QLD::problem ( int  nrvar,
int  nreq,
int  nrineq,
int  ldq = -1 
)

Allocate the memory necessary to solve a problem with the given dimensions.

Parameters
nrvarSize of the variable vector \(x\).
nreqNumber of equality constraints, i.e. the row size of \(A_{eq}\).
nrineqNumber of inequality constraints, i.e. the row size of \(A_{ineq}\).
ldqLeading dimension of the matrix \(Q\) (see QLDDirect::problem). If smaller than nrvar, nrvar will be used instead (which is the default case).

◆ solve()

template<typename MatObj , typename VecObj , typename MatEq , typename VecEq , typename MatIneq , typename VecIneq , typename VecVar >
bool Eigen::QLD::solve ( const MatrixBase< MatObj > &  Q,
const MatrixBase< VecObj > &  c,
const MatrixBase< MatEq > &  Aeq,
const MatrixBase< VecEq > &  beq,
const MatrixBase< MatIneq > &  Aineq,
const MatrixBase< VecIneq > &  bineq,
const MatrixBase< VecVar > &  xl,
const MatrixBase< VecVar > &  xu,
bool  isDecomp = false,
double  eps = 1e-12 
)
inline

Solve the problem

\begin{align} \underset{{x} \in \mathbb{R}^n}{\text{minimize}} & \ \frac{1}{2}{x^TQx} + {c^Tx} \nonumber \\ \text{subject to} & \ A_{eq} x = b_{eq} \\ & \ A_{ineq} x \leq b_{ineq} \\ & \ x_l \leq x \leq x_u \end{align}

For a positive definite matrix \(Q\), we have the Cholesky decomposition \(Q = R^T R\) with \(R\) upper triangular. If isDecomp is true, \(R\) should be given instead of \(Q\). Only the upper triangular part of \(R\) will be used.

Parameters
QThe matrix \(Q\), or its Cholesky factor \(R\) if isDecomp is true. \(Q\) should be symmetric and positive definite, \(R\) should be upper triangular.
cThe vector \(c\).
AeqThe matrix \(A_{eq}\).
beqThe vector \(b_{eq}\).
AineqThe matrix \(A_{ineq}\).
bineqThe vector \(b_{ineq}\).
xlThe vector \(x_{l}\).
xuThe vector \(x_{u}\).
isDecompspecify if the Cholesky decomposition of \(Q\) is used or not.
epsDesired final accuracy.

This is a wrapper around QLDDirect::solver.


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