.. _program_listing_file__tmp_ws_src_proxsuite_include_proxsuite_proxqp_dense_wrapper.hpp: Program Listing for File wrapper.hpp ==================================== |exhale_lsh| :ref:`Return to documentation for file ` (``/tmp/ws/src/proxsuite/include/proxsuite/proxqp/dense/wrapper.hpp``) .. |exhale_lsh| unicode:: U+021B0 .. UPWARDS ARROW WITH TIP LEFTWARDS .. code-block:: cpp // // Copyright (c) 2022 INRIA // #ifndef PROXSUITE_PROXQP_DENSE_WRAPPER_HPP #define PROXSUITE_PROXQP_DENSE_WRAPPER_HPP #include #include #include #include namespace proxsuite { namespace proxqp { namespace dense { template struct QP { Results results; Settings settings; Model model; Workspace work; preconditioner::RuizEquilibration ruiz; QP(isize _dim, isize _n_eq, isize _n_in) : results(_dim, _n_eq, _n_in) , settings() , model(_dim, _n_eq, _n_in) , work(_dim, _n_eq, _n_in) , ruiz(preconditioner::RuizEquilibration{ _dim, _n_eq + _n_in }) { work.timer.stop(); } void init(optional> H, optional> g, optional> A, optional> b, optional> C, optional> l, optional> u, bool compute_preconditioner = true, optional rho = nullopt, optional mu_eq = nullopt, optional mu_in = nullopt) { // dense case if (settings.compute_timings) { work.timer.stop(); work.timer.start(); } // check the model is valid if (g != nullopt && g.value().size() != 0) { PROXSUITE_CHECK_ARGUMENT_SIZE( g.value().size(), model.dim, "the dimension wrt the primal variable x variable for initializing g " "is not valid."); } else { g.reset(); } if (b != nullopt && b.value().size() != 0) { PROXSUITE_CHECK_ARGUMENT_SIZE( b.value().size(), model.n_eq, "the dimension wrt equality constrained variables for initializing b " "is not valid."); } else { b.reset(); } if (u != nullopt && u.value().size() != 0) { PROXSUITE_CHECK_ARGUMENT_SIZE( u.value().size(), model.n_in, "the dimension wrt inequality constrained variables for initializing u " "is not valid."); } else { u.reset(); } if (l != nullopt && l.value().size() != 0) { PROXSUITE_CHECK_ARGUMENT_SIZE( l.value().size(), model.n_in, "the dimension wrt inequality constrained variables for initializing l " "is not valid."); } else { l.reset(); } if (H != nullopt && H.value().size() != 0) { PROXSUITE_CHECK_ARGUMENT_SIZE( H.value().rows(), model.dim, "the row dimension for initializing H is not valid."); PROXSUITE_CHECK_ARGUMENT_SIZE( H.value().cols(), model.dim, "the column dimension for initializing H is not valid."); } else { H.reset(); } if (A != nullopt && A.value().size() != 0) { PROXSUITE_CHECK_ARGUMENT_SIZE( A.value().rows(), model.n_eq, "the row dimension for initializing A is not valid."); PROXSUITE_CHECK_ARGUMENT_SIZE( A.value().cols(), model.dim, "the column dimension for initializing A is not valid."); } else { A.reset(); } if (C != nullopt && C.value().size() != 0) { PROXSUITE_CHECK_ARGUMENT_SIZE( C.value().rows(), model.n_in, "the row dimension for initializing C is not valid."); PROXSUITE_CHECK_ARGUMENT_SIZE( C.value().cols(), model.dim, "the column dimension for initializing C is not valid."); } else { C.reset(); } if (settings.initial_guess == InitialGuessStatus::WARM_START_WITH_PREVIOUS_RESULT) { work.refactorize = true; // necessary for the first solve (then refactorize only if there // is an update of the matrices) } else { work.refactorize = false; } work.proximal_parameter_update = false; if (settings.compute_timings) { work.timer.stop(); work.timer.start(); } PreconditionerStatus preconditioner_status; if (compute_preconditioner) { preconditioner_status = proxsuite::proxqp::PreconditionerStatus::EXECUTE; } else { preconditioner_status = proxsuite::proxqp::PreconditionerStatus::IDENTITY; } proxsuite::proxqp::dense::update_proximal_parameters( settings, results, work, rho, mu_eq, mu_in); proxsuite::proxqp::dense::setup(H, g, A, b, C, l, u, settings, model, work, results, ruiz, preconditioner_status); work.is_initialized = true; if (settings.compute_timings) { results.info.setup_time = work.timer.elapsed().user; // in microseconds } }; void update(optional> H, optional> g, optional> A, optional> b, optional> C, optional> l, optional> u, bool update_preconditioner = true, optional rho = nullopt, optional mu_eq = nullopt, optional mu_in = nullopt) { if (!work.is_initialized) { init(H, g, A, b, C, l, u, update_preconditioner, rho, mu_eq, mu_in); return; } // dense case work.refactorize = false; work.proximal_parameter_update = false; if (settings.compute_timings) { work.timer.stop(); work.timer.start(); } PreconditionerStatus preconditioner_status; if (update_preconditioner) { preconditioner_status = proxsuite::proxqp::PreconditionerStatus::EXECUTE; } else { preconditioner_status = proxsuite::proxqp::PreconditionerStatus::KEEP; } const bool matrix_update = !(H == nullopt && g == nullopt && A == nullopt && b == nullopt && C == nullopt && u == nullopt && l == nullopt); if (matrix_update) { proxsuite::proxqp::dense::update(H, g, A, b, C, l, u, model, work); } proxsuite::proxqp::dense::update_proximal_parameters( settings, results, work, rho, mu_eq, mu_in); typedef optional> optional_MatRef; typedef optional> optional_VecRef; proxsuite::proxqp::dense::setup(/* avoid double assignation */ optional_MatRef(nullopt), optional_VecRef(nullopt), optional_MatRef(nullopt), optional_VecRef(nullopt), optional_MatRef(nullopt), optional_VecRef(nullopt), optional_VecRef(nullopt), settings, model, work, results, ruiz, preconditioner_status); if (settings.compute_timings) { results.info.setup_time = work.timer.elapsed().user; // in microseconds } }; void solve() { qp_solve( // settings, model, results, work, ruiz); }; void solve(optional> x, optional> y, optional> z) { proxsuite::proxqp::dense::warm_start(x, y, z, results, settings, model); qp_solve( // settings, model, results, work, ruiz); }; void cleanup() { results.cleanup(settings); work.cleanup(); } }; template proxqp::Results solve( optional> H, optional> g, optional> A, optional> b, optional> C, optional> l, optional> u, optional> x = nullopt, optional> y = nullopt, optional> z = nullopt, optional eps_abs = nullopt, optional eps_rel = nullopt, optional rho = nullopt, optional mu_eq = nullopt, optional mu_in = nullopt, optional verbose = nullopt, bool compute_preconditioner = true, bool compute_timings = false, optional max_iter = nullopt, proxsuite::proxqp::InitialGuessStatus initial_guess = proxsuite::proxqp::InitialGuessStatus::EQUALITY_CONSTRAINED_INITIAL_GUESS, bool check_duality_gap = false, optional eps_duality_gap_abs = nullopt, optional eps_duality_gap_rel = nullopt) { isize n(0); isize n_eq(0); isize n_in(0); if (H != nullopt) { n = H.value().rows(); } if (A != nullopt) { n_eq = A.value().rows(); } if (C != nullopt) { n_in = C.value().rows(); } QP Qp(n, n_eq, n_in); Qp.settings.initial_guess = initial_guess; Qp.settings.check_duality_gap = check_duality_gap; if (eps_abs != nullopt) { Qp.settings.eps_abs = eps_abs.value(); } if (eps_rel != nullopt) { Qp.settings.eps_rel = eps_rel.value(); } if (verbose != nullopt) { Qp.settings.verbose = verbose.value(); } if (max_iter != nullopt) { Qp.settings.max_iter = max_iter.value(); } if (eps_duality_gap_abs != nullopt) { Qp.settings.eps_duality_gap_abs = eps_duality_gap_abs.value(); } if (eps_duality_gap_rel != nullopt) { Qp.settings.eps_duality_gap_rel = eps_duality_gap_rel.value(); } Qp.settings.compute_timings = compute_timings; Qp.init(H, g, A, b, C, l, u, compute_preconditioner, rho, mu_eq, mu_in); Qp.solve(x, y, z); return Qp.results; } template bool operator==(const QP& qp1, const QP& qp2) { bool value = qp1.model == qp2.model && qp1.settings == qp2.settings && qp1.results == qp2.results; return value; } template bool operator!=(const QP& qp1, const QP& qp2) { return !(qp1 == qp2); } template struct BatchQP { std::vector> vector_qp; dense::isize m_size; explicit BatchQP(size_t batch_size) { if (vector_qp.max_size() != batch_size) { vector_qp.clear(); vector_qp.reserve(batch_size); } m_size = 0; } QP& init_qp_in_place(dense::isize dim, dense::isize n_eq, dense::isize n_in) { vector_qp.emplace_back(dim, n_eq, n_in); auto& qp = vector_qp.back(); m_size++; return qp; }; void insert(QP& qp) { vector_qp.emplace_back(qp); }; QP& get(isize i) { return vector_qp.at(i); }; QP& operator[](isize i) { return vector_qp.at(i); }; dense::isize size() { return m_size; }; }; } // namespace dense } // namespace proxqp } // namespace proxsuite #endif /* end of include guard PROXSUITE_PROXQP_DENSE_WRAPPER_HPP */