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//
7#ifndef PROXSUITE_PROXQP_DENSE_UTILS_HPP
8#define PROXSUITE_PROXQP_DENSE_UTILS_HPP
9
10#include <iostream>
11#include <fstream>
12#include <cmath>
13#include <type_traits>
14
23
24// #include <fmt/format.h>
25// #include <fmt/ostream.h>
26
27namespace proxsuite {
28namespace proxqp {
29namespace dense {
30
31template<typename T>
32void
34 const Results<T>& results,
35 const Model<T>& model,
36 const bool box_constraints,
37 const DenseBackend& dense_backend,
38 const HessianType& hessian_type)
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 << std::endl;
49
50 // Print Settings
51 std::cout << "settings: " << std::endl;
52 std::cout << " backend = dense," << std::endl;
53 std::cout << " eps_abs = " << settings.eps_abs
54 << " eps_rel = " << settings.eps_rel << std::endl;
55 std::cout << " eps_prim_inf = " << settings.eps_primal_inf
56 << ", eps_dual_inf = " << settings.eps_dual_inf << "," << std::endl;
57
58 std::cout << " rho = " << results.info.rho
59 << ", mu_eq = " << results.info.mu_eq
60 << ", mu_in = " << results.info.mu_in << "," << std::endl;
61 std::cout << " max_iter = " << settings.max_iter
62 << ", max_iter_in = " << settings.max_iter_in << "," << std::endl;
63 if (box_constraints) {
64 std::cout << " box constraints: on, " << std::endl;
65 } else {
66 std::cout << " box constraints: off, " << std::endl;
67 }
68 switch (dense_backend) {
70 std::cout << " dense backend: PrimalDualLDLT, " << std::endl;
71 break;
73 std::cout << " dense backend: PrimalLDLT, " << std::endl;
74 break;
76 break;
77 }
78 switch (hessian_type) {
80 std::cout << " problem type: Quadratic Program, " << std::endl;
81 break;
83 std::cout << " problem type: Linear Program, " << std::endl;
84 break;
86 std::cout
87 << " problem type: Quadratic Program with diagonal Hessian, "
88 << std::endl;
89 break;
90 }
91 if (settings.compute_preconditioner) {
92 std::cout << " scaling: on, " << std::endl;
93 } else {
94 std::cout << " scaling: off, " << std::endl;
95 }
96 if (settings.compute_timings) {
97 std::cout << " timings: on, " << std::endl;
98 } else {
99 std::cout << " timings: off, " << std::endl;
100 }
101 switch (settings.initial_guess) {
103 std::cout << " initial guess: warm start. \n" << std::endl;
104 break;
106 std::cout << " initial guess: no initial guess. \n" << std::endl;
107 break;
109 std::cout
110 << " initial guess: warm start with previous result. \n"
111 << std::endl;
112 break;
114 std::cout
115 << " initial guess: cold start with previous result. \n"
116 << std::endl;
117 break;
119 std::cout
120 << " initial guess: equality constrained initial guess. \n"
121 << std::endl;
122 }
123}
124
131template<typename Derived>
132void
133save_data(const std::string& filename, const ::Eigen::MatrixBase<Derived>& mat)
134{
135 // https://eigen.tuxfamily.org/dox/structEigen_1_1IOFormat.html
136 const static Eigen::IOFormat CSVFormat(
137 Eigen::FullPrecision, Eigen::DontAlignCols, ", ", "\n");
138
139 std::ofstream file(filename);
140 if (file.is_open()) {
141 file << mat.format(CSVFormat);
142 file.close();
143 }
144}
145
164template<typename T>
165void
167 Results<T>& qpresults,
168 const Settings<T>& qpsettings,
169 Workspace<T>& qpwork,
171 const bool box_constraints,
172 T& primal_feasibility_lhs,
173 T& primal_feasibility_eq_rhs_0,
174 T& primal_feasibility_in_rhs_0,
175 T& primal_feasibility_eq_lhs,
176 T& primal_feasibility_in_lhs)
177{
178 // COMPUTES:
179 // primal_residual_eq_scaled = scaled(Ax - b)
180 //
181 // primal_feasibility_lhs = max(norm(unscaled(Ax - b)),
182 // norm(unscaled([Cx - u]+ + [Cx - l]-)))
183 // primal_feasibility_eq_rhs_0 = norm(unscaled(Ax))
184 // primal_feasibility_in_rhs_0 = norm(unscaled(Cx))
185 //
186 // MAY_ALIAS[primal_residual_in_scaled_u, primal_residual_in_scaled_l]
187 //
188 // INDETERMINATE:
189 // primal_residual_in_scaled_u = unscaled(Cx)
190 // primal_residual_in_scaled_l = unscaled([Cx - u]+ + [Cx - l]-)
191 // qpwork.primal_residual_eq_scaled.noalias() = qpwork.A_scaled * qpresults.x;
192 qpresults.se.noalias() = qpwork.A_scaled * qpresults.x;
193 qpwork.primal_residual_in_scaled_up.head(qpmodel.n_in).noalias() =
194 qpwork.C_scaled * qpresults.x;
195 if (box_constraints) {
196 qpwork.primal_residual_in_scaled_up.tail(qpmodel.dim) = qpresults.x;
197 // qpwork.primal_residual_in_scaled_up.tail(qpmodel.dim).array() *=
198 // qpwork.i_scaled.array();
200 from_eigen, qpwork.primal_residual_in_scaled_up.tail(qpmodel.dim) });
201 // ruiz.unscale_box_primal_residual_in_place_in(VectorViewMut<T>{from_eigen,
202 // qpwork.primal_residual_in_scaled_up.tail(qpmodel.dim)});
203 }
204 // ruiz.unscale_primal_residual_in_place_eq(
205 // VectorViewMut<T>{ from_eigen, qpwork.primal_residual_eq_scaled });
206 // primal_feasibility_eq_rhs_0 = infty_norm(qpwork.primal_residual_eq_scaled);
208 VectorViewMut<T>{ from_eigen, qpresults.se });
209 primal_feasibility_eq_rhs_0 = infty_norm(qpresults.se);
211 from_eigen, qpwork.primal_residual_in_scaled_up.head(qpmodel.n_in) });
212 primal_feasibility_in_rhs_0 =
213 infty_norm(qpwork.primal_residual_in_scaled_up.head(qpmodel.n_in));
214 qpresults.si.head(qpmodel.n_in) =
216 qpwork.primal_residual_in_scaled_up.head(qpmodel.n_in) - qpmodel.u) +
218 qpwork.primal_residual_in_scaled_up.head(qpmodel.n_in) - qpmodel.l);
219 if (box_constraints) {
220 qpresults.si.tail(qpmodel.dim) =
222 qpwork.primal_residual_in_scaled_up.tail(qpmodel.dim) - qpmodel.u_box) +
224 qpwork.primal_residual_in_scaled_up.tail(qpmodel.dim) - qpmodel.l_box);
225 qpwork.active_part_z.tail(qpmodel.dim) =
226 qpresults.x - qpresults.si.tail(qpmodel.dim);
227 primal_feasibility_in_rhs_0 =
228 std::max(primal_feasibility_in_rhs_0,
229 infty_norm(qpwork.active_part_z.tail(qpmodel.dim)));
230 primal_feasibility_in_rhs_0 =
231 std::max(primal_feasibility_in_rhs_0, infty_norm(qpresults.x));
232 }
233 // qpwork.primal_residual_eq_scaled -= qpmodel.b;
234 qpresults.se -= qpmodel.b;
235
236 primal_feasibility_in_lhs = infty_norm(qpresults.si);
237 // primal_feasibility_eq_lhs = infty_norm(qpwork.primal_residual_eq_scaled);
238 primal_feasibility_eq_lhs = infty_norm(qpresults.se);
239 primal_feasibility_lhs =
240 std::max(primal_feasibility_eq_lhs, primal_feasibility_in_lhs);
241 if (qpsettings.primal_infeasibility_solving &&
243 qpwork.rhs.head(qpmodel.dim).noalias() =
244 qpmodel.A.transpose() * qpresults.se;
245 qpwork.rhs.head(qpmodel.dim).noalias() +=
246 qpmodel.C.transpose() * qpresults.si;
247 primal_feasibility_lhs = infty_norm(qpwork.rhs.head(qpmodel.dim));
248 }
250 VectorViewMut<T>{ from_eigen, qpresults.se });
251 // VectorViewMut<T>{ from_eigen, qpwork.primal_residual_eq_scaled });
252}
253
269template<typename T>
270bool
272 VectorViewMut<T> ATdy,
273 VectorViewMut<T> CTdz,
276 Workspace<T>& qpwork,
277 const Model<T>& qpmodel,
278 const Settings<T>& qpsettings,
279 const bool box_constraints,
281{
282
283 // The problem is primal infeasible if the following four conditions hold:
284 //
285 // ||unscaled(A^Tdy)|| <= eps_p_inf ||unscaled(dy)||
286 // b^T dy <= -eps_p_inf ||unscaled(dy)||
287 // ||unscaled(C^Tdz)|| <= eps_p_inf ||unscaled(dz)||
288 // u^T [dz]_+ - l^T[-dz]_+ <= -eps_p_inf ||unscaled(dz)||
289 //
290 // the variables in entry are changed in place
291
292 bool res = infty_norm(dy.to_eigen()) != 0 || infty_norm(dz.to_eigen()) != 0;
293 if (!res) {
294 return res;
295 }
298 T lower_bound_1 = dy.to_eigen().dot(qpwork.b_scaled) +
299 helpers::positive_part(dz.to_eigen().head(qpmodel.n_in))
300 .dot(qpwork.u_scaled) -
301 helpers::negative_part(dz.to_eigen().head(qpmodel.n_in))
302 .dot(qpwork.l_scaled);
305 VectorViewMut<T>{ from_eigen, dz.to_eigen().head(qpmodel.n_in) });
306 if (box_constraints) {
307 // in_inf+=
308 // helpers::positive_part(dz.to_eigen().tail(qpmodel.dim)).dot(qpmodel.u) -
309 // helpers::negative_part(dz.to_eigen().tail(qpmodel.dim)).dot(qpmodel.l);
310 lower_bound_1 += helpers::positive_part(dz.to_eigen().tail(qpmodel.dim))
311 .dot(qpwork.u_box_scaled) -
312 helpers::negative_part(dz.to_eigen().tail(qpmodel.dim))
313 .dot(qpwork.l_box_scaled);
315 VectorViewMut<T>{ from_eigen, dz.to_eigen().tail(qpmodel.dim) });
316 }
317 T upper_bound =
318 qpsettings.eps_primal_inf *
319 std::max(infty_norm(dy.to_eigen()), infty_norm(dz.to_eigen()));
320 T lower_bound_2 = infty_norm(ATdy.to_eigen() + CTdz.to_eigen());
321
322 res = lower_bound_2 <= upper_bound && lower_bound_1 <= -upper_bound;
323 return res;
324}
325
343template<typename T>
344bool
350 Workspace<T>& qpwork,
351 const Settings<T>& qpsettings,
352 const Model<T>& qpmodel,
353 const bool box_constraints,
355{
356
357 // The problem is dual infeasible the two following conditions hold:
358 //
359 // FIRST
360 // ||unscaled(Adx)|| <= eps_d_inf ||unscaled(dx)||
361 // unscaled(Cdx)_i \in [-eps_d_inf,eps_d_inf] ||unscaled(dx)|| if u_i and l_i
362 // are finite or >=
363 // -eps_d_inf||unscaled(dx)|| if u_i = +inf or <= eps_d_inf||unscaled(dx)|| if
364 // l_i = -inf
365 //
366 // SECOND
367 //
368 // ||unscaled(Hdx)|| <= c eps_d_inf * ||unscaled(dx)|| and q^Tdx <= -c
369 // eps_d_inf ||unscaled(dx)|| the variables in
370 // entry are changed in place
374 VectorViewMut<T>{ from_eigen, Cdx.to_eigen().head(qpmodel.n_in) });
375 if (box_constraints) {
377 VectorViewMut<T>{ from_eigen, Cdx.to_eigen().tail(qpmodel.dim) });
378 }
379 T gdx = (dx.to_eigen()).dot(qpwork.g_scaled);
381
382 T bound = infty_norm(dx.to_eigen()) * qpsettings.eps_dual_inf;
383 T bound_neg = -bound;
384
385 bool first_cond = infty_norm(Adx.to_eigen()) <= bound;
386
387 for (i64 iter = 0; iter < qpmodel.n_in; ++iter) {
388 T Cdx_i = Cdx.to_eigen()[iter];
389 if (qpwork.u_scaled[iter] <= 1.E20 && qpwork.l_scaled[iter] >= -1.E20) {
390 first_cond = first_cond && Cdx_i <= bound && Cdx_i >= bound_neg;
391 } else if (qpwork.u_scaled[iter] > 1.E20) {
392 first_cond = first_cond && Cdx_i >= bound_neg;
393 } else if (qpwork.l_scaled[iter] < -1.E20) {
394 first_cond = first_cond && Cdx_i <= bound;
395 }
396 }
397 if (box_constraints) {
398 for (i64 iter = 0; iter < qpmodel.dim; ++iter) {
399 T dx_i = dx.to_eigen()[iter];
400 // if (qpmodel.u_box[iter] <= 1.E20 && qpmodel.l_box[iter] >= -1.E20) {
401 if (qpwork.u_box_scaled[iter] <= 1.E20 &&
402 qpwork.l_box_scaled[iter] >= -1.E20) {
403 first_cond = first_cond && dx_i <= bound && dx_i >= bound_neg;
404 } else if (qpwork.u_box_scaled[iter] > 1.E20) {
405 first_cond = first_cond && dx_i >= bound_neg;
406 } else if (qpwork.l_box_scaled[iter] < -1.E20) {
407 first_cond = first_cond && dx_i <= bound;
408 }
409 }
410 }
411 bound *= ruiz.c;
412 bound_neg *= ruiz.c;
413 bool second_cond_alt1 =
414 infty_norm(Hdx.to_eigen()) <= bound && gdx <= bound_neg;
415 bound_neg *= qpsettings.eps_dual_inf;
416
417 bool res = first_cond && second_cond_alt1 && infty_norm(dx.to_eigen()) != 0;
418 return res;
419}
420
437template<typename T>
438void
440 Workspace<T>& qpwork,
441 const Model<T>& qpmodel,
442 const bool box_constraints,
444 T& dual_feasibility_lhs,
445 T& dual_feasibility_rhs_0,
446 T& dual_feasibility_rhs_1,
447 T& dual_feasibility_rhs_3,
448 T& rhs_duality_gap,
449 T& duality_gap,
450 const HessianType& hessian_type)
451{
452 // dual_feasibility_lhs = norm(dual_residual_scaled)
453 // dual_feasibility_rhs_0 = norm(unscaled(Hx))
454 // dual_feasibility_rhs_1 = norm(unscaled(ATy))
455 // dual_feasibility_rhs_3 = norm(unscaled(CTz))
456 //
457 // dual_residual_scaled = scaled(Hx + g + ATy + CTz)
458
459 qpwork.dual_residual_scaled = qpwork.g_scaled;
460
461 switch (hessian_type) {
463 dual_feasibility_rhs_0 = 0;
464 break;
466 qpwork.CTz.noalias() =
467 qpwork.H_scaled.template selfadjointView<Eigen::Lower>() * qpresults.x;
468 qpwork.dual_residual_scaled += qpwork.CTz;
470 VectorViewMut<T>{ from_eigen, qpwork.CTz }); // contains unscaled Hx
471 dual_feasibility_rhs_0 = infty_norm(qpwork.CTz);
472 break;
474 qpwork.CTz.array() =
475 qpwork.H_scaled.diagonal().array() * qpresults.x.array();
476 qpwork.dual_residual_scaled += qpwork.CTz;
478 VectorViewMut<T>{ from_eigen, qpwork.CTz }); // contains unscaled Hx
479 dual_feasibility_rhs_0 = infty_norm(qpwork.CTz);
480 break;
481 }
482 ruiz.unscale_primal_in_place(VectorViewMut<T>{ from_eigen, qpresults.x });
483 duality_gap = (qpmodel.g).dot(qpresults.x);
484 rhs_duality_gap = std::fabs(duality_gap);
485 T xHx(0);
486 switch (hessian_type) {
488 break;
490 xHx = (qpwork.CTz).dot(qpresults.x);
491 duality_gap += xHx; // contains now xHx+g.Tx
492 rhs_duality_gap = std::max(rhs_duality_gap, std::abs(xHx));
493 break;
495 xHx = (qpwork.CTz).dot(qpresults.x);
496 duality_gap += xHx; // contains now xHx+g.Tx
497 rhs_duality_gap = std::max(rhs_duality_gap, std::abs(xHx));
498 break;
499 }
500 ruiz.scale_primal_in_place(VectorViewMut<T>{ from_eigen, qpresults.x });
501
502 qpwork.CTz.noalias() = qpwork.A_scaled.transpose() * qpresults.y;
503 qpwork.dual_residual_scaled += qpwork.CTz;
505 VectorViewMut<T>{ from_eigen, qpwork.CTz });
506 dual_feasibility_rhs_1 = infty_norm(qpwork.CTz);
507
508 qpwork.CTz.noalias() =
509 qpwork.C_scaled.transpose() * qpresults.z.head(qpmodel.n_in);
510 qpwork.dual_residual_scaled += qpwork.CTz;
512 VectorViewMut<T>{ from_eigen, qpwork.CTz });
513 dual_feasibility_rhs_3 = infty_norm(qpwork.CTz);
514 if (box_constraints) {
515 qpwork.CTz.noalias() = qpresults.z.tail(qpmodel.dim);
516 // FALSE : we should add I_scaled * z_scaled and I_scaled = D, so we
517 // unscaled z_scaled to unscale then the dual residual
518 // ruiz.unscale_primal_in_place(VectorViewMut<T>{from_eigen,qpwork.CTz});
519 qpwork.CTz.array() *= qpwork.i_scaled.array();
520
521 qpwork.dual_residual_scaled += qpwork.CTz;
523 VectorViewMut<T>{ from_eigen, qpwork.CTz });
524 dual_feasibility_rhs_3 =
525 std::max(infty_norm(qpwork.CTz), dual_feasibility_rhs_3);
526 }
527
529 VectorViewMut<T>{ from_eigen, qpwork.dual_residual_scaled });
530
531 dual_feasibility_lhs = infty_norm(qpwork.dual_residual_scaled);
532
534 VectorViewMut<T>{ from_eigen, qpwork.dual_residual_scaled });
535
536 ruiz.unscale_dual_in_place_eq(VectorViewMut<T>{ from_eigen, qpresults.y });
537 const T by = (qpmodel.b).dot(qpresults.y);
538 rhs_duality_gap = std::max(rhs_duality_gap, std::abs(by));
539 duality_gap += by;
540 ruiz.scale_dual_in_place_eq(VectorViewMut<T>{ from_eigen, qpresults.y });
541
543 VectorViewMut<T>{ from_eigen, qpresults.z.head(qpmodel.n_in) });
544
545 T zu =
546 helpers::select(qpwork.active_set_up.head(qpmodel.n_in),
547 qpresults.z.head(qpmodel.n_in),
548 0)
550 rhs_duality_gap = std::max(rhs_duality_gap, std::abs(zu));
551 duality_gap += zu;
552
553 T zl =
554 helpers::select(qpwork.active_set_low.head(qpmodel.n_in),
555 qpresults.z.head(qpmodel.n_in),
556 0)
558 rhs_duality_gap = std::max(rhs_duality_gap, std::abs(zl));
559 duality_gap += zl;
560
562 VectorViewMut<T>{ from_eigen, qpresults.z.head(qpmodel.n_in) });
563
564 if (box_constraints) {
566 VectorViewMut<T>{ from_eigen, qpresults.z.tail(qpmodel.dim) });
567
568 zu = helpers::select(qpwork.active_set_up.tail(qpmodel.dim),
569 qpresults.z.tail(qpmodel.dim),
570 0)
571 .dot(helpers::at_most(qpmodel.u_box,
573 rhs_duality_gap = std::max(rhs_duality_gap, std::abs(zu));
574 duality_gap += zu;
575
576 zl = helpers::select(qpwork.active_set_low.tail(qpmodel.dim),
577 qpresults.z.tail(qpmodel.dim),
578 0)
579 .dot(helpers::at_least(qpmodel.l_box,
581 rhs_duality_gap = std::max(rhs_duality_gap, std::abs(zl));
582 duality_gap += zl;
583
585 VectorViewMut<T>{ from_eigen, qpresults.z.tail(qpmodel.dim) });
586 }
587}
588
589} // namespace dense
590} // namespace proxqp
591} // namespace proxsuite
592
593#endif /* end of include guard PROXSUITE_PROXQP_DENSE_UTILS_HPP */
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.
void global_dual_residual(Results< T > &qpresults, Workspace< T > &qpwork, const Model< T > &qpmodel, const bool box_constraints, const preconditioner::RuizEquilibration< T > &ruiz, T &dual_feasibility_lhs, T &dual_feasibility_rhs_0, T &dual_feasibility_rhs_1, T &dual_feasibility_rhs_3, T &rhs_duality_gap, T &duality_gap, const HessianType &hessian_type)
Definition utils.hpp:439
bool global_dual_residual_infeasibility(VectorViewMut< T > Adx, VectorViewMut< T > Cdx, VectorViewMut< T > Hdx, VectorViewMut< T > dx, Workspace< T > &qpwork, const Settings< T > &qpsettings, const Model< T > &qpmodel, const bool box_constraints, const preconditioner::RuizEquilibration< T > &ruiz)
Definition utils.hpp:345
bool global_primal_residual_infeasibility(VectorViewMut< T > ATdy, VectorViewMut< T > CTdz, VectorViewMut< T > dy, VectorViewMut< T > dz, Workspace< T > &qpwork, const Model< T > &qpmodel, const Settings< T > &qpsettings, const bool box_constraints, const preconditioner::RuizEquilibration< T > &ruiz)
Definition utils.hpp:271
void global_primal_residual(const Model< T > &qpmodel, Results< T > &qpresults, const Settings< T > &qpsettings, Workspace< T > &qpwork, const preconditioner::RuizEquilibration< T > &ruiz, const bool box_constraints, T &primal_feasibility_lhs, T &primal_feasibility_eq_rhs_0, T &primal_feasibility_in_rhs_0, T &primal_feasibility_eq_lhs, T &primal_feasibility_in_lhs)
Definition utils.hpp:166
void print_setup_header(const Settings< T > &settings, const Results< T > &results, const Model< T > &model, const bool box_constraints, const DenseBackend &dense_backend, const HessianType &hessian_type)
Definition utils.hpp:33
void save_data(const std::string &filename, const ::Eigen::MatrixBase< Derived > &mat)
Definition utils.hpp:133
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
sparse::Vec< T > z
Definition results.hpp:74
sparse::Vec< T > si
Definition results.hpp:77
sparse::Vec< T > y
Definition results.hpp:73
sparse::Vec< T > se
Definition results.hpp:75
sparse::Vec< T > x
Definition results.hpp:72
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:24
This class defines the workspace of the dense solver.
Definition workspace.hpp:26
void unscale_primal_in_place(VectorViewMut< T > primal) const
Definition ruiz.hpp:554
void scale_primal_residual_in_place_eq(VectorViewMut< T > primal_eq) const
Definition ruiz.hpp:618
void scale_box_dual_in_place_in(VectorViewMut< T > dual) const
Definition ruiz.hpp:580
void unscale_dual_in_place_in(VectorViewMut< T > dual) const
Definition ruiz.hpp:598
void unscale_primal_residual_in_place_in(VectorViewMut< T > primal_in) const
Definition ruiz.hpp:675
void scale_primal_in_place(VectorViewMut< T > primal) const
Definition ruiz.hpp:518
void unscale_dual_in_place_eq(VectorViewMut< T > dual) const
Definition ruiz.hpp:589
void unscale_primal_residual_in_place_eq(VectorViewMut< T > primal_eq) const
Definition ruiz.hpp:667
void unscale_dual_residual_in_place(VectorViewMut< T > dual) const
Definition ruiz.hpp:691
void scale_dual_residual_in_place(VectorViewMut< T > dual) const
Definition ruiz.hpp:642
void scale_dual_in_place_in(VectorViewMut< T > dual) const
Definition ruiz.hpp:545
void scale_dual_in_place_eq(VectorViewMut< T > dual) const
Definition ruiz.hpp:536
void unscale_box_dual_in_place_in(VectorViewMut< T > dual) const
Definition ruiz.hpp:571
void unscale_box_primal_residual_in_place_in(VectorViewMut< T > primal_in) const
Definition ruiz.hpp:683