Template Struct QP

Struct Documentation

template<typename T>
struct QP

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;
}

Public Functions

inline QP(isize _dim, isize _n_eq, isize _n_in)

Default constructor using QP model dimensions.

Parameters:
  • _dim – primal variable dimension.

  • _n_eq – number of equality constraints.

  • _n_in – number of inequality constraints.

inline void init(optional<MatRef<T>> H, optional<VecRef<T>> g, optional<MatRef<T>> A, optional<VecRef<T>> b, optional<MatRef<T>> 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)

Setups the QP model (with dense matrix format) and equilibrates it if specified by the user.

Parameters:
  • H – quadratic cost input defining the QP model.

  • g – linear cost input defining the QP model.

  • A – equality constraint matrix input defining the QP model.

  • b – equality constraint vector input defining the QP model.

  • C – inequality constraint matrix input defining the QP model.

  • l – lower inequality constraint vector input defining the QP model.

  • u – upper inequality constraint vector input defining the QP model.

  • compute_preconditioner – boolean parameter for executing or not the preconditioner.

  • rho – proximal step size wrt primal variable.

  • mu_eq – proximal step size wrt equality constrained multiplier.

  • mu_in – proximal step size wrt inequality constrained multiplier.

inline void update(optional<MatRef<T>> H, optional<VecRef<T>> g, optional<MatRef<T>> A, optional<VecRef<T>> b, optional<MatRef<T>> C, optional<VecRef<T>> l, optional<VecRef<T>> u, bool update_preconditioner = true, optional<T> rho = nullopt, optional<T> mu_eq = nullopt, optional<T> mu_in = nullopt)

Updates the QP model (with dense matrix format) and re-equilibrates it if specified by the user.

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.

Parameters:
  • H – quadratic cost input defining the QP model.

  • g – linear cost input defining the QP model.

  • A – equality constraint matrix input defining the QP model.

  • b – equality constraint vector input defining the QP model.

  • C – inequality constraint matrix input defining the QP model.

  • l – lower inequality constraint vector input defining the QP model.

  • u – upper inequality constraint vector input defining the QP model.

  • update_preconditioner – bool parameter for updating or not the preconditioner and the associated scaled model.

  • rho – proximal step size wrt primal variable.

  • mu_eq – proximal step size wrt equality constrained multiplier.

  • mu_in – proximal step size wrt inequality constrained multiplier.

inline void solve()

Solves the QP problem using PRXOQP algorithm.

inline void solve(optional<VecRef<T>> x, optional<VecRef<T>> y, optional<VecRef<T>> z)

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

Parameters:
  • x – primal warm start.

  • y – dual equality warm start.

  • z – dual inequality warm start.

inline void cleanup()

Clean-ups solver’s results and workspace.

Public Members

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