proxsuite 0.6.7
The Advanced Proximal Optimization Toolbox
Loading...
Searching...
No Matches
helpers.hpp
Go to the documentation of this file.
1//
2// Copyright (c) 2022-2023 INRIA
3//
8#ifndef PROXSUITE_PROXQP_DENSE_HELPERS_HPP
9#define PROXSUITE_PROXQP_DENSE_HELPERS_HPP
10
16#include <chrono>
18#include <Eigen/Eigenvalues>
19
20namespace proxsuite {
21namespace proxqp {
22namespace dense {
23
24template<typename T,
25 typename MatIn,
26 typename VecIn1,
27 typename VecIn2,
28 typename VecIn3>
29T
30power_iteration(const Eigen::MatrixBase<MatIn>& H,
31 const Eigen::MatrixBase<VecIn1>& dw,
32 const Eigen::MatrixBase<VecIn2>& rhs,
33 const Eigen::MatrixBase<VecIn3>& err_v,
34 T power_iteration_accuracy,
35 isize nb_power_iteration)
36{
37 auto& dw_cc = dw.const_cast_derived();
38 auto& rhs_cc = rhs.const_cast_derived();
39 auto& err_v_cc = err_v.const_cast_derived();
40 // computes maximal eigen value of the bottom right matrix of the LDLT
41 isize dim = H.rows();
42 rhs_cc.setZero();
43 // stores eigenvector
44 rhs_cc.array() += 1. / std::sqrt(dim);
45 // stores Hx
46 dw_cc.noalias() = H.template selfadjointView<Eigen::Lower>() * rhs_cc; // Hx
47 T eig = 0;
48 for (isize i = 0; i < nb_power_iteration; i++) {
49
50 rhs_cc = dw_cc / dw_cc.norm();
51 dw_cc.noalias() = (H.template selfadjointView<Eigen::Lower>() * rhs_cc);
52 // calculate associated eigenvalue
53 eig = rhs.dot(dw_cc);
54 // calculate associated error
55 err_v_cc = dw_cc - eig * rhs_cc;
56 T err = proxsuite::proxqp::dense::infty_norm(err_v_cc);
57 // std::cout << "power iteration max: i " << i << " err " << err <<
58 // std::endl;
59 if (err <= power_iteration_accuracy) {
60 break;
61 }
62 }
63 return eig;
64}
65template<typename T,
66 typename MatIn,
67 typename VecIn1,
68 typename VecIn2,
69 typename VecIn3>
70T
72 const Eigen::MatrixBase<MatIn>& H,
73 const Eigen::MatrixBase<VecIn1>& dw,
74 const Eigen::MatrixBase<VecIn2>& rhs,
75 const Eigen::MatrixBase<VecIn3>& err_v,
76 T max_eigen_value,
77 T power_iteration_accuracy,
78 isize nb_power_iteration)
79{
80 // performs power iteration on the matrix: max_eigen_value I - H
81 // estimates then the minimal eigenvalue with: minimal_eigenvalue =
82 // max_eigen_value - eig
83 auto& dw_cc = dw.const_cast_derived();
84 auto& rhs_cc = rhs.const_cast_derived();
85 auto& err_v_cc = err_v.const_cast_derived();
86 isize dim = H.rows();
87 rhs_cc.setZero();
88 // stores eigenvector
89 rhs_cc.array() += 1. / std::sqrt(dim);
90 // stores Hx
91 dw_cc.noalias() =
92 -(H.template selfadjointView<Eigen::Lower>() * rhs_cc); // Hx
93 dw_cc += max_eigen_value * rhs_cc;
94 T eig = 0;
95 for (isize i = 0; i < nb_power_iteration; i++) {
96
97 rhs_cc = dw_cc / dw_cc.norm();
98 dw_cc.noalias() = -(H.template selfadjointView<Eigen::Lower>() * rhs_cc);
99 dw_cc += max_eigen_value * rhs_cc;
100 // calculate associated eigenvalue
101 eig = rhs_cc.dot(dw_cc);
102 // calculate associated error
103 err_v_cc = dw_cc - eig * rhs_cc;
104 T err = proxsuite::proxqp::dense::infty_norm(err_v_cc);
105 // std::cout << "power iteration min: i " << i << " err " << err <<
106 // std::endl;
107 if (err <= power_iteration_accuracy) {
108 break;
109 }
110 }
111 T minimal_eigenvalue = max_eigen_value - eig;
112 return minimal_eigenvalue;
113}
115
123template<typename T, typename MatIn>
124T
126 const Eigen::MatrixBase<MatIn>& H,
127 EigenValueEstimateMethodOption estimate_method_option,
128 T power_iteration_accuracy,
129 isize nb_power_iteration)
130{
132 (!H.isApprox(H.transpose(), std::numeric_limits<T>::epsilon())),
133 std::invalid_argument,
134 "H is not symmetric.");
135 if (H.size()) {
137 H.rows(),
138 H.cols(),
139 "H has a number of rows different of the number of columns.");
140 }
141 isize dim = H.rows();
142 T res(0.);
143 switch (estimate_method_option) {
145 Vec<T> dw(dim);
146 Vec<T> rhs(dim);
147 Vec<T> err(dim);
148 T dominant_eigen_value = power_iteration(
149 H, dw, rhs, err, power_iteration_accuracy, nb_power_iteration);
150 T min_eigenvalue =
152 dw,
153 rhs,
154 err,
155 dominant_eigen_value,
156 power_iteration_accuracy,
157 nb_power_iteration);
158 res = std::min(min_eigenvalue, dominant_eigen_value);
159 } break;
161 Eigen::SelfAdjointEigenSolver<Mat<T>> es(H, Eigen::EigenvaluesOnly);
162 res = T(es.eigenvalues()[0]);
163 } break;
164 }
165 return res;
166}
168
174template<typename T>
175void
177 optional<T> manual_minimal_H_eigenvalue,
178 Results<T>& results,
179 Settings<T>& settings)
180{
181 if (manual_minimal_H_eigenvalue != nullopt) {
183 manual_minimal_H_eigenvalue.value();
184 results.info.minimal_H_eigenvalue_estimate =
186 }
187 settings.default_rho += std::abs(results.info.minimal_H_eigenvalue_estimate);
188 results.info.rho = settings.default_rho;
189}
199template<typename T>
200void
202 const Settings<T>& qpsettings,
203 const Model<T>& qpmodel,
204 const isize n_constraints,
205 const DenseBackend& dense_backend,
206 const HessianType& hessian_type,
207 Results<T>& qpresults)
208{
209
210 qpwork.rhs.setZero();
211 qpwork.rhs.head(qpmodel.dim) = -qpwork.g_scaled;
212 qpwork.rhs.segment(qpmodel.dim, qpmodel.n_eq) = qpwork.b_scaled;
214 qpsettings,
215 qpmodel,
216 qpresults,
217 qpwork,
218 n_constraints,
219 dense_backend,
220 hessian_type,
221 T(1),
222 qpmodel.dim + qpmodel.n_eq);
223
224 qpresults.x = qpwork.dw_aug.head(qpmodel.dim);
225 qpresults.y = qpwork.dw_aug.segment(qpmodel.dim, qpmodel.n_eq);
226 qpwork.dw_aug.setZero();
227 qpwork.rhs.setZero();
228}
229
239template<typename T>
240void
242 const Model<T>& qpmodel,
243 Results<T>& qpresults,
244 const DenseBackend& dense_backend,
245 const HessianType& hessian_type)
246{
247
250 qpwork.ldl_stack.as_mut(),
251 };
252 switch (hessian_type) {
254 qpwork.kkt.topLeftCorner(qpmodel.dim, qpmodel.dim) = qpwork.H_scaled;
255 break;
257 qpwork.kkt.topLeftCorner(qpmodel.dim, qpmodel.dim).setZero();
258 break;
260 qpwork.kkt.topLeftCorner(qpmodel.dim, qpmodel.dim) = qpwork.H_scaled;
261 break;
262 }
263 qpwork.kkt.topLeftCorner(qpmodel.dim, qpmodel.dim).diagonal().array() +=
264 qpresults.info.rho;
265 switch (dense_backend) {
267 qpwork.kkt.block(0, qpmodel.dim, qpmodel.dim, qpmodel.n_eq) =
268 qpwork.A_scaled.transpose();
269 qpwork.kkt.block(qpmodel.dim, 0, qpmodel.n_eq, qpmodel.dim) =
270 qpwork.A_scaled;
271 qpwork.kkt.bottomRightCorner(qpmodel.n_eq, qpmodel.n_eq).setZero();
272 qpwork.kkt.diagonal()
273 .segment(qpmodel.dim, qpmodel.n_eq)
274 .setConstant(-qpresults.info.mu_eq);
275 qpwork.ldl.factorize(qpwork.kkt.transpose(), stack);
276 break;
278 qpwork.kkt.noalias() += qpresults.info.mu_eq_inv *
279 (qpwork.A_scaled.transpose() * qpwork.A_scaled);
280 qpwork.ldl.factorize(qpwork.kkt.transpose(), stack);
281 break;
283 break;
284 }
285}
298template<typename T>
299void
301 const Settings<T>& qpsettings,
302 const bool box_constraints,
303 const HessianType hessian_type,
305 bool execute_preconditioner)
306{
307
308 QpViewBoxMut<T> qp_scaled{
309 { from_eigen, qpwork.H_scaled }, { from_eigen, qpwork.g_scaled },
310 { from_eigen, qpwork.A_scaled }, { from_eigen, qpwork.b_scaled },
311 { from_eigen, qpwork.C_scaled }, { from_eigen, qpwork.u_scaled },
312 { from_eigen, qpwork.l_scaled }, { from_eigen, qpwork.i_scaled },
313 { from_eigen, qpwork.l_box_scaled }, { from_eigen, qpwork.u_box_scaled },
314 };
315
318 qpwork.ldl_stack.as_mut(),
319 };
320 ruiz.scale_qp_in_place(qp_scaled,
321 execute_preconditioner,
323 qpsettings.preconditioner_max_iter,
324 qpsettings.preconditioner_accuracy,
325 hessian_type,
326 box_constraints,
327 stack);
328 qpwork.correction_guess_rhs_g = infty_norm(qpwork.g_scaled);
329}
330
340template<typename T>
341void
343 Settings<T>& qpsettings,
344 Model<T>& qpmodel,
345 Results<T>& qpresults)
346{
347
348 switch (qpsettings.initial_guess) {
351 qpwork, qpsettings, qpmodel, qpresults);
352 break;
353 }
354 }
355}
372template<typename T>
373void
381 optional<VecRef<T>> l_box,
382 optional<VecRef<T>> u_box,
383 Model<T>& model,
384 Workspace<T>& work,
385 const bool box_constraints)
386{
387 // check the model is valid
388 if (g != nullopt) {
389 PROXSUITE_CHECK_ARGUMENT_SIZE(g.value().size(),
390 model.dim,
391 "the dimension wrt the primal variable x "
392 "variable for updating g is not valid.");
393 }
394 if (b != nullopt) {
395 PROXSUITE_CHECK_ARGUMENT_SIZE(b.value().size(),
396 model.n_eq,
397 "the dimension wrt equality constrained "
398 "variables for updating b is not valid.");
399 }
400 if (u != nullopt) {
401 PROXSUITE_CHECK_ARGUMENT_SIZE(u.value().size(),
402 model.n_in,
403 "the dimension wrt inequality constrained "
404 "variables for updating u is not valid.");
405 }
406 if (l != nullopt) {
407 PROXSUITE_CHECK_ARGUMENT_SIZE(l.value().size(),
408 model.n_in,
409 "the dimension wrt inequality constrained "
410 "variables for updating l is not valid.");
411 }
412 if (H != nullopt) {
414 H.value().rows(),
415 model.dim,
416 "the row dimension for updating H is not valid.");
418 H.value().cols(),
419 model.dim,
420 "the column dimension for updating H is not valid.");
421 }
422 if (A != nullopt) {
424 A.value().rows(),
425 model.n_eq,
426 "the row dimension for updating A is not valid.");
428 A.value().cols(),
429 model.dim,
430 "the column dimension for updating A is not valid.");
431 }
432 if (C != nullopt) {
434 C.value().rows(),
435 model.n_in,
436 "the row dimension for updating C is not valid.");
438 C.value().cols(),
439 model.dim,
440 "the column dimension for updating C is not valid.");
441 }
442
443 // update the model
444 if (g != nullopt) {
445 model.g = g.value().eval();
446 }
447 if (b != nullopt) {
448 model.b = b.value().eval();
449 }
450 if (u != nullopt) {
451 model.u = u.value().eval();
452 }
453 if (l != nullopt) {
454 model.l = l.value().eval();
455 }
456 if (u_box != nullopt && box_constraints) {
457 model.u_box = u_box.value();
458 } // else qpmodel.u_box remains initialized to a matrix with zero elements or
459 // zero shape
460
461 if (l_box != nullopt && box_constraints) {
462 model.l_box = l_box.value();
463 } // else qpmodel.l_box remains initialized to a matrix with zero elements or
464 // zero shape
465
466 if (H != nullopt || A != nullopt || C != nullopt) {
467 work.refactorize = true;
468 }
469
470 if (H != nullopt) {
471 model.H = H.value();
472 }
473 if (A != nullopt) {
474 model.A = A.value();
475 }
476 if (C != nullopt) {
477 model.C = C.value();
478 }
479 assert(model.is_valid(box_constraints));
480}
500template<typename T>
501void
510 optional<VecRef<T>> l_box,
511 optional<VecRef<T>> u_box,
512 Settings<T>& qpsettings,
513 Model<T>& qpmodel,
514 Workspace<T>& qpwork,
515 Results<T>& qpresults,
516 const bool box_constraints,
518 PreconditionerStatus preconditioner_status,
519 const HessianType hessian_type)
520{
521
522 switch (qpsettings.initial_guess) {
524 if (qpwork.proximal_parameter_update) {
526 } else {
527 qpresults.cleanup(qpsettings);
528 }
529 qpwork.cleanup(box_constraints);
530 break;
531 }
533 // keep solutions but restart workspace and results
534 if (qpwork.proximal_parameter_update) {
535 qpresults.cleanup_statistics();
536 } else {
537 qpresults.cold_start(qpsettings);
538 }
539 qpwork.cleanup(box_constraints);
540 break;
541 }
543 if (qpwork.proximal_parameter_update) {
545 } else {
546 qpresults.cleanup(qpsettings);
547 }
548 qpwork.cleanup(box_constraints);
549 break;
550 }
552 if (qpwork.proximal_parameter_update) {
553 qpresults
554 .cleanup_all_except_prox_parameters(); // the warm start is given at
555 // the solve function
556 } else {
557 qpresults.cleanup(qpsettings);
558 }
559 qpwork.cleanup(box_constraints);
560 break;
561 }
563 if (qpwork.refactorize || qpwork.proximal_parameter_update) {
564 qpwork.cleanup(box_constraints); // meaningful for when there is an
565 // upate of the model and one wants to
566 // warm start with previous result
567 qpwork.refactorize = true;
568 }
569 qpresults.cleanup_statistics();
570 break;
571 }
572 }
573 if (H != nullopt) {
574 qpmodel.H = H.value();
575 } // else qpmodel.H remains initialzed to a matrix with zero elements
576 if (g != nullopt) {
577 qpmodel.g = g.value();
578 }
579
580 if (A != nullopt) {
581 qpmodel.A = A.value();
582 } // else qpmodel.A remains initialized to a matrix with zero elements or zero
583 // shape
584
585 if (b != nullopt) {
586 qpmodel.b = b.value();
587 } // else qpmodel.b remains initialized to a matrix with zero elements or zero
588 // shape
589
590 if (C != nullopt) {
591 qpmodel.C = C.value();
592 } // else qpmodel.C remains initialized to a matrix with zero elements or zero
593 // shape
594
595 if (u != nullopt) {
596 qpmodel.u = u.value();
597 } // else qpmodel.u remains initialized to a matrix with zero elements or zero
598 // shape
599
600 if (l != nullopt) {
601 qpmodel.l = l.value();
602 } // else qpmodel.l remains initialized to a matrix with zero elements or zero
603 // shape
604 if (u_box != nullopt) {
605 qpmodel.u_box = u_box.value();
606 } // else qpmodel.u_box remains initialized to a matrix with zero elements or
607 // zero shape
608
609 if (l_box != nullopt) {
610 qpmodel.l_box = l_box.value();
611 } // else qpmodel.l_box remains initialized to a matrix with zero elements or
612 // zero shape
613 assert(qpmodel.is_valid(box_constraints));
614 switch (hessian_type) {
616 break;
618 qpwork.H_scaled = qpmodel.H;
619 break;
621 qpwork.H_scaled = qpmodel.H;
622 break;
623 }
624 qpwork.g_scaled = qpmodel.g;
625 qpwork.A_scaled = qpmodel.A;
626 qpwork.b_scaled = qpmodel.b;
627 qpwork.C_scaled = qpmodel.C;
628 qpwork.u_scaled =
629 (qpmodel.u.array() <= T(1.E20))
630 .select(qpmodel.u,
631 Eigen::Matrix<T, Eigen::Dynamic, 1>::Zero(qpmodel.n_in).array() +
632 T(1.E20));
633 qpwork.l_scaled =
634 (qpmodel.l.array() >= T(-1.E20))
635 .select(qpmodel.l,
636 Eigen::Matrix<T, Eigen::Dynamic, 1>::Zero(qpmodel.n_in).array() -
637 T(1.E20));
638 if (box_constraints) {
639 qpwork.u_box_scaled =
640 (qpmodel.u_box.array() <= T(1.E20))
641 .select(qpmodel.u_box,
642 Eigen::Matrix<T, Eigen::Dynamic, 1>::Zero(qpmodel.dim).array() +
643 T(1.E20));
644 qpwork.l_box_scaled =
645 (qpmodel.l_box.array() >= T(-1.E20))
646 .select(qpmodel.l_box,
647 Eigen::Matrix<T, Eigen::Dynamic, 1>::Zero(qpmodel.dim).array() -
648 T(1.E20));
649 }
650
651 qpwork.dual_feasibility_rhs_2 = infty_norm(qpmodel.g);
652 switch (preconditioner_status) {
655 qpwork, qpsettings, box_constraints, hessian_type, ruiz, true);
656 break;
659 qpwork, qpsettings, box_constraints, hessian_type, ruiz, false);
660 break;
662 // keep previous one
664 qpwork, qpsettings, box_constraints, hessian_type, ruiz, false);
665 break;
666 }
667}
669
678template<typename T>
679void
681 Results<T>& results,
682 Workspace<T>& work,
683 optional<T> rho_new,
684 optional<T> mu_eq_new,
685 optional<T> mu_in_new)
686{
687
688 if (rho_new != nullopt) {
689 settings.default_rho = rho_new.value();
690 results.info.rho = rho_new.value();
691 work.proximal_parameter_update = true;
692 }
693 if (mu_eq_new != nullopt) {
694 settings.default_mu_eq = mu_eq_new.value();
695 results.info.mu_eq = mu_eq_new.value();
696 results.info.mu_eq_inv = T(1) / results.info.mu_eq;
697 work.proximal_parameter_update = true;
698 }
699 if (mu_in_new != nullopt) {
700 settings.default_mu_in = mu_in_new.value();
701 results.info.mu_in = mu_in_new.value();
702 results.info.mu_in_inv = T(1) / results.info.mu_in;
703 work.proximal_parameter_update = true;
704 }
705}
715template<typename T>
716void
718 optional<VecRef<T>> y_wm,
719 optional<VecRef<T>> z_wm,
720 Results<T>& results,
721 Settings<T>& settings,
722 Model<T>& model)
723{
724 if (x_wm == nullopt && y_wm == nullopt && z_wm == nullopt)
725 return;
726
728
729 // first check problem dimensions
730 if (x_wm != nullopt) {
732 x_wm.value().rows(),
733 model.dim,
734 "the dimension wrt primal variable x for warm start is not valid.");
735 }
736
737 if (y_wm != nullopt) {
738 PROXSUITE_CHECK_ARGUMENT_SIZE(y_wm.value().rows(),
739 model.n_eq,
740 "the dimension wrt equality constrained "
741 "variables for warm start is not valid.");
742 }
743
744 if (z_wm != nullopt) {
746 z_wm.value().rows(),
747 model.n_in,
748 "the dimension wrt inequality constrained variables for warm start "
749 "is not valid.");
750 }
751
752 if (x_wm != nullopt) {
753 results.x = x_wm.value().eval();
754 }
755
756 if (y_wm != nullopt) {
757 results.y = y_wm.value().eval();
758 }
759
760 if (z_wm != nullopt) {
761 results.z = z_wm.value().eval();
762 }
763}
764} // namespace dense
765} // namespace proxqp
766} // namespace proxsuite
767
768#endif /* end of include guard PROXSUITE_PROXQP_DENSE_HELPERS_HPP */
TL_OPTIONAL_11_CONSTEXPR T & value() &
#define PROXSUITE_CHECK_ARGUMENT_SIZE(size, expected_size, message)
Definition macros.hpp:28
#define PROXSUITE_THROW_PRETTY(condition, exception, message)
Definition macros.hpp:18
void setup_factorization(Workspace< T > &qpwork, const Model< T > &qpmodel, Results< T > &qpresults, const DenseBackend &dense_backend, const HessianType &hessian_type)
Definition helpers.hpp:241
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 > &model)
Definition helpers.hpp:717
void initial_guess(Workspace< T > &qpwork, Settings< T > &qpsettings, Model< T > &qpmodel, Results< T > &qpresults)
Definition helpers.hpp:342
void iterative_solve_with_permut_fact(const Settings< T > &qpsettings, const Model< T > &qpmodel, Results< T > &qpresults, Workspace< T > &qpwork, const isize n_constraints, const DenseBackend &dense_backend, const HessianType &hessian_type, T eps, isize inner_pb_dim)
Definition solver.hpp:408
void setup_equilibration(Workspace< T > &qpwork, const Settings< T > &qpsettings, const bool box_constraints, const HessianType hessian_type, preconditioner::RuizEquilibration< T > &ruiz, bool execute_preconditioner)
Definition helpers.hpp:300
Eigen::Ref< Mat< T, l > const > MatRef
Definition fwd.hpp:33
T estimate_minimal_eigen_value_of_symmetric_matrix(const Eigen::MatrixBase< MatIn > &H, EigenValueEstimateMethodOption estimate_method_option, T power_iteration_accuracy, isize nb_power_iteration)
Definition helpers.hpp:125
void update_proximal_parameters(Settings< T > &settings, Results< T > &results, Workspace< T > &work, optional< T > rho_new, optional< T > mu_eq_new, optional< T > mu_in_new)
Definition helpers.hpp:680
void update_default_rho_with_minimal_Hessian_eigen_value(optional< T > manual_minimal_H_eigenvalue, Results< T > &results, Settings< T > &settings)
Definition helpers.hpp:176
Eigen::Matrix< T, DYN, 1 > Vec
Definition fwd.hpp:24
void setup(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, optional< VecRef< T > > l_box, optional< VecRef< T > > u_box, Settings< T > &qpsettings, Model< T > &qpmodel, Workspace< T > &qpwork, Results< T > &qpresults, const bool box_constraints, preconditioner::RuizEquilibration< T > &ruiz, PreconditionerStatus preconditioner_status, const HessianType hessian_type)
Definition helpers.hpp:502
T power_iteration(const Eigen::MatrixBase< MatIn > &H, const Eigen::MatrixBase< VecIn1 > &dw, const Eigen::MatrixBase< VecIn2 > &rhs, const Eigen::MatrixBase< VecIn3 > &err_v, T power_iteration_accuracy, isize nb_power_iteration)
Definition helpers.hpp:30
T min_eigen_value_via_modified_power_iteration(const Eigen::MatrixBase< MatIn > &H, const Eigen::MatrixBase< VecIn1 > &dw, const Eigen::MatrixBase< VecIn2 > &rhs, const Eigen::MatrixBase< VecIn3 > &err_v, T max_eigen_value, T power_iteration_accuracy, isize nb_power_iteration)
Definition helpers.hpp:71
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, optional< VecRef< T > > l_box, optional< VecRef< T > > u_box, Model< T > &model, Workspace< T > &work, const bool box_constraints)
Definition helpers.hpp:374
void compute_equality_constrained_initial_guess(Workspace< T > &qpwork, const Settings< T > &qpsettings, const Model< T > &qpmodel, const isize n_constraints, const DenseBackend &dense_backend, const HessianType &hessian_type, Results< T > &qpresults)
Definition helpers.hpp:201
Eigen::Ref< Vec< T > const > VecRef
Definition fwd.hpp:26
constexpr nullopt_t nullopt
Definition optional.hpp:42
VEG_NODISCARD VEG_INLINE auto as_mut() VEG_NOEXCEPT -> SliceMut< T >
Definition vec.hpp:957
This class stores all the results of PROXQP solvers with sparse and dense backends.
Definition results.hpp:68
sparse::Vec< T > z
Definition results.hpp:74
void cold_start(optional< Settings< T > > settings=nullopt)
Definition results.hpp:175
void cleanup(optional< Settings< T > > settings=nullopt)
Definition results.hpp:149
sparse::Vec< T > y
Definition results.hpp:73
sparse::Vec< T > x
Definition results.hpp:72
void cleanup_all_except_prox_parameters()
Definition results.hpp:195
This class defines the settings of PROXQP solvers with sparse and dense backends.
Definition settings.hpp:89
InitialGuessStatus initial_guess
Definition settings.hpp:123
This class stores the model of the QP problem.
Definition model.hpp:24
bool is_valid(const bool box_constraints)
Definition model.hpp:104
This class defines the workspace of the dense solver.
Definition workspace.hpp:26
proxsuite::linalg::veg::Vec< unsigned char > ldl_stack
Definition workspace.hpp:30
void cleanup(const bool box_constraints)
proxsuite::linalg::dense::Ldlt< T > ldl
Definition workspace.hpp:29
void scale_qp_in_place(QpViewBoxMut< T > qp, bool execute_preconditioner, bool preconditioning_for_infeasible_problems, const isize max_iter, const T epsilon, const HessianType &HessianType, const bool box_constraints, proxsuite::linalg::veg::dynstack::DynStackMut stack)
Definition ruiz.hpp:403