proxsuite 0.6.7
The Advanced Proximal Optimization Toolbox
Loading...
Searching...
No Matches
proxsuite::proxqp::sparse::QP< T, I > Struct Template Reference

This class defines the API of PROXQP solver with sparse backend. More...

#include <proxsuite/proxqp/sparse/wrapper.hpp>

Public Member Functions

 QP (isize dim, isize n_eq, isize n_in)
 
 QP (const SparseMat< bool, I > &H, const SparseMat< bool, I > &A, const SparseMat< bool, I > &C)
 
void init (optional< SparseMat< T, I > > H, optional< VecRef< T > > g, optional< SparseMat< T, I > > A, optional< VecRef< T > > b, optional< SparseMat< T, I > > C, optional< VecRef< T > > l, optional< VecRef< T > > u, bool compute_preconditioner_=true, optional< T > rho=nullopt, optional< T > mu_eq=nullopt, optional< T > mu_in=nullopt, optional< T > manual_minimal_H_eigenvalue=nullopt)
 
void update (const optional< SparseMat< T, I > > H, optional< VecRef< T > > g, const optional< SparseMat< T, I > > A, optional< VecRef< T > > b, const optional< SparseMat< T, I > > C, optional< VecRef< T > > l, optional< VecRef< T > > u, bool update_preconditioner=false, optional< T > rho=nullopt, optional< T > mu_eq=nullopt, optional< T > mu_in=nullopt, optional< T > manual_minimal_H_eigenvalue=nullopt)
 
void solve ()
 
void solve (optional< VecRef< T > > x, optional< VecRef< T > > y, optional< VecRef< T > > z)
 
void cleanup ()
 

Public Attributes

Results< T > results
 
Settings< T > settings
 
Model< T, I > model
 
Workspace< T, I > work
 
preconditioner::RuizEquilibration< T, I > ruiz
 

Detailed Description

template<typename T, typename I>
struct proxsuite::proxqp::sparse::QP< T, I >

This class defines the API of PROXQP solver with sparse backend.

Wrapper class for using proxsuite API with dense backend for solving linearly constrained convex QP problem using ProxQp algorithm.

Example usage:

