proxsuite 0.7.1
The Advanced Proximal Optimization Toolbox
Loading...
Searching...
No Matches
wrapper.hpp
Go to the documentation of this file.
1//
2// Copyright (c) 2022 INRIA
3//
8#ifndef PROXSUITE_PROXQP_SPARSE_WRAPPER_HPP
9#define PROXSUITE_PROXQP_SPARSE_WRAPPER_HPP
14
15namespace proxsuite {
16namespace proxqp {
17namespace sparse {
21
89template<typename T, typename I>
90struct QP
91{
103 QP(isize dim, isize n_eq, isize n_in)
104 : results(dim, n_eq, n_in)
105 , settings()
106 , model(dim, n_eq, n_in)
107 , work()
108 , ruiz(dim, n_eq + n_in, 1e-3, 10, preconditioner::Symmetry::UPPER)
109 {
110 work.timer.stop();
111 work.internal.do_symbolic_fact = true;
112 work.internal.is_initialized = false;
113 }
123 const SparseMat<bool, I>& A,
124 const SparseMat<bool, I>& C)
125 : QP(H.rows(), A.rows(), C.rows())
126 {
127 if (settings.compute_timings) {
128 work.timer.stop();
129 work.timer.start();
130 }
131 SparseMat<bool, I> H_triu = H.template triangularView<Eigen::Upper>();
132 SparseMat<bool, I> AT = A.transpose();
133 SparseMat<bool, I> CT = C.transpose();
135 proxsuite::linalg::sparse::from_eigen, H_triu
136 };
138 proxsuite::linalg::sparse::from_eigen, AT
139 };
141 proxsuite::linalg::sparse::from_eigen, CT
142 };
143 work.setup_symbolic_factorizaton(
144 model, Href.symbolic(), ATref.symbolic(), CTref.symbolic());
145 if (settings.compute_timings) {
146 results.info.setup_time = work.timer.elapsed().user; // in microseconds
147 }
148 }
149
150 QP(const QP&) = delete;
151 QP& operator=(const QP&) = delete;
152 QP(QP&&) = default;
153 QP& operator=(QP&&) = default;
154
177 bool compute_preconditioner_ = true,
178 optional<T> rho = nullopt,
179 optional<T> mu_eq = nullopt,
180 optional<T> mu_in = nullopt,
181 optional<T> manual_minimal_H_eigenvalue = nullopt)
182 {
183 if (settings.compute_timings) {
184 work.timer.stop();
185 work.timer.start();
186 }
187 if (g != nullopt && g.value().size() != 0) {
189 g.value().size(),
190 model.dim,
191 "the dimension wrt the primal variable x variable for initializing g "
192 "is not valid.");
193 } else {
194 g.reset();
195 }
196 if (b != nullopt && b.value().size() != 0) {
198 b.value().size(),
199 model.n_eq,
200 "the dimension wrt equality constrained variables for initializing b "
201 "is not valid.");
202 } else {
203 b.reset();
204 }
205 if (u != nullopt && u.value().size() != 0) {
207 u.value().size(),
208 model.n_in,
209 "the dimension wrt inequality constrained variables for initializing u "
210 "is not valid.");
211 } else {
212 u.reset();
213 }
214 if (l != nullopt && l.value().size() != 0) {
216 l.value().size(),
217 model.n_in,
218 "the dimension wrt inequality constrained variables for initializing l "
219 "is not valid.");
220 } else {
221 l.reset();
222 }
223 if (H != nullopt && H.value().size() != 0) {
225 H.value().rows(),
226 model.dim,
227 "the row dimension for initializing H is not valid.");
229 H.value().cols(),
230 model.dim,
231 "the column dimension for initializing H is not valid.");
232 } else {
233 H.reset();
234 }
235 if (A != nullopt && A.value().size() != 0) {
237 A.value().rows(),
238 model.n_eq,
239 "the row dimension for initializing A is not valid.");
241 A.value().cols(),
242 model.dim,
243 "the column dimension for initializing A is not valid.");
244 } else {
245 A.reset();
246 }
247 if (C != nullopt && C.value().size() != 0) {
249 C.value().rows(),
250 model.n_in,
251 "the row dimension for initializing C is not valid.");
253 C.value().cols(),
254 model.dim,
255 "the column dimension for initializing C is not valid.");
256 } else {
257 C.reset();
258 }
259 work.internal.proximal_parameter_update = false;
260 PreconditionerStatus preconditioner_status;
261 if (compute_preconditioner_) {
263 } else {
265 }
267 settings, results, work, rho, mu_eq, mu_in);
268
269 if (g != nullopt) {
270 model.g = g.value();
271 } // else qpmodel.g remains initialzed to a matrix with zero elements or
272 // zero shape
273 if (b != nullopt) {
274 model.b = b.value();
275 } // else qpmodel.b remains initialzed to a matrix with zero elements or
276 // zero shape
277 if (u != nullopt) {
278 model.u = u.value();
279 } // else qpmodel.u remains initialzed to a matrix with zero elements or
280 // zero shape
281 if (l != nullopt) {
282 model.l = l.value();
283 } // else qpmodel.l remains initialzed to a matrix with zero elements or
284 // zero shape
285
286 // avoid allocations when H is not nullopt
287 SparseMat<T, I> AT(model.dim, model.n_eq);
288 if (A != nullopt) {
289 AT = (A.value()).transpose();
290 } else {
291 AT.setZero();
292 }
293 SparseMat<T, I> CT(model.dim, model.n_in);
294 if (C != nullopt) {
295 CT = (C.value()).transpose();
296 } else {
297 CT.setZero();
298 }
299 if (H != nullopt) {
300 SparseMat<T, I> H_triu =
301 (H.value()).template triangularView<Eigen::Upper>();
302
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 }
311 };
312 qp_setup(qp, results, model, work, settings, ruiz, preconditioner_status);
315 manual_minimal_H_eigenvalue, results, settings);
316 } else {
317 SparseMat<T, I> H_triu(model.dim, model.dim);
318 H_triu.setZero();
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 }
328 };
329 qp_setup(qp, results, model, work, settings, ruiz, preconditioner_status);
330 }
331 work.internal.is_initialized = true;
332
333 if (settings.compute_timings) {
334 results.info.setup_time += work.timer.elapsed().user; // in microseconds
335 }
336 };
359 const optional<SparseMat<T, I>> A,
361 const optional<SparseMat<T, I>> C,
364 bool update_preconditioner = false,
365 optional<T> rho = nullopt,
366 optional<T> mu_eq = nullopt,
367 optional<T> mu_in = nullopt,
368 optional<T> manual_minimal_H_eigenvalue = nullopt)
369 {
370 if (!work.internal.is_initialized) {
371 init(H, g, A, b, C, l, u, update_preconditioner, rho, mu_eq, mu_in);
372 return;
373 }
374 if (settings.compute_timings) {
375 work.timer.stop();
376 work.timer.start();
377 }
378 work.internal.dirty = false;
379 work.internal.proximal_parameter_update = false;
380 PreconditionerStatus preconditioner_status;
381 if (update_preconditioner) {
383 } else {
384 preconditioner_status = proxsuite::proxqp::PreconditionerStatus::KEEP;
385 }
386 isize n = model.dim;
387 isize n_eq = model.n_eq;
388 isize n_in = model.n_in;
390 model.kkt_mut_unscaled();
391
392 auto kkt_top_n_rows = detail::top_rows_mut_unchecked(
393 proxsuite::linalg::veg::unsafe, kkt_unscaled, n);
394
396 detail::middle_cols_mut(kkt_top_n_rows, 0, n, model.H_nnz);
397
399 detail::middle_cols_mut(kkt_top_n_rows, n, n_eq, model.A_nnz);
400
402 detail::middle_cols_mut(kkt_top_n_rows, n + n_eq, n_in, model.C_nnz);
403
404 // check the model is valid
405 if (g != nullopt) {
406 PROXSUITE_CHECK_ARGUMENT_SIZE(g.value().size(),
407 model.dim,
408 "the dimension wrt the primal variable x "
409 "variable for updating g is not valid.");
410 }
411 if (b != nullopt) {
412 PROXSUITE_CHECK_ARGUMENT_SIZE(b.value().size(),
413 model.n_eq,
414 "the dimension wrt equality constrained "
415 "variables for updating b is not valid.");
416 }
417 if (u != nullopt) {
418 PROXSUITE_CHECK_ARGUMENT_SIZE(u.value().size(),
419 model.n_in,
420 "the dimension wrt inequality constrained "
421 "variables for updating u is not valid.");
422 }
423 if (l != nullopt) {
424 PROXSUITE_CHECK_ARGUMENT_SIZE(l.value().size(),
425 model.n_in,
426 "the dimension wrt inequality constrained "
427 "variables for updating l is not valid.");
428 }
429 if (H != nullopt) {
431 H.value().rows(),
432 model.dim,
433 "the row dimension for updating H is not valid.");
435 H.value().cols(),
436 model.dim,
437 "the column dimension for updating H is not valid.");
438 }
439 if (A != nullopt) {
441 A.value().rows(),
442 model.n_eq,
443 "the row dimension for updating A is not valid.");
445 A.value().cols(),
446 model.dim,
447 "the column dimension for updating A is not valid.");
448 }
449 if (C != nullopt) {
451 C.value().rows(),
452 model.n_in,
453 "the row dimension for updating C is not valid.");
455 C.value().cols(),
456 model.dim,
457 "the column dimension for updating C is not valid.");
458 }
459
460 // update the model
461
462 if (g != nullopt) {
463 model.g = g.value();
464 }
465 if (b != nullopt) {
466 model.b = b.value();
467 }
468 if (u != nullopt) {
469 model.u = u.value();
470 }
471 if (l != nullopt) {
472 model.l = l.value();
473 }
474 if (H != nullopt) {
475 SparseMat<T, I> H_triu =
476 H.value().template triangularView<Eigen::Upper>();
477 if (A != nullopt) {
478 if (C != nullopt) {
479 bool res =
481 H_unscaled.as_const(),
482 { proxsuite::linalg::sparse::from_eigen, H_triu }) &&
483 have_same_structure(AT_unscaled.as_const(),
484 { proxsuite::linalg::sparse::from_eigen,
485 SparseMat<T, I>(A.value().transpose()) }) &&
486 have_same_structure(CT_unscaled.as_const(),
487 { proxsuite::linalg::sparse::from_eigen,
488 SparseMat<T, I>(C.value().transpose()) });
489 /* TO PUT IN DEBUG MODE
490 std::cout << "have same structure = " << res << std::endl;
491 */
492 if (res) {
493 copy(H_unscaled,
494 { proxsuite::linalg::sparse::from_eigen,
495 H_triu }); // copy rhs into lhs
496 copy(
497 AT_unscaled,
498 { proxsuite::linalg::sparse::from_eigen,
499 SparseMat<T, I>(A.value().transpose()) }); // copy rhs into lhs
500 copy(
501 CT_unscaled,
502 { proxsuite::linalg::sparse::from_eigen,
503 SparseMat<T, I>(C.value().transpose()) }); // copy rhs into lhs
504 }
505 } else {
506 bool res =
508 H_unscaled.as_const(),
509 { proxsuite::linalg::sparse::from_eigen, H_triu }) &&
510 have_same_structure(AT_unscaled.as_const(),
511 { proxsuite::linalg::sparse::from_eigen,
512 SparseMat<T, I>(A.value().transpose()) });
513 /* TO PUT IN DEBUG MODE
514 std::cout << "have same structure = " << res << std::endl;
515 */
516 if (res) {
517 copy(H_unscaled,
518 { proxsuite::linalg::sparse::from_eigen,
519 H_triu }); // copy rhs into lhs
520 copy(
521 AT_unscaled,
522 { proxsuite::linalg::sparse::from_eigen,
523 SparseMat<T, I>(A.value().transpose()) }); // copy rhs into lhs
524 }
525 }
526 } else if (C != nullopt) {
527 bool res =
529 H_unscaled.as_const(),
530 { proxsuite::linalg::sparse::from_eigen, H_triu }) &&
531 have_same_structure(CT_unscaled.as_const(),
532 { proxsuite::linalg::sparse::from_eigen,
533 SparseMat<T, I>(C.value().transpose()) });
534 /* TO PUT IN DEBUG MODE
535 std::cout << "have same structure = " << res << std::endl;
536 */
537 if (res) {
538 copy(H_unscaled,
539 { proxsuite::linalg::sparse::from_eigen,
540 H_triu }); // copy rhs into lhs
541 copy(CT_unscaled,
542 { proxsuite::linalg::sparse::from_eigen,
543 SparseMat<T, I>(C.value().transpose()) }); // copy rhs into lhs
544 }
545 } else {
546
547 bool res = have_same_structure(
548 H_unscaled.as_const(),
549 { proxsuite::linalg::sparse::from_eigen, H_triu });
550 /* TO PUT IN DEBUG MODE
551 std::cout << "have same structure = " << res << std::endl;
552 */
553 if (res) {
554 copy(H_unscaled,
555 { proxsuite::linalg::sparse::from_eigen,
556 H.value() }); // copy rhs into lhs
557 }
558 }
559 } else if (A != nullopt) {
560 if (C != nullopt) {
561 bool res =
562 have_same_structure(AT_unscaled.as_const(),
563 { proxsuite::linalg::sparse::from_eigen,
564 SparseMat<T, I>(A.value().transpose()) }) &&
565 have_same_structure(CT_unscaled.as_const(),
566 { proxsuite::linalg::sparse::from_eigen,
567 SparseMat<T, I>(C.value().transpose()) });
568 /* TO PUT IN DEBUG MODE
569 std::cout << "have same structure = " << res << std::endl;
570 */
571 if (res) {
572 copy(AT_unscaled,
573 { proxsuite::linalg::sparse::from_eigen,
574 SparseMat<T, I>(A.value().transpose()) }); // copy rhs into lhs
575 copy(CT_unscaled,
576 { proxsuite::linalg::sparse::from_eigen,
577 SparseMat<T, I>(C.value().transpose()) }); // copy rhs into lhs
578 }
579 } else {
580 bool res =
581 have_same_structure(AT_unscaled.as_const(),
582 { proxsuite::linalg::sparse::from_eigen,
583 SparseMat<T, I>(A.value().transpose()) });
584 /* TO PUT IN DEBUG MODE
585 std::cout << "have same structure = " << res << std::endl;
586 */
587 if (res) {
588 copy(AT_unscaled,
589 { proxsuite::linalg::sparse::from_eigen,
590 SparseMat<T, I>(A.value().transpose()) }); // copy rhs into lhs
591 }
592 }
593 } else if (C != nullopt) {
594 bool res =
595 have_same_structure(CT_unscaled.as_const(),
596 { proxsuite::linalg::sparse::from_eigen,
597 SparseMat<T, I>(C.value().transpose()) });
598 /* TO PUT IN DEBUG MODE
599 std::cout << "have same structure = " << res << std::endl;
600 */
601 if (res) {
602 copy(CT_unscaled,
603 { proxsuite::linalg::sparse::from_eigen,
604 SparseMat<T, I>(C.value().transpose()) }); // copy rhs into lhs
605 }
606 }
607
608 SparseMat<T, I> H_triu =
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 }
618 };
620 settings, results, work, rho, mu_eq, mu_in);
621 qp_setup(qp,
622 results,
623 model,
624 work,
625 settings,
626 ruiz,
627 preconditioner_status); // store model value + performs scaling
628 // according to chosen options
629 if (H != nullopt) {
632 manual_minimal_H_eigenvalue, results, settings);
633 }
634 if (settings.compute_timings) {
635 results.info.setup_time = work.timer.elapsed().user; // in microseconds
636 }
637 };
638
642 void solve()
643 {
644 qp_solve( //
645 results,
646 model,
647 settings,
648 work,
649 ruiz);
650 };
660 {
662 qp_solve( //
663 results,
664 model,
665 settings,
666 work,
667 ruiz);
668 };
672 void cleanup() { results.cleanup(settings); }
673};
710template<typename T, typename I>
723 optional<T> eps_abs = nullopt,
724 optional<T> eps_rel = nullopt,
725 optional<T> rho = nullopt,
726 optional<T> mu_eq = nullopt,
727 optional<T> mu_in = nullopt,
728 optional<bool> verbose = nullopt,
729 bool compute_preconditioner = true,
730 bool compute_timings = false,
731 optional<isize> max_iter = nullopt,
734 proxsuite::proxqp::SparseBackend sparse_backend =
736 bool check_duality_gap = false,
737 optional<T> eps_duality_gap_abs = nullopt,
738 optional<T> eps_duality_gap_rel = nullopt,
739 bool primal_infeasibility_solving = false,
740 optional<T> manual_minimal_H_eigenvalue = nullopt)
741{
742
743 isize n(0);
744 isize n_eq(0);
745 isize n_in(0);
746 if (H != nullopt) {
747 n = H.value().rows();
748 }
749 if (A != nullopt) {
750 n_eq = A.value().rows();
751 }
752 if (C != nullopt) {
753 n_in = C.value().rows();
754 }
755
756 proxqp::sparse::QP<T, I> Qp(n, n_eq, n_in);
757 Qp.settings.initial_guess = initial_guess;
758 Qp.settings.check_duality_gap = check_duality_gap;
759
760 if (eps_abs != nullopt) {
761 Qp.settings.eps_abs = eps_abs.value();
762 }
763 if (eps_rel != nullopt) {
764 Qp.settings.eps_rel = eps_rel.value();
765 }
766 if (verbose != nullopt) {
767 Qp.settings.verbose = verbose.value();
768 }
769 if (max_iter != nullopt) {
770 Qp.settings.max_iter = max_iter.value();
771 }
772 if (eps_duality_gap_abs != nullopt) {
773 Qp.settings.eps_duality_gap_abs = eps_duality_gap_abs.value();
774 }
775 if (eps_duality_gap_rel != nullopt) {
776 Qp.settings.eps_duality_gap_rel = eps_duality_gap_rel.value();
777 }
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) {
782 Qp.init(H,
783 g,
784 A,
785 b,
786 C,
787 l,
788 u,
789 compute_preconditioner,
790 rho,
791 mu_eq,
792 mu_in,
793 manual_minimal_H_eigenvalue.value());
794 } else {
795 Qp.init(
796 H, g, A, b, C, l, u, compute_preconditioner, rho, mu_eq, mu_in, nullopt);
797 }
798 Qp.solve(x, y, z);
799
800 return Qp.results;
801}
802
804template<typename T, typename I>
806{
811 std::vector<QP<T, I>> qp_vector;
813
814 BatchQP(long unsigned int batchSize)
815 {
816 if (qp_vector.max_size() != batchSize) {
817 qp_vector.clear();
818 qp_vector.reserve(batchSize);
819 }
820 m_size = 0;
821 }
822 BatchQP(const BatchQP&) = delete;
823 BatchQP& operator=(const BatchQP&) = delete;
824 BatchQP(BatchQP&&) = default;
825 BatchQP& operator=(BatchQP&&) = default;
826
831 sparse::isize n_eq,
832 sparse::isize n_in)
833 {
834 qp_vector.emplace_back(dim, n_eq, n_in);
835 auto& qp = qp_vector.back();
836 m_size++;
837 return qp;
838 };
839
840 // /*!
841 // * Init a QP in place and return a reference to it
842 // */
843 // QP<T, I>& init_qp_in_place(const sparse::SparseMat<bool, I>& H,
844 // const sparse::SparseMat<bool, I>& A,
845 // const sparse::SparseMat<bool, I>& C)
846 // {
847 // qp_vector.emplace_back(H.rows(), A.rows(), C.rows());
848 // auto& qp = qp_vector.back();
849 // m_size++;
850 // return qp;
851 // };
852
856 void insert(QP<T, I>& qp) { qp_vector.emplace_back(qp); };
857
861 QP<T, I>& get(isize i) { return qp_vector.at(size_t(i)); };
862
866 const QP<T, I>& get(isize i) const { return qp_vector.at(size_t(i)); };
867
871 QP<T, I>& operator[](isize i) { return get(i); };
872
876 const QP<T, I>& operator[](isize i) const { return get(i); };
877
878 sparse::isize size() { return m_size; };
879};
880
881} // namespace sparse
882} // namespace proxqp
883} // namespace proxsuite
884
885#endif /* end of include guard PROXSUITE_PROXQP_SPARSE_WRAPPER_HPP */
#define PROXSUITE_CHECK_ARGUMENT_SIZE(size, expected_size, message)
Definition macros.hpp:28
_detail::_meta::make_signed< usize >::Type isize
Definition typedefs.hpp:43
auto middle_cols_mut(proxsuite::linalg::sparse::MatMut< T, I > mat, isize start, isize ncols, isize nnz) -> proxsuite::linalg::sparse::MatMut< T, I >
Definition utils.hpp:400
auto top_rows_mut_unchecked(proxsuite::linalg::veg::Unsafe, proxsuite::linalg::sparse::MatMut< T, I > mat, isize nrows) -> proxsuite::linalg::sparse::MatMut< T, I >
Definition utils.hpp:440
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)
Definition helpers.hpp:219
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)
Definition helpers.hpp:183
void qp_solve(Results< T > &results, Model< T, I > &data, const Settings< T > &settings, Workspace< T, I > &work, P &precond)
Definition solver.hpp:344
Eigen::SparseMatrix< T, Eigen::ColMajor, I > SparseMat
Definition fwd.hpp:31
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)
Definition wrapper.hpp:712
void copy(proxsuite::linalg::sparse::MatMut< T, I > a, proxsuite::linalg::sparse::MatRef< T, I > b)
Definition helpers.hpp:443
Eigen::Ref< Eigen::Matrix< T, DYN, 1 > const > VecRef
Definition fwd.hpp:34
auto have_same_structure(proxsuite::linalg::sparse::MatRef< T, I > a, proxsuite::linalg::sparse::MatRef< T, I > b) -> bool
Definition helpers.hpp:414
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)
Definition helpers.hpp:282
void update_default_rho_with_minimal_Hessian_eigen_value(optional< T > manual_minimal_H_eigenvalue, Results< T > &results, Settings< T > &settings)
Definition helpers.hpp:159
constexpr nullopt_t nullopt
Definition optional.hpp:42
auto as_const() const noexcept -> MatRef< T, I >
Definition core.hpp:484
auto to_eigen() const noexcept -> Eigen::Map< Eigen::SparseMatrix< T, Eigen::ColMajor, I > >
Definition core.hpp:508
This class stores all the results of PROXQP solvers with sparse and dense backends.
Definition results.hpp:68
This class defines the settings of PROXQP solvers with sparse and dense backends.
Definition settings.hpp:89
std::vector< QP< T, I > > qp_vector
Definition wrapper.hpp:811
void insert(QP< T, I > &qp)
Definition wrapper.hpp:856
QP< T, I > & operator[](isize i)
Definition wrapper.hpp:871
BatchQP & operator=(const BatchQP &)=delete
QP< T, I > & init_qp_in_place(sparse::isize dim, sparse::isize n_eq, sparse::isize n_in)
Definition wrapper.hpp:830
QP< T, I > & get(isize i)
Definition wrapper.hpp:861
const QP< T, I > & operator[](isize i) const
Definition wrapper.hpp:876
BatchQP & operator=(BatchQP &&)=default
const QP< T, I > & get(isize i) const
Definition wrapper.hpp:866
BatchQP(long unsigned int batchSize)
Definition wrapper.hpp:814
BatchQP(const BatchQP &)=delete
This class stores the model of the QP problem.
Definition model.hpp:23
This class defines the API of PROXQP solver with sparse backend.
Definition wrapper.hpp:91
QP & operator=(const QP &)=delete
Workspace< T, I > work
Definition wrapper.hpp:95
QP(const SparseMat< bool, I > &H, const SparseMat< bool, I > &A, const SparseMat< bool, I > &C)
Definition wrapper.hpp:122
void solve(optional< VecRef< T > > x, optional< VecRef< T > > y, optional< VecRef< T > > z)
Definition wrapper.hpp:657
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)
Definition wrapper.hpp:357
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)
Definition wrapper.hpp:170
QP(isize dim, isize n_eq, isize n_in)
Definition wrapper.hpp:103
preconditioner::RuizEquilibration< T, I > ruiz
Definition wrapper.hpp:96
This class defines the workspace of the sparse solver.