proxsuite 0.6.7
The Advanced Proximal Optimization Toolbox
Loading...
Searching...
No Matches
utils.hpp
Go to the documentation of this file.
1//
2// Copyright (c) 2022-2024 INRIA
3//
6#ifndef PROXSUITE_PROXQP_SPARSE_UTILS_HPP
7#define PROXSUITE_PROXQP_SPARSE_UTILS_HPP
8
9#include <iostream>
10#include <Eigen/IterativeLinearSolvers>
11#include <unsupported/Eigen/IterativeSolvers>
12
29
30namespace proxsuite {
31namespace proxqp {
32namespace sparse {
33
34template<typename T, typename I>
35void
37 Results<T>& results,
38 const Model<T, I>& model)
39{
40
42
43 // Print variables and constraints
44 std::cout << "problem: " << std::noshowpos << std::endl;
45 std::cout << " variables n = " << model.dim
46 << ", equality constraints n_eq = " << model.n_eq << ",\n"
47 << " inequality constraints n_in = " << model.n_in
48 << ", nnz = " << model.H_nnz + model.A_nnz + model.C_nnz << ",\n"
49 << std::endl;
50
51 // Print Settings
52 std::cout << "settings: " << std::endl;
53 std::cout << " backend = sparse," << std::endl;
54 std::cout << " sparse_backend = " << settings.sparse_backend;
56 std::cout << " -> " << results.info.sparse_backend;
57 }
58 std::cout << "," << std::endl;
59 std::cout << " eps_abs = " << settings.eps_abs
60 << ", eps_rel = " << settings.eps_rel << std::endl;
61 std::cout << " eps_prim_inf = " << settings.eps_primal_inf
62 << ", eps_dual_inf = " << settings.eps_dual_inf << "," << std::endl;
63
64 std::cout << " rho = " << results.info.rho
65 << ", mu_eq = " << results.info.mu_eq
66 << ", mu_in = " << results.info.mu_in << "," << std::endl;
67 std::cout << " max_iter = " << settings.max_iter
68 << ", max_iter_in = " << settings.max_iter_in << "," << std::endl;
69
70 if (settings.compute_preconditioner) {
71 std::cout << " scaling: on, " << std::endl;
72 } else {
73 std::cout << " scaling: off, " << std::endl;
74 }
75 if (settings.compute_timings) {
76 std::cout << " timings: on, " << std::endl;
77 } else {
78 std::cout << " timings: off, " << std::endl;
79 }
80 switch (settings.initial_guess) {
82 std::cout << " initial guess: warm start. \n" << std::endl;
83 break;
85 std::cout << " initial guess: initial guess. \n" << std::endl;
86 break;
88 std::cout
89 << " initial guess: warm start with previous result. \n"
90 << std::endl;
91 break;
93 std::cout
94 << " initial guess: cold start with previous result. \n"
95 << std::endl;
96 break;
98 std::cout
99 << " initial guess: equality constrained initial guess. \n"
100 << std::endl;
101 }
102}
103
104namespace detail {
105
106template<typename T, typename I>
107VEG_NO_INLINE void
109 VectorViewMut<T> out_l,
110 VectorViewMut<T> out_r,
112 VectorView<T> in_l,
113 VectorView<T> in_r)
114{
115 VEG_ASSERT_ALL_OF /* NOLINT */ (a.nrows() == out_r.dim,
116 a.ncols() == in_r.dim,
117 a.ncols() == out_l.dim,
118 a.nrows() == in_l.dim);
119 // equivalent to
120 // out_r.to_eigen().noalias() += a.to_eigen() * in_r.to_eigen();
121 // out_l.to_eigen().noalias() += a.to_eigen().transpose() * in_l.to_eigen();
122
123 auto* ai = a.row_indices();
124 auto* ax = a.values();
125 auto n = a.ncols();
126
127 for (usize j = 0; j < usize(n); ++j) {
128 usize col_start = a.col_start(j);
129 usize col_end = a.col_end(j);
130
131 T acc0 = 0;
132 T acc1 = 0;
133 T acc2 = 0;
134 T acc3 = 0;
135
136 T in_rj = in_r(isize(j));
137
138 usize pcount = col_end - col_start;
139
140 usize p = col_start;
141
142 auto zx = proxsuite::linalg::sparse::util::zero_extend;
143
144 for (; p < col_start + pcount / 4 * 4; p += 4) {
145 auto i0 = isize(zx(ai[p + 0]));
146 auto i1 = isize(zx(ai[p + 1]));
147 auto i2 = isize(zx(ai[p + 2]));
148 auto i3 = isize(zx(ai[p + 3]));
149
150 T ai0j = ax[p + 0];
151 T ai1j = ax[p + 1];
152 T ai2j = ax[p + 2];
153 T ai3j = ax[p + 3];
154
155 out_r(i0) += ai0j * in_rj;
156 out_r(i1) += ai1j * in_rj;
157 out_r(i2) += ai2j * in_rj;
158 out_r(i3) += ai3j * in_rj;
159
160 acc0 += ai0j * in_l(i0);
161 acc1 += ai1j * in_l(i1);
162 acc2 += ai2j * in_l(i2);
163 acc3 += ai3j * in_l(i3);
164 }
165
166 for (; p < col_end; ++p) {
167 auto i = isize(zx(ai[p]));
168
169 T aij = ax[p];
170 out_r(i) += aij * in_rj;
171 acc0 += aij * in_l(i);
172 }
173
174 acc0 = ((acc0 + acc1) + (acc2 + acc3));
175 out_l(isize(j)) += acc0;
176 }
177}
178
179template<typename T, typename I>
180VEG_NO_INLINE void
184 VectorView<T> in)
185{
186 VEG_ASSERT_ALL_OF /* NOLINT */ ( //
187 a.nrows() == a.ncols(),
188 a.nrows() == out.dim,
189 a.ncols() == in.dim);
190 // equivalent to
191 // out.to_eigen().noalias() +=
192 // a.to_eigen().template selfadjointView<Eigen::Upper>() *
193 // in.to_eigen();
194
195 auto* ai = a.row_indices();
196 auto* ax = a.values();
197 auto n = a.ncols();
198
199 for (usize j = 0; j < usize(n); ++j) {
200 usize col_start = a.col_start(j);
201 usize col_end = a.col_end(j);
202
203 if (col_start == col_end) {
204 continue;
205 }
206
207 T acc0 = 0;
208 T acc1 = 0;
209 T acc2 = 0;
210 T acc3 = 0;
211
212 T in_j = in(isize(j));
213
214 usize pcount = col_end - col_start;
215
216 auto zx = proxsuite::linalg::sparse::util::zero_extend;
217
218 if (zx(ai[col_end - 1]) == j) {
219 T ajj = ax[col_end - 1];
220 out(isize(j)) += ajj * in_j;
221 pcount -= 1;
222 }
223
224 usize p = col_start;
225
226 for (; p < col_start + pcount / 4 * 4; p += 4) {
227 auto i0 = isize(zx(ai[p + 0]));
228 auto i1 = isize(zx(ai[p + 1]));
229 auto i2 = isize(zx(ai[p + 2]));
230 auto i3 = isize(zx(ai[p + 3]));
231
232 T ai0j = ax[p + 0];
233 T ai1j = ax[p + 1];
234 T ai2j = ax[p + 2];
235 T ai3j = ax[p + 3];
236
237 out(i0) += ai0j * in_j;
238 out(i1) += ai1j * in_j;
239 out(i2) += ai2j * in_j;
240 out(i3) += ai3j * in_j;
241
242 acc0 += ai0j * in(i0);
243 acc1 += ai1j * in(i1);
244 acc2 += ai2j * in(i2);
245 acc3 += ai3j * in(i3);
246 }
247 for (; p < col_start + pcount; ++p) {
248 auto i = isize(zx(ai[p]));
249
250 T aij = ax[p];
251 out(i) += aij * in_j;
252 acc0 += aij * in(i);
253 }
254 acc0 = ((acc0 + acc1) + (acc2 + acc3));
255 out(isize(j)) += acc0;
256 }
257}
258
259template<typename OutL, typename OutR, typename A, typename InL, typename InR>
260void
262 OutR&& out_r,
263 A const& a,
264 InL const& in_l,
265 InR const& in_r)
266{
267 // noalias general vector matrix matrix vector add
268 noalias_gevmmv_add_impl<typename A::Scalar, typename A::StorageIndex>(
269 { proxqp::from_eigen, out_l },
270 { proxqp::from_eigen, out_r },
272 { proxqp::from_eigen, in_l },
273 { proxqp::from_eigen, in_r });
274}
275
276template<typename Out, typename A, typename In>
277void
278noalias_symhiv_add(Out&& out, A const& a, In const& in)
279{
280 // noalias symmetric (hi) matrix vector add
281 noalias_symhiv_add_impl<typename A::Scalar, typename A::StorageIndex>(
282 { proxqp::from_eigen, out },
284 { proxqp::from_eigen, in });
285}
286
287template<typename T, typename I>
288struct AugmentedKkt : Eigen::EigenBase<AugmentedKkt<T, I>>
289{
301
302 AugmentedKkt /* NOLINT */ (Raw raw) noexcept
303 : _{ raw }
304 {
305 }
306
307 using Scalar = T;
308 using RealScalar = T;
309 using StorageIndex = I;
310 enum
311 {
312 ColsAtCompileTime = Eigen::Dynamic,
313 MaxColsAtCompileTime = Eigen::Dynamic,
314 IsRowMajor = false,
315 };
316
317 auto rows() const noexcept -> isize { return _.n + _.n_eq + _.n_in; }
318 auto cols() const noexcept -> isize { return rows(); }
319 template<typename Rhs>
320 auto operator*(Eigen::MatrixBase<Rhs> const& x) const
321 -> Eigen::Product<AugmentedKkt, Rhs, Eigen::AliasFreeProduct>
322 {
323 return Eigen::Product< //
325 Rhs,
326 Eigen::AliasFreeProduct>{
327 *this,
328 x.derived(),
329 };
330 }
331};
332
333template<typename T>
334using VecMapMut = Eigen::Map<Eigen::Matrix<T, Eigen::Dynamic, 1>,
335 Eigen::Unaligned,
336 Eigen::Stride<Eigen::Dynamic, 1>>;
337template<typename T>
338using VecMap = Eigen::Map<Eigen::Matrix<T, Eigen::Dynamic, 1> const,
339 Eigen::Unaligned,
340 Eigen::Stride<Eigen::Dynamic, 1>>;
341
342template<typename V>
343auto
345{
346 static_assert(V::InnerStrideAtCompileTime == 1, ".");
347 return {
348 v.data(),
349 v.rows(),
350 v.cols(),
351 Eigen::Stride<Eigen::Dynamic, 1>{
352 v.outerStride(),
353 v.innerStride(),
354 },
355 };
356}
357
358template<typename V>
359auto
362{
363 static_assert(
364 proxsuite::linalg::veg::uncvref_t<V>::InnerStrideAtCompileTime == 1, ".");
365 return {
366 v.data(),
367 v.rows(),
368 v.cols(),
369 Eigen::Stride<Eigen::Dynamic, 1>{
370 v.outerStride(),
371 v.innerStride(),
372 },
373 };
374}
375
376template<typename T, typename I>
377auto
379 isize start,
380 isize ncols,
382{
383 VEG_ASSERT(start <= mat.ncols());
384 VEG_ASSERT(ncols <= mat.ncols() - start);
385
386 return {
387 proxsuite::linalg::sparse::from_raw_parts,
388 mat.nrows(),
389 ncols,
390 nnz,
391 mat.col_ptrs() + start,
392 mat.is_compressed() ? nullptr : (mat.nnz_per_col() + start),
393 mat.row_indices(),
394 mat.values(),
395 };
396}
397
398template<typename T, typename I>
399auto
401 isize start,
402 isize ncols,
404{
405 VEG_ASSERT(start <= mat.ncols());
406 VEG_ASSERT(ncols <= mat.ncols() - start);
407 return {
408 proxsuite::linalg::sparse::from_raw_parts,
409 mat.nrows(),
410 ncols,
411 nnz,
412 mat.col_ptrs_mut() + start,
413 mat.is_compressed() ? nullptr : (mat.nnz_per_col_mut() + start),
414 mat.row_indices_mut(),
415 mat.values_mut(),
416 };
417}
418
419template<typename T, typename I>
420auto
421top_rows_unchecked(proxsuite::linalg::veg::Unsafe /*unsafe*/,
424{
425 VEG_ASSERT(nrows <= mat.nrows());
426 return {
427 proxsuite::linalg::sparse::from_raw_parts,
428 nrows,
429 mat.ncols(),
430 mat.nnz(),
431 mat.col_ptrs(),
432 mat.nnz_per_col(),
433 mat.row_indices(),
434 mat.values(),
435 };
436}
437
438template<typename T, typename I>
439auto
440top_rows_mut_unchecked(proxsuite::linalg::veg::Unsafe /*unsafe*/,
443{
444 VEG_ASSERT(nrows <= mat.nrows());
445 return {
446 proxsuite::linalg::sparse::from_raw_parts,
447 nrows,
448 mat.ncols(),
449 mat.nnz(),
450 mat.col_ptrs_mut(),
451 mat.nnz_per_col_mut(),
452 mat.row_indices_mut(),
453 mat.values_mut(),
454 };
455}
471template<typename T, typename I, typename P>
472bool
474 VectorViewMut<T> CTdz,
477 const QpView<T, I> qp_scaled,
478 const Settings<T>& qpsettings,
479 const P& ruiz)
480{
481
482 // The problem is primal infeasible if the following four conditions hold:
483 //
484 // ||unscaled(A^Tdy)|| <= eps_p_inf ||unscaled(dy)||
485 // b^T dy <= -eps_p_inf ||unscaled(dy)||
486 // ||unscaled(C^Tdz)|| <= eps_p_inf ||unscaled(dz)||
487 // u^T [dz]_+ - l^T[-dz]_+ <= -eps_p_inf ||unscaled(dz)||
488 //
489 // the variables in entry are changed in place
490
491 bool res = infty_norm(dy.to_eigen()) != 0 || infty_norm(dz.to_eigen()) != 0;
492 if (!res) {
493 return res;
494 }
495 ruiz.unscale_dual_residual_in_place(ATdy);
496 ruiz.unscale_dual_residual_in_place(CTdz);
497 T lower_bound_1 =
498 dy.to_eigen().dot(qp_scaled.b.to_eigen()) +
499 helpers::positive_part(dz.to_eigen()).dot(qp_scaled.u.to_eigen()) -
500 helpers::negative_part(dz.to_eigen()).dot(qp_scaled.l.to_eigen());
501 ruiz.unscale_dual_in_place_eq(dy);
502 ruiz.unscale_dual_in_place_in(dz);
503 T lower_bound_2 = infty_norm(ATdy.to_eigen() + CTdz.to_eigen());
504 T upper_bound =
505 qpsettings.eps_primal_inf *
506 std::max(infty_norm(dy.to_eigen()), infty_norm(dz.to_eigen()));
507 res = lower_bound_2 <= upper_bound && lower_bound_1 <= -upper_bound;
508 return res;
509}
527template<typename T, typename I, typename P>
528bool
533 const QpView<T, I> qp_scaled,
534 const Settings<T>& qpsettings,
535 const Model<T, I>& qpmodel,
536 const P& ruiz)
537{
538
539 // The problem is dual infeasible if one of the conditions hold:
540 //
541 // FIRST
542 // ||unscaled(Adx)|| <= eps_d_inf ||unscaled(dx)||
543 // unscaled(Cdx)_i \in [-eps_d_inf,eps_d_inf] ||unscaled(dx)|| if u_i and l_i
544 // are finite or >=
545 // -eps_d_inf||unscaled(dx)|| if u_i = +inf or <= eps_d_inf||unscaled(dx)|| if
546 // l_i = -inf
547 //
548 // SECOND
549 //
550 // ||unscaled(Hdx)|| <= c eps_d_inf * ||unscaled(dx)|| and q^Tdx <= -c
551 // eps_d_inf ||unscaled(dx)||
552 // the variables in entry are changed in place
553 ruiz.unscale_dual_residual_in_place(Hdx);
554 ruiz.unscale_primal_residual_in_place_eq(Adx);
555 ruiz.unscale_primal_residual_in_place_in(Cdx);
556 T gdx = (dx.to_eigen()).dot(qp_scaled.g.to_eigen());
557 ruiz.unscale_primal_in_place(dx);
558
559 T bound = infty_norm(dx.to_eigen()) * qpsettings.eps_dual_inf;
560 T bound_neg = -bound;
561
562 bool first_cond = infty_norm(Adx.to_eigen()) <= bound;
563
564 for (i64 iter = 0; iter < qpmodel.n_in; ++iter) {
565 T Cdx_i = Cdx.to_eigen()[iter];
566 if (qp_scaled.u.to_eigen()[iter] <= T(1.E20) &&
567 qp_scaled.l.to_eigen()[iter] >= T(-1.E20)) {
568 first_cond = first_cond && Cdx_i <= bound && Cdx_i >= bound_neg;
569 } else if (qp_scaled.u.to_eigen()[iter] > T(1.E20)) {
570 first_cond = first_cond && Cdx_i >= bound_neg;
571 } else if (qp_scaled.l.to_eigen()[iter] < T(-1.E20)) {
572 first_cond = first_cond && Cdx_i <= bound;
573 }
574 }
575
576 bound *= ruiz.c;
577 bound_neg *= ruiz.c;
578 bool second_cond_alt1 =
579 infty_norm(Hdx.to_eigen()) <= bound && gdx <= bound_neg;
580 bound_neg *= qpsettings.eps_dual_inf;
581
582 bool res = first_cond && second_cond_alt1 && infty_norm(dx.to_eigen()) != 0;
583 return res;
584}
585
614template<typename T, typename I, typename P>
615auto
617 Results<T>& results,
618 const Settings<T>& settings,
619 VecMapMut<T> primal_residual_eq_scaled,
620 VecMapMut<T> primal_residual_in_scaled_lo,
621 VecMapMut<T> primal_residual_in_scaled_up,
622 VecMapMut<T> dual_residual_scaled,
623 T& primal_feasibility_eq_rhs_0,
624 T& primal_feasibility_in_rhs_0,
625 T& dual_feasibility_rhs_0,
626 T& dual_feasibility_rhs_1,
627 T& dual_feasibility_rhs_3,
628 T& rhs_duality_gap,
629 const P& precond,
630 Model<T, I> const& data,
631 const QpView<T, I> qp_scaled,
632 VecMapMut<T> x_e,
633 VecMapMut<T> y_e,
634 VecMapMut<T> z_e,
637{
638 isize n = x_e.rows();
639
640 LDLT_TEMP_VEC_UNINIT(T, tmp, n, stack);
641 dual_residual_scaled = qp_scaled.g.to_eigen();
642 {
643 tmp.setZero();
644 noalias_symhiv_add(tmp, qp_scaled.H.to_eigen(), x_e);
645 dual_residual_scaled += tmp;
646
647 precond.unscale_dual_residual_in_place(
648 { proxqp::from_eigen, tmp }); // contains unscaled Hx
649 dual_feasibility_rhs_0 = infty_norm(tmp);
650 precond.unscale_primal_in_place({ proxqp::from_eigen, x_e });
651 results.info.duality_gap = x_e.dot(data.g); // contains gTx
652 rhs_duality_gap = std::fabs(results.info.duality_gap);
653
654 const T xHx = (tmp).dot(x_e);
655 results.info.duality_gap += xHx;
656 rhs_duality_gap = std::max(rhs_duality_gap, std::abs(xHx));
657 tmp += data.g; // contains now Hx+g
658 precond.scale_primal_in_place({ proxqp::from_eigen, x_e });
659
660 precond.unscale_dual_in_place_eq({ proxsuite::proxqp::from_eigen, y_e });
661 const T by = (data.b).dot(y_e);
662 results.info.duality_gap += by;
663 rhs_duality_gap = std::max(rhs_duality_gap, std::abs(by));
664 precond.scale_dual_in_place_eq({ proxsuite::proxqp::from_eigen, y_e });
665
666 precond.unscale_dual_in_place_in({ proxsuite::proxqp::from_eigen, z_e });
667
668 const T zl =
669 helpers::select(work.active_set_low, results.z, 0)
671 results.info.duality_gap += zl;
672 rhs_duality_gap = std::max(rhs_duality_gap, std::abs(zl));
673
674 const T zu =
675 helpers::select(work.active_set_up, results.z, 0)
677 results.info.duality_gap += zu;
678 rhs_duality_gap = std::max(rhs_duality_gap, std::abs(zu));
679
680 precond.scale_dual_in_place_in({ proxsuite::proxqp::from_eigen, z_e });
681 }
682
683 {
684 auto ATy = tmp;
685 ATy.setZero();
686 primal_residual_eq_scaled.setZero();
687
689 primal_residual_eq_scaled, ATy, qp_scaled.AT.to_eigen(), x_e, y_e);
690
691 dual_residual_scaled += ATy;
692
693 precond.unscale_dual_residual_in_place({ proxqp::from_eigen, ATy });
694 dual_feasibility_rhs_1 = infty_norm(ATy);
695 }
696
697 {
698 auto CTz = tmp;
699 CTz.setZero();
700 primal_residual_in_scaled_up.setZero();
701
703 primal_residual_in_scaled_up, CTz, qp_scaled.CT.to_eigen(), x_e, z_e);
704
705 dual_residual_scaled += CTz;
706
707 precond.unscale_dual_residual_in_place({ proxqp::from_eigen, CTz });
708 dual_feasibility_rhs_3 = infty_norm(CTz);
709 }
710 precond.unscale_primal_residual_in_place_eq(
711 { proxqp::from_eigen, primal_residual_eq_scaled });
712 primal_feasibility_eq_rhs_0 = infty_norm(primal_residual_eq_scaled);
713
714 precond.unscale_primal_residual_in_place_in(
715 { proxqp::from_eigen, primal_residual_in_scaled_up });
716 primal_feasibility_in_rhs_0 = infty_norm(primal_residual_in_scaled_up);
717
718 auto b = data.b;
719 auto l = data.l;
720 auto u = data.u;
721 primal_residual_in_scaled_lo =
722 helpers::positive_part(primal_residual_in_scaled_up - u) +
723 helpers::negative_part(primal_residual_in_scaled_up - l);
724
725 primal_residual_eq_scaled -= b;
726 T primal_feasibility_eq_lhs = infty_norm(primal_residual_eq_scaled);
727 T primal_feasibility_in_lhs = infty_norm(primal_residual_in_scaled_lo);
728 T primal_feasibility_lhs =
729 std::max(primal_feasibility_eq_lhs, primal_feasibility_in_lhs);
730
731 if ((settings.primal_infeasibility_solving &&
732 results.info.status == QPSolverOutput::PROXQP_PRIMAL_INFEASIBLE)) {
733 tmp.setZero();
734 {
735 results.se = primal_residual_eq_scaled;
736 results.si = primal_residual_in_scaled_lo;
737 precond.unscale_primal_residual_in_place_eq(
738 { proxqp::from_eigen,
739 primal_residual_eq_scaled }); // E^{-1}(unscaled Ax-b)
740 tmp.noalias() = qp_scaled.AT.to_eigen() * primal_residual_eq_scaled;
741 }
742
743 {
744 precond.unscale_primal_residual_in_place_in(
745 { proxqp::from_eigen,
746 primal_residual_in_scaled_lo }); // E^{-1}(unscaled Ax-b)
747 tmp.noalias() += qp_scaled.CT.to_eigen() * primal_residual_in_scaled_lo;
748 }
749 precond.unscale_dual_residual_in_place({ proxqp::from_eigen, tmp });
750
751 primal_feasibility_lhs = infty_norm(tmp);
752 precond.scale_primal_residual_in_place_eq(
753 { proxqp::from_eigen, primal_residual_eq_scaled });
754 }
755
756 // scaled Ax - b
757 precond.scale_primal_residual_in_place_eq(
758 { proxqp::from_eigen, primal_residual_eq_scaled });
759 // scaled Cx
760 precond.scale_primal_residual_in_place_in(
761 { proxqp::from_eigen, primal_residual_in_scaled_up });
762
763 precond.unscale_dual_residual_in_place(
764 { proxqp::from_eigen, dual_residual_scaled });
765 T dual_feasibility_lhs = infty_norm(dual_residual_scaled);
766 precond.scale_dual_residual_in_place(
767 { proxqp::from_eigen, dual_residual_scaled });
768
769 return proxsuite::linalg::veg::tuplify(primal_feasibility_lhs,
770 dual_feasibility_lhs);
771}
772
773} // namespace detail
774} // namespace sparse
775} // namespace proxqp
776} // namespace proxsuite
777
778namespace Eigen {
779namespace internal {
780template<typename T, typename I>
781struct traits<proxsuite::proxqp::sparse::detail::AugmentedKkt<T, I>>
782 : Eigen::internal::traits<Eigen::SparseMatrix<T, Eigen::ColMajor, I>>
783{};
784
785template<typename Rhs, typename T, typename I>
786struct generic_product_impl<
787 proxsuite::proxqp::sparse::detail::AugmentedKkt<T, I>,
788 Rhs,
789 SparseShape,
790 DenseShape,
791 GemvProduct>
792 : generic_product_impl_base<
793 proxsuite::proxqp::sparse::detail::AugmentedKkt<T, I>,
794 Rhs,
795 generic_product_impl<
796 proxsuite::proxqp::sparse::detail::AugmentedKkt<T, I>,
797 Rhs>>
798{
800
801 using Scalar = typename Product<Mat_, Rhs>::Scalar;
802
803 template<typename Dst>
804 static void scaleAndAddTo(Dst& dst,
805 Mat_ const& lhs,
806 Rhs const& rhs,
807 PROXSUITE_MAYBE_UNUSED Scalar const& alpha)
808 {
810
811 VEG_ASSERT(alpha == Scalar(1));
813 dst, lhs._.kkt_active.to_eigen(), rhs);
814
815 {
816 isize n = lhs._.n;
817 isize n_eq = lhs._.n_eq;
818 isize n_in = lhs._.n_in;
819
820 auto dst_x = dst.head(n);
821 auto dst_y = dst.segment(n, n_eq);
822 auto dst_z = dst.tail(n_in);
823
824 auto rhs_x = rhs.head(n);
825 auto rhs_y = rhs.segment(n, n_eq);
826 auto rhs_z = rhs.tail(n_in);
827
828 dst_x += lhs._.rho * rhs_x;
829 dst_y += (-1 / lhs._.mu_eq) * rhs_y;
830 for (isize i = 0; i < n_in; ++i) {
831 dst_z[i] +=
832 (lhs._.active_constraints[i] ? -1 / lhs._.mu_in : T(1)) * rhs_z[i];
833 }
834 }
835 }
836};
837} // namespace internal
838} // namespace Eigen
839
840#endif /* end of include guard PROXSUITE_PROXQP_SPARSE_UTILS_HPP */
#define VEG_ASSERT_ALL_OF(...)
#define VEG_ASSERT(...)
#define LDLT_TEMP_VEC_UNINIT(Type, Name, Rows, Stack)
Definition core.hpp:76
#define PROXSUITE_MAYBE_UNUSED
Definition fwd.hpp:20
#define VEG_NO_INLINE
Definition macros.hpp:123
auto select(Condition const &condition, T const &expr, const Scalar value) PROXSUITE_DEDUCE_RET((condition).select(expr
Select the components of the expression if the condition is fullfiled. Otherwise, set the component t...
auto positive_part(T const &expr) PROXSUITE_DEDUCE_RET((expr.array() > 0).select(expr
Returns the positive part of an expression.
auto at_least(T const &expr, const Scalar value) PROXSUITE_DEDUCE_RET((expr.array() > value).select(expr
Returns the part of the expression which is greater than value.
auto at_most(T const &expr, const Scalar value) PROXSUITE_DEDUCE_RET((expr.array()< value).select(expr
Returns the part of the expression which is lower than value.
auto negative_part(T const &expr) PROXSUITE_DEDUCE_RET((expr.array()< 0).select(expr
Returns the negative part of an expression.
_detail::_meta::make_signed< usize >::Type isize
Definition typedefs.hpp:43
void noalias_gevmmv_add(OutL &&out_l, OutR &&out_r, A const &a, InL const &in_l, InR const &in_r)
Definition utils.hpp:261
VEG_NO_INLINE void noalias_symhiv_add_impl(VectorViewMut< T > out, proxsuite::linalg::sparse::MatRef< T, I > a, VectorView< T > in)
Definition utils.hpp:181
bool global_primal_residual_infeasibility(VectorViewMut< T > ATdy, VectorViewMut< T > CTdz, VectorViewMut< T > dy, VectorViewMut< T > dz, const QpView< T, I > qp_scaled, const Settings< T > &qpsettings, const P &ruiz)
Definition utils.hpp:473
auto middle_cols(proxsuite::linalg::sparse::MatRef< T, I > mat, isize start, isize ncols, isize nnz) -> proxsuite::linalg::sparse::MatRef< T, I >
Definition utils.hpp:378
void noalias_symhiv_add(Out &&out, A const &a, In const &in)
Definition utils.hpp:278
auto top_rows_unchecked(proxsuite::linalg::veg::Unsafe, proxsuite::linalg::sparse::MatRef< T, I > mat, isize nrows) -> proxsuite::linalg::sparse::MatRef< T, I >
Definition utils.hpp:421
VEG_NO_INLINE void noalias_gevmmv_add_impl(VectorViewMut< T > out_l, VectorViewMut< T > out_r, proxsuite::linalg::sparse::MatRef< T, I > a, VectorView< T > in_l, VectorView< T > in_r)
Definition utils.hpp:108
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
auto vec_mut(V &&v) -> VecMapMut< typename proxsuite::linalg::veg::uncvref_t< V >::Scalar >
Definition utils.hpp:360
auto vec(V const &v) -> VecMap< typename V::Scalar >
Definition utils.hpp:344
bool global_dual_residual_infeasibility(VectorViewMut< T > Adx, VectorViewMut< T > Cdx, VectorViewMut< T > Hdx, VectorViewMut< T > dx, const QpView< T, I > qp_scaled, const Settings< T > &qpsettings, const Model< T, I > &qpmodel, const P &ruiz)
Definition utils.hpp:529
Eigen::Map< Eigen::Matrix< T, Eigen::Dynamic, 1 >, Eigen::Unaligned, Eigen::Stride< Eigen::Dynamic, 1 > > VecMapMut
Definition utils.hpp:334
auto unscaled_primal_dual_residual(Workspace< T, I > &work, Results< T > &results, const Settings< T > &settings, VecMapMut< T > primal_residual_eq_scaled, VecMapMut< T > primal_residual_in_scaled_lo, VecMapMut< T > primal_residual_in_scaled_up, VecMapMut< T > dual_residual_scaled, T &primal_feasibility_eq_rhs_0, T &primal_feasibility_in_rhs_0, T &dual_feasibility_rhs_0, T &dual_feasibility_rhs_1, T &dual_feasibility_rhs_3, T &rhs_duality_gap, const P &precond, Model< T, I > const &data, const QpView< T, I > qp_scaled, VecMapMut< T > x_e, VecMapMut< T > y_e, VecMapMut< T > z_e, proxsuite::linalg::veg::dynstack::DynStackMut stack) -> proxsuite::linalg::veg::Tuple< T, T >
Definition utils.hpp:616
Eigen::Map< Eigen::Matrix< T, Eigen::Dynamic, 1 > const, Eigen::Unaligned, Eigen::Stride< Eigen::Dynamic, 1 > > VecMap
Definition utils.hpp:338
void print_setup_header(const Settings< T > &settings, Results< T > &results, const Model< T, I > &model)
Definition utils.hpp:36
decltype(sizeof(0)) usize
Definition views.hpp:24
void print_preambule()
Definition prints.hpp:30
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
InitialGuessStatus initial_guess
Definition settings.hpp:123
VEG_INLINE auto to_eigen() const -> detail::VecMapMut< T >
Definition views.hpp:677
This class stores the model of the QP problem.
Definition model.hpp:23
proxsuite::linalg::sparse::DenseVecRef< T > u
Definition views.hpp:36
proxsuite::linalg::sparse::DenseVecRef< T > b
Definition views.hpp:33
proxsuite::linalg::sparse::DenseVecRef< T > g
Definition views.hpp:31
proxsuite::linalg::sparse::DenseVecRef< T > l
Definition views.hpp:35
This class defines the workspace of the sparse solver.
proxsuite::linalg::sparse::MatRef< T, I > kkt_active
Definition utils.hpp:292
proxsuite::linalg::veg::Slice< bool > active_constraints
Definition utils.hpp:293
auto cols() const noexcept -> isize
Definition utils.hpp:318
auto operator*(Eigen::MatrixBase< Rhs > const &x) const -> Eigen::Product< AugmentedKkt, Rhs, Eigen::AliasFreeProduct >
Definition utils.hpp:320
auto rows() const noexcept -> isize
Definition utils.hpp:317
struct proxsuite::proxqp::sparse::detail::AugmentedKkt::Raw _