Template Function proxsuite::proxqp::dense::dense_backend_choice
Defined in File wrapper.hpp
Function Documentation
-
template<typename T>
DenseBackend proxsuite::proxqp::dense::dense_backend_choice(DenseBackend _dense_backend, isize dim, isize n_eq, isize n_in, bool box_constraints) This class defines the API of PROXQP solver with dense 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 <proxsuite/proxqp/dense/dense.hpp> #include <proxsuite/linalg/veg/util/dbg.hpp> #include <util.hpp> using T = double; 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); Qp<T> qp{ random_with_dim_and_neq_and_n_in, dim, n_eq, n_in, sparsity_factor, strong_convexity_factor}; // Solve the problem proxqp::dense::QP<T> Qp{dim, n_eq, n_in}; // creating QP object Qp.settings.eps_abs = eps_abs; // choose accuracy needed Qp.init(qp.H, qp.g, qp.A, qp.b, qp.C, qp.u, qp.l); // setup the QP object Qp.solve(); // solve the problem // 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; }