#include <Eigen/Core>
#include <Eigen/Cholesky>
#include <util.hpp>
using T = double;
using I = c_int;
auto main() -> int {
// Generate a random QP problem with primal variable dimension of size
dim; n_eq equality constraints and n_in inequality constraints
::proxsuite::proxqp::test::rand::set_seed(1);
proxqp::isize dim = 10;
proxqp::isize n_eq(dim / 4);
proxqp::isize n_in(dim / 4);
T strong_convexity_factor(1.e-2);
T sparsity_factor = 0.15; // controls the sparsity of each matrix of the
problem generated T eps_abs = T(1e-9); double p = 1.0; T conditioning(10.0);
auto H =
::proxsuite::proxqp::test::rand::sparse_positive_definite_rand(n, conditioning,
p); auto g = ::proxsuite::proxqp::test::rand::vector_rand<T>(n); auto A =
::proxsuite::proxqp::test::rand::sparse_matrix_rand<T>(n_eq,n, p); auto b =
::proxsuite::proxqp::test::rand::vector_rand<T>(n_eq); auto C =
::proxsuite::proxqp::test::rand::sparse_matrix_rand<T>(n_in,n, p); auto l =
::proxsuite::proxqp::test::rand::vector_rand<T>(n_in); auto u = (l.array() +
1).matrix().eval();
proxqp::sparse::QP<T,I> Qp(n, n_eq, n_in);
Qp.settings.eps_abs = 1.E-9;
Qp.settings.verbose = true;
Qp.setup_sparse_matrices(H,g,A,b,C,u,l);
Qp.solve();
// Solve the problem
proxqp::sparse::QP<T,I> Qp(n, n_eq, n_in);
Qp.settings.eps_abs = 1.E-9;
Qp.settings.verbose = true;
Qp.setup_sparse_matrices(H,g,A,b,C,u,l);
Qp.solve();
// Verify solution accuracy
T pri_res = std::max(
(qp.A * Qp.results.x - qp.b).lpNorm<Eigen::Infinity>(),
(helpers::positive_part(qp.C * Qp.results.x -
qp.u) + helpers::negative_part(qp.C * Qp.results.x - qp.l))
.lpNorm<Eigen::Infinity>());
T dua_res = (qp.H * Qp.results.x + qp.g + qp.A.transpose() *
Qp.results.y + qp.C.transpose() * Qp.results.z) .lpNorm<Eigen::Infinity>();
VEG_ASSERT(pri_res <= eps_abs);
VEG_ASSERT(dua_res <= eps_abs);
// Some solver statistics
std::cout << "------solving qp with dim: " << dim
<< " neq: " << n_eq << " nin: "
<< n_in << std::endl; std::cout << "primal residual: " << pri_res << std::endl;
std::cout << "dual residual: " << dua_res << std::endl;
std::cout << "total number of iteration: " << Qp.results.info.iter
<< std::endl;
}
#define VEG_ASSERT(...)
auto positive_part(T const &expr) PROXSUITE_DEDUCE_RET((expr.array() > 0).select(expr
Returns the positive part of an expression.
auto negative_part(T const &expr) PROXSUITE_DEDUCE_RET((expr.array()< 0).select(expr
Returns the negative part of an expression.
This class defines the API of PROXQP solver with sparse backend.
Definition wrapper.hpp:91

Definition at line 90 of file wrapper.hpp.

Constructor & Destructor Documentation

◆ QP() [1/2]

template<typename T , typename I >
proxsuite::proxqp::sparse::QP< T, I >::QP ( isize dim,
isize n_eq,
isize n_in )
inline

Default constructor using the dimension of the matrices in entry.

Parameters
dimprimal variable dimension.
n_eqnumber of equality constraints.
n_innumber of inequality constraints.

Definition at line 103 of file wrapper.hpp.

◆ QP() [2/2]

template<typename T , typename I >
proxsuite::proxqp::sparse::QP< T, I >::QP ( const SparseMat< bool, I > & H,
const SparseMat< bool, I > & A,
const SparseMat< bool, I > & C )
inline

Default constructor using the sparsity structure of the matrices in entry.

Parameters
Hboolean mask of the quadratic cost input defining the QP model.
Aboolean mask of the equality constraint matrix input defining the QP model.
Cboolean mask of the inequality constraint matrix input defining the QP model.

Definition at line 122 of file wrapper.hpp.

Member Function Documentation

◆ init()

template<typename T , typename I >
void proxsuite::proxqp::sparse::QP< T, I >::init ( optional< SparseMat< T, I > > H,
optional< VecRef< T > > g,
optional< SparseMat< T, I > > A,
optional< VecRef< T > > b,
optional< SparseMat< T, I > > C,
optional< VecRef< T > > l,
optional< VecRef< T > > u,
bool compute_preconditioner_ = true,
optional< T > rho = nullopt,
optional< T > mu_eq = nullopt,
optional< T > mu_in = nullopt,
optional< T > manual_minimal_H_eigenvalue = nullopt )
inline

Setups the QP model (with sparse matrix format) and equilibrates it.

Parameters
Hquadratic cost input defining the QP model.
glinear cost input defining the QP model.
Aequality constraint matrix input defining the QP model.
bequality constraint vector input defining the QP model.
Cinequality constraint matrix input defining the QP model.
llower inequality constraint vector input defining the QP model.
uupper inequality constraint vector input defining the QP model.
compute_preconditionerboolean parameter for executing or not the preconditioner.
rhoproximal step size wrt primal variable.
mu_eqproximal step size wrt equality constrained multiplier.
mu_inproximal step size wrt inequality constrained multiplier.

Definition at line 165 of file wrapper.hpp.

◆ update()

template<typename T , typename I >
void proxsuite::proxqp::sparse::QP< T, I >::update ( const optional< SparseMat< T, I > > H,
optional< VecRef< T > > g,
const optional< SparseMat< T, I > > A,
optional< VecRef< T > > b,
const optional< SparseMat< T, I > > C,
optional< VecRef< T > > l,
optional< VecRef< T > > u,
bool update_preconditioner = false,
optional< T > rho = nullopt,
optional< T > mu_eq = nullopt,
optional< T > mu_in = nullopt,
optional< T > manual_minimal_H_eigenvalue = nullopt )
inline

Updates the QP model (with sparse matrix format) and re-equilibrates it if specified by the user. If matrices in entry are not null, the update is effective only if the sparsity structure of entry is the same as the one used for the initialization.

Parameters
Hquadratic cost input defining the QP model.
glinear cost input defining the QP model.
Aequality constraint matrix input defining the QP model.
bequality constraint vector input defining the QP model.
Cinequality constraint matrix input defining the QP model.
llower inequality constraint vector input defining the QP model.
ulower inequality constraint vector input defining the QP model.
update_preconditionerbool parameter for updating or not the preconditioner and the associated scaled model.
rhoproximal step size wrt primal variable.
mu_eqproximal step size wrt equality constrained multiplier.
mu_inproximal step size wrt inequality constrained multiplier.
Note
The init method should be called before update. If it has not been done before, init is called depending on the is_initialized flag.

Definition at line 352 of file wrapper.hpp.

◆ solve() [1/2]

template<typename T , typename I >
void proxsuite::proxqp::sparse::QP< T, I >::solve ( )
inline

Solves the QP problem using PRXOQP algorithm.

Definition at line 637 of file wrapper.hpp.

◆ solve() [2/2]

template<typename T , typename I >
void proxsuite::proxqp::sparse::QP< T, I >::solve ( optional< VecRef< T > > x,
optional< VecRef< T > > y,
optional< VecRef< T > > z )
inline

Solves the QP problem using PROXQP algorithm and a warm start.

Parameters
xprimal warm start.
ydual equality warm start.
zdual inequality warm start.

Definition at line 652 of file wrapper.hpp.

◆ cleanup()

template<typename T , typename I >
void proxsuite::proxqp::sparse::QP< T, I >::cleanup ( )
inline

Clean-ups solver's results.

Definition at line 667 of file wrapper.hpp.

Member Data Documentation

◆ results

template<typename T , typename I >
Results<T> proxsuite::proxqp::sparse::QP< T, I >::results

Definition at line 92 of file wrapper.hpp.

◆ settings

template<typename T , typename I >
Settings<T> proxsuite::proxqp::sparse::QP< T, I >::settings

Definition at line 93 of file wrapper.hpp.

◆ model

template<typename T , typename I >
Model<T, I> proxsuite::proxqp::sparse::QP< T, I >::model

Definition at line 94 of file wrapper.hpp.

◆ work

template<typename T , typename I >
Workspace<T, I> proxsuite::proxqp::sparse::QP< T, I >::work

Definition at line 95 of file wrapper.hpp.

◆ ruiz

template<typename T , typename I >
preconditioner::RuizEquilibration<T, I> proxsuite::proxqp::sparse::QP< T, I >::ruiz

Definition at line 96 of file wrapper.hpp.


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