8#ifndef PROXSUITE_PROXQP_SPARSE_WRAPPER_HPP
9#define PROXSUITE_PROXQP_SPARSE_WRAPPER_HPP
89template<
typename T,
typename I>
106 ,
model(dim, n_eq, n_in)
108 ,
ruiz(dim, n_eq + n_in, 1e-3, 10, preconditioner::
Symmetry::UPPER)
111 work.internal.do_symbolic_fact =
true;
112 work.internal.is_initialized =
false;
125 :
QP(H.rows(), A.rows(), C.rows())
135 proxsuite::linalg::sparse::from_eigen, H_triu
138 proxsuite::linalg::sparse::from_eigen, AT
141 proxsuite::linalg::sparse::from_eigen, CT
143 work.setup_symbolic_factorizaton(
144 model, Href.symbolic(), ATref.symbolic(), CTref.symbolic());
146 results.info.setup_time = work.timer.elapsed().user;
177 bool compute_preconditioner_ =
true,
187 if (g !=
nullopt && g.value().size() != 0) {
191 "the dimension wrt the primal variable x variable for initializing g "
196 if (b !=
nullopt && b.value().size() != 0) {
200 "the dimension wrt equality constrained variables for initializing b "
205 if (u !=
nullopt && u.value().size() != 0) {
209 "the dimension wrt inequality constrained variables for initializing u "
214 if (l !=
nullopt && l.value().size() != 0) {
218 "the dimension wrt inequality constrained variables for initializing l "
223 if (H !=
nullopt && H.value().size() != 0) {
227 "the row dimension for initializing H is not valid.");
231 "the column dimension for initializing H is not valid.");
235 if (A !=
nullopt && A.value().size() != 0) {
239 "the row dimension for initializing A is not valid.");
243 "the column dimension for initializing A is not valid.");
247 if (C !=
nullopt && C.value().size() != 0) {
251 "the row dimension for initializing C is not valid.");
255 "the column dimension for initializing C is not valid.");
259 work.internal.proximal_parameter_update =
false;
261 if (compute_preconditioner_) {
289 AT = (A.value()).transpose();
295 CT = (C.value()).transpose();
301 (H.value()).template triangularView<Eigen::Upper>();
304 { proxsuite::linalg::sparse::from_eigen, H_triu },
305 { proxsuite::linalg::sparse::from_eigen,
model.g },
306 { proxsuite::linalg::sparse::from_eigen, AT },
307 { proxsuite::linalg::sparse::from_eigen,
model.b },
308 { proxsuite::linalg::sparse::from_eigen, CT },
309 { proxsuite::linalg::sparse::from_eigen,
model.l },
310 { proxsuite::linalg::sparse::from_eigen,
model.u }
319 H_triu = (H.value()).template triangularView<Eigen::Upper>();
321 { proxsuite::linalg::sparse::from_eigen, H_triu },
322 { proxsuite::linalg::sparse::from_eigen,
model.g },
323 { proxsuite::linalg::sparse::from_eigen, AT },
324 { proxsuite::linalg::sparse::from_eigen,
model.b },
325 { proxsuite::linalg::sparse::from_eigen, CT },
326 { proxsuite::linalg::sparse::from_eigen,
model.l },
327 { proxsuite::linalg::sparse::from_eigen,
model.u }
331 work.internal.is_initialized =
true;
334 results.info.setup_time +=
work.timer.elapsed().user;
364 bool update_preconditioner =
false,
370 if (!
work.internal.is_initialized) {
371 init(H, g, A, b, C, l, u, update_preconditioner, rho, mu_eq, mu_in);
378 work.internal.dirty =
false;
379 work.internal.proximal_parameter_update =
false;
381 if (update_preconditioner) {
390 model.kkt_mut_unscaled();
393 proxsuite::linalg::veg::unsafe, kkt_unscaled, n);
408 "the dimension wrt the primal variable x "
409 "variable for updating g is not valid.");
414 "the dimension wrt equality constrained "
415 "variables for updating b is not valid.");
420 "the dimension wrt inequality constrained "
421 "variables for updating u is not valid.");
426 "the dimension wrt inequality constrained "
427 "variables for updating l is not valid.");
433 "the row dimension for updating H is not valid.");
437 "the column dimension for updating H is not valid.");
443 "the row dimension for updating A is not valid.");
447 "the column dimension for updating A is not valid.");
453 "the row dimension for updating C is not valid.");
457 "the column dimension for updating C is not valid.");
476 H.value().template triangularView<Eigen::Upper>();
482 { proxsuite::linalg::sparse::from_eigen, H_triu }) &&
484 { proxsuite::linalg::sparse::from_eigen,
485 SparseMat<T, I>(A.value().transpose()) }) &&
487 { proxsuite::linalg::sparse::from_eigen,
488 SparseMat<T, I>(C.value().transpose()) });
494 { proxsuite::linalg::sparse::from_eigen,
498 { proxsuite::linalg::sparse::from_eigen,
502 { proxsuite::linalg::sparse::from_eigen,
509 { proxsuite::linalg::sparse::from_eigen, H_triu }) &&
511 { proxsuite::linalg::sparse::from_eigen,
512 SparseMat<T, I>(A.value().transpose()) });
518 { proxsuite::linalg::sparse::from_eigen,
522 { proxsuite::linalg::sparse::from_eigen,
530 { proxsuite::linalg::sparse::from_eigen, H_triu }) &&
532 { proxsuite::linalg::sparse::from_eigen,
533 SparseMat<T, I>(C.value().transpose()) });
539 { proxsuite::linalg::sparse::from_eigen,
542 { proxsuite::linalg::sparse::from_eigen,
549 { proxsuite::linalg::sparse::from_eigen, H_triu });
555 { proxsuite::linalg::sparse::from_eigen,
563 { proxsuite::linalg::sparse::from_eigen,
564 SparseMat<T, I>(A.value().transpose()) }) &&
566 { proxsuite::linalg::sparse::from_eigen,
567 SparseMat<T, I>(C.value().transpose()) });
573 { proxsuite::linalg::sparse::from_eigen,
576 { proxsuite::linalg::sparse::from_eigen,
582 { proxsuite::linalg::sparse::from_eigen,
583 SparseMat<T, I>(A.value().transpose()) });
589 { proxsuite::linalg::sparse::from_eigen,
596 { proxsuite::linalg::sparse::from_eigen,
597 SparseMat<T, I>(C.value().transpose()) });
603 { proxsuite::linalg::sparse::from_eigen,
609 H_unscaled.
to_eigen().template triangularView<Eigen::Upper>();
611 { proxsuite::linalg::sparse::from_eigen, H_triu },
612 { proxsuite::linalg::sparse::from_eigen,
model.g },
613 { proxsuite::linalg::sparse::from_eigen, AT_unscaled.
to_eigen() },
614 { proxsuite::linalg::sparse::from_eigen,
model.b },
615 { proxsuite::linalg::sparse::from_eigen, CT_unscaled.
to_eigen() },
616 { proxsuite::linalg::sparse::from_eigen,
model.l },
617 { proxsuite::linalg::sparse::from_eigen,
model.u }
627 preconditioner_status);
635 results.info.setup_time =
work.timer.elapsed().user;
710template<
typename T,
typename I>
729 bool compute_preconditioner =
true,
730 bool compute_timings =
false,
736 bool check_duality_gap =
false,
739 bool primal_infeasibility_solving =
false,
747 n = H.value().rows();
750 n_eq = A.value().rows();
753 n_in = C.value().rows();
757 Qp.
settings.initial_guess = initial_guess;
758 Qp.
settings.check_duality_gap = check_duality_gap;
761 Qp.
settings.eps_abs = eps_abs.value();
764 Qp.
settings.eps_rel = eps_rel.value();
767 Qp.
settings.verbose = verbose.value();
770 Qp.
settings.max_iter = max_iter.value();
772 if (eps_duality_gap_abs !=
nullopt) {
773 Qp.
settings.eps_duality_gap_abs = eps_duality_gap_abs.value();
775 if (eps_duality_gap_rel !=
nullopt) {
776 Qp.
settings.eps_duality_gap_rel = eps_duality_gap_rel.value();
778 Qp.
settings.compute_timings = compute_timings;
779 Qp.
settings.sparse_backend = sparse_backend;
780 Qp.
settings.primal_infeasibility_solving = primal_infeasibility_solving;
781 if (manual_minimal_H_eigenvalue !=
nullopt) {
789 compute_preconditioner,
793 manual_minimal_H_eigenvalue.value());
796 H, g, A, b, C, l, u, compute_preconditioner, rho, mu_eq, mu_in,
nullopt);
804template<
typename T,
typename I>
#define PROXSUITE_CHECK_ARGUMENT_SIZE(size, expected_size, message)
_detail::_meta::make_signed< usize >::Type isize
auto middle_cols_mut(proxsuite::linalg::sparse::MatMut< T, I > mat, isize start, isize ncols, isize nnz) -> proxsuite::linalg::sparse::MatMut< T, I >
auto top_rows_mut_unchecked(proxsuite::linalg::veg::Unsafe, proxsuite::linalg::sparse::MatMut< T, I > mat, isize nrows) -> proxsuite::linalg::sparse::MatMut< T, I >
void warm_start(optional< VecRef< T > > x_wm, optional< VecRef< T > > y_wm, optional< VecRef< T > > z_wm, Results< T > &results, Settings< T > &settings, Model< T, I > &model)
void update_proximal_parameters(Settings< T > &settings, Results< T > &results, Workspace< T, I > &work, optional< T > rho_new, optional< T > mu_eq_new, optional< T > mu_in_new)
void qp_solve(Results< T > &results, Model< T, I > &data, const Settings< T > &settings, Workspace< T, I > &work, P &precond)
Eigen::SparseMatrix< T, Eigen::ColMajor, I > SparseMat
proxqp::Results< T > solve(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, optional< VecRef< T > > x=nullopt, optional< VecRef< T > > y=nullopt, optional< VecRef< T > > z=nullopt, optional< T > eps_abs=nullopt, optional< T > eps_rel=nullopt, optional< T > rho=nullopt, optional< T > mu_eq=nullopt, optional< T > mu_in=nullopt, optional< bool > verbose=nullopt, bool compute_preconditioner=true, bool compute_timings=false, optional< isize > max_iter=nullopt, proxsuite::proxqp::InitialGuessStatus initial_guess=proxsuite::proxqp::InitialGuessStatus::EQUALITY_CONSTRAINED_INITIAL_GUESS, proxsuite::proxqp::SparseBackend sparse_backend=proxsuite::proxqp::SparseBackend::Automatic, bool check_duality_gap=false, optional< T > eps_duality_gap_abs=nullopt, optional< T > eps_duality_gap_rel=nullopt, bool primal_infeasibility_solving=false, optional< T > manual_minimal_H_eigenvalue=nullopt)
void copy(proxsuite::linalg::sparse::MatMut< T, I > a, proxsuite::linalg::sparse::MatRef< T, I > b)
Eigen::Ref< Eigen::Matrix< T, DYN, 1 > const > VecRef
auto have_same_structure(proxsuite::linalg::sparse::MatRef< T, I > a, proxsuite::linalg::sparse::MatRef< T, I > b) -> bool
void qp_setup(QpView< T, I > qp, Results< T > &results, Model< T, I > &data, Workspace< T, I > &work, Settings< T > &settings, P &precond, PreconditionerStatus &preconditioner_status)
void update_default_rho_with_minimal_Hessian_eigen_value(optional< T > manual_minimal_H_eigenvalue, Results< T > &results, Settings< T > &settings)
@ EQUALITY_CONSTRAINED_INITIAL_GUESS
constexpr nullopt_t nullopt
auto as_const() const noexcept -> MatRef< T, I >
auto to_eigen() const noexcept -> Eigen::Map< Eigen::SparseMatrix< T, Eigen::ColMajor, I > >
This class stores all the results of PROXQP solvers with sparse and dense backends.
This class defines the settings of PROXQP solvers with sparse and dense backends.
std::vector< QP< T, I > > qp_vector
void insert(QP< T, I > &qp)
BatchQP(BatchQP &&)=default
QP< T, I > & operator[](isize i)
BatchQP & operator=(const BatchQP &)=delete
QP< T, I > & init_qp_in_place(sparse::isize dim, sparse::isize n_eq, sparse::isize n_in)
QP< T, I > & get(isize i)
const QP< T, I > & operator[](isize i) const
BatchQP & operator=(BatchQP &&)=default
const QP< T, I > & get(isize i) const
BatchQP(long unsigned int batchSize)
BatchQP(const BatchQP &)=delete
This class stores the model of the QP problem.
This class defines the API of PROXQP solver with sparse backend.
QP & operator=(const QP &)=delete
QP(const SparseMat< bool, I > &H, const SparseMat< bool, I > &A, const SparseMat< bool, I > &C)
void solve(optional< VecRef< T > > x, optional< VecRef< T > > y, optional< VecRef< T > > z)
QP & operator=(QP &&)=default
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 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)
QP(isize dim, isize n_eq, isize n_in)
preconditioner::RuizEquilibration< T, I > ruiz
This class defines the workspace of the sparse solver.