WPILibC++ 2027.0.0-alpha-5
Loading...
Searching...
No Matches
feasibility_restoration.hpp
Go to the documentation of this file.
1// Copyright (c) Sleipnir contributors
2
3#pragma once
4
5#include <algorithm>
6#include <cmath>
7#include <span>
8#include <tuple>
9#include <utility>
10
11#include <Eigen/Core>
12#include <gch/small_vector.hpp>
13
21
22namespace slp {
23
24template <typename Scalar>
26 const InteriorPointMatrixCallbacks<Scalar>& matrix_callbacks,
27 std::span<std::function<bool(const IterationInfo<Scalar>& info)>>
28 iteration_callbacks,
29 const Options& options, bool in_feasibility_restoration,
30#ifdef SLEIPNIR_ENABLE_BOUND_PROJECTION
31 const Eigen::ArrayX<bool>& bound_constraint_mask,
32#endif
33 Eigen::Vector<Scalar, Eigen::Dynamic>& x,
34 Eigen::Vector<Scalar, Eigen::Dynamic>& s,
35 Eigen::Vector<Scalar, Eigen::Dynamic>& y,
36 Eigen::Vector<Scalar, Eigen::Dynamic>& z, Scalar& μ);
37
38/// Computes initial values for p and n in feasibility restoration.
39///
40/// @tparam Scalar Scalar type.
41/// @param[in] c The constraint violation.
42/// @param[in] ρ Scaling parameter.
43/// @param[in] μ Barrier parameter.
44/// @return Tuple of positive and negative constraint relaxation slack variables
45/// respectively.
46template <typename Scalar>
47std::tuple<Eigen::Vector<Scalar, Eigen::Dynamic>,
48 Eigen::Vector<Scalar, Eigen::Dynamic>>
49compute_p_n(const Eigen::Vector<Scalar, Eigen::Dynamic>& c, Scalar ρ,
50 Scalar μ) {
51 // From equation (33) of [2]:
52 // ______________________
53 // μ − ρ c(x) /(μ − ρ c(x))² μ c(x)
54 // n = −−−−−−−−−− + / (−−−−−−−−−−) + −−−−−− (1)
55 // 2ρ √ ( 2ρ ) 2ρ
56 //
57 // The quadratic formula:
58 // ________
59 // -b + √b² - 4ac
60 // x = −−−−−−−−−−−−−− (2)
61 // 2a
62 //
63 // Rearrange (1) to fit the quadratic formula better:
64 // _________________________
65 // μ - ρ c(x) + √(μ - ρ c(x))² + 2ρ μ c(x)
66 // n = −−−−−−−−−−−−−−−−−−−−−−−−−−−−−−−−−−−−−−−
67 // 2ρ
68 //
69 // Solve for coefficients:
70 //
71 // a = ρ (3)
72 // b = ρ c(x) - μ (4)
73 //
74 // -4ac = 2ρ μ c(x)
75 // -4(ρ)c = 2ρ μ c(x)
76 // -4c = 2μ c(x)
77 // c = -μ c(x)/2 (5)
78 //
79 // p = c(x) + n (6)
80
81 using DenseVector = Eigen::Vector<Scalar, Eigen::Dynamic>;
82
83 using std::sqrt;
84
85 DenseVector p{c.rows()};
86 DenseVector n{c.rows()};
87 for (int row = 0; row < p.rows(); ++row) {
88 Scalar _a = ρ;
89 Scalar _b = ρ * c[row] - μ;
90 Scalar _c = -μ * c[row] / Scalar(2);
91
92 n[row] = (-_b + sqrt(_b * _b - Scalar(4) * _a * _c)) / (Scalar(2) * _a);
93 p[row] = c[row] + n[row];
94 }
95
96 return {std::move(p), std::move(n)};
97}
98
99// @cond Suppress Doxygen
100
101/// Finds the iterate that minimizes the constraint violation while not
102/// deviating too far from the starting point. This is a fallback procedure when
103/// the normal Sequential Quadratic Programming method fails to converge to a
104/// feasible point.
105///
106/// @tparam Scalar Scalar type.
107/// @param[in] matrix_callbacks Matrix callbacks.
108/// @param[in] iteration_callbacks The list of callbacks to call at the
109/// beginning of each iteration.
110/// @param[in] options Solver options.
111/// @param[in,out] x The decision variables from the normal solve.
112/// @param[in,out] y The equality constraint dual variables from the normal
113/// solve.
114/// @return The exit status.
115template <typename Scalar>
116ExitStatus feasibility_restoration(
117 const SQPMatrixCallbacks<Scalar>& matrix_callbacks,
118 std::span<std::function<bool(const IterationInfo<Scalar>& info)>>
119 iteration_callbacks,
120 const Options& options, Eigen::Vector<Scalar, Eigen::Dynamic>& x,
121 Eigen::Vector<Scalar, Eigen::Dynamic>& y) {
122 // Feasibility restoration
123 //
124 // min ρ Σ (pₑ + nₑ) + ζ/2 (x - xᵣ)ᵀDᵣ(x - xᵣ)
125 // x
126 // pₑ,nₑ
127 //
128 // s.t. cₑ(x) - pₑ + nₑ = 0
129 // pₑ ≥ 0
130 // nₑ ≥ 0
131 //
132 // where ρ = 1000, ζ = √μ where μ is the barrier parameter, xᵣ is original
133 // iterate before feasibility restoration, and Dᵣ is a scaling matrix defined
134 // by
135 //
136 // Dᵣ = diag(min(1, 1/xᵣ[i]²) for i in x.rows())
137
138 using DenseVector = Eigen::Vector<Scalar, Eigen::Dynamic>;
139 using DiagonalMatrix = Eigen::DiagonalMatrix<Scalar, Eigen::Dynamic>;
140 using SparseMatrix = Eigen::SparseMatrix<Scalar>;
141 using SparseVector = Eigen::SparseVector<Scalar>;
142
143 using std::sqrt;
144
145 const auto& matrices = matrix_callbacks;
146 const auto& num_vars = matrices.num_decision_variables;
147 const auto& num_eq = matrices.num_equality_constraints;
148
149 constexpr Scalar ρ(1e3);
150 const Scalar μ(options.tolerance / 10.0);
151 const Scalar ζ = sqrt(μ);
152
153 const DenseVector c_e = matrices.c_e(x);
154
155 const auto& x_r = x;
156 const auto [p_e_0, n_e_0] = compute_p_n(c_e, ρ, μ);
157
158 // Dᵣ = diag(min(1, 1/xᵣ[i]²) for i in x.rows())
159 const DiagonalMatrix D_r =
160 x.cwiseSquare().cwiseInverse().cwiseMin(Scalar(1)).asDiagonal();
161
162 DenseVector fr_x{num_vars + 2 * num_eq};
163 fr_x << x, p_e_0, n_e_0;
164
165 DenseVector fr_s = DenseVector::Ones(2 * num_eq);
166
167 DenseVector fr_y = DenseVector::Zero(num_eq);
168
169 DenseVector fr_z{2 * num_eq};
170 fr_z << μ * p_e_0.cwiseInverse(), μ * n_e_0.cwiseInverse();
171
172 Scalar fr_μ = std::max(μ, c_e.template lpNorm<Eigen::Infinity>());
173
174 InteriorPointMatrixCallbacks<Scalar> fr_matrix_callbacks{
175 static_cast<int>(fr_x.rows()),
176 static_cast<int>(fr_y.rows()),
177 static_cast<int>(fr_z.rows()),
178 [&](const DenseVector& x_p) -> Scalar {
179 auto x = x_p.segment(0, num_vars);
180
181 // Cost function
182 //
183 // ρ Σ (pₑ + nₑ) + ζ/2 (x - xᵣ)ᵀDᵣ(x - xᵣ)
184
185 auto diff = x - x_r;
186 return ρ * x_p.segment(num_vars, 2 * num_eq).array().sum() +
187 ζ / Scalar(2) * diff.transpose() * D_r * diff;
188 },
189 [&](const DenseVector& x_p) -> SparseVector {
190 auto x = x_p.segment(0, num_vars);
191
192 // Cost function gradient
193 //
194 // [ζDᵣ(x − xᵣ)]
195 // [ ρ ]
196 // [ ρ ]
197 DenseVector g{x_p.rows()};
198 g.segment(0, num_vars) = ζ * D_r * (x - x_r);
199 g.segment(num_vars, 2 * num_eq).setConstant(ρ);
200 return g.sparseView();
201 },
202 [&](const DenseVector& x_p, const DenseVector& y_p,
203 [[maybe_unused]] const DenseVector& z_p) -> SparseMatrix {
204 auto x = x_p.segment(0, num_vars);
205 const auto& y = y_p;
206
207 // Cost function Hessian
208 //
209 // [ζDᵣ 0 0]
210 // [ 0 0 0]
211 // [ 0 0 0]
213 triplets.reserve(x_p.rows());
214 append_as_triplets(triplets, 0, 0, {SparseMatrix{ζ * D_r}});
215 SparseMatrix d2f_dx2{x_p.rows(), x_p.rows()};
216 d2f_dx2.setFromSortedTriplets(triplets.begin(), triplets.end());
217
218 // Constraint part of original problem's Lagrangian Hessian
219 //
220 // −∇ₓₓ²yᵀcₑ(x)
221 auto H_c = matrices.H_c(x, y);
222 H_c.resize(x_p.rows(), x_p.rows());
223
224 // Lagrangian Hessian
225 //
226 // [ζDᵣ 0 0]
227 // [ 0 0 0] − ∇ₓₓ²yᵀcₑ(x)
228 // [ 0 0 0]
229 return d2f_dx2 + H_c;
230 },
231 [&](const DenseVector& x_p, [[maybe_unused]] const DenseVector& y_p,
232 [[maybe_unused]] const DenseVector& z_p) -> SparseMatrix {
233 return SparseMatrix{x_p.rows(), x_p.rows()};
234 },
235 [&](const DenseVector& x_p) -> DenseVector {
236 auto x = x_p.segment(0, num_vars);
237 auto p_e = x_p.segment(num_vars, num_eq);
238 auto n_e = x_p.segment(num_vars + num_eq, num_eq);
239
240 // Equality constraints
241 //
242 // cₑ(x) - pₑ + nₑ = 0
243 return matrices.c_e(x) - p_e + n_e;
244 },
245 [&](const DenseVector& x_p) -> SparseMatrix {
246 auto x = x_p.segment(0, num_vars);
247
248 // Equality constraint Jacobian
249 //
250 // [Aₑ −I I]
251
252 SparseMatrix A_e = matrices.A_e(x);
253
255 triplets.reserve(A_e.nonZeros() + 2 * num_eq);
256
257 append_as_triplets(triplets, 0, 0, {A_e});
259 triplets, 0, num_vars,
260 DenseVector::Constant(num_eq, Scalar(-1)).eval());
262 triplets, 0, num_vars + num_eq,
263 DenseVector::Constant(num_eq, Scalar(1)).eval());
264
265 SparseMatrix A_e_p{A_e.rows(), x_p.rows()};
266 A_e_p.setFromSortedTriplets(triplets.begin(), triplets.end());
267 return A_e_p;
268 },
269 [&](const DenseVector& x_p) -> DenseVector {
270 // Inequality constraints
271 //
272 // pₑ ≥ 0
273 // nₑ ≥ 0
274 return x_p.segment(num_vars, 2 * num_eq);
275 },
276 [&](const DenseVector& x_p) -> SparseMatrix {
277 // Inequality constraint Jacobian
278 //
279 // [0 I 0]
280 // [0 0 I]
281
283 triplets.reserve(2 * num_eq);
284
286 triplets, 0, num_vars,
287 DenseVector::Constant(2 * num_eq, Scalar(1)).eval());
288
289 SparseMatrix A_i_p{2 * num_eq, x_p.rows()};
290 A_i_p.setFromSortedTriplets(triplets.begin(), triplets.end());
291 return A_i_p;
292 }};
293
294 auto status = interior_point<Scalar>(fr_matrix_callbacks, iteration_callbacks,
295 options, true,
296#ifdef SLEIPNIR_ENABLE_BOUND_PROJECTION
297 {},
298#endif
299 fr_x, fr_s, fr_y, fr_z, fr_μ);
300
301 x = fr_x.segment(0, x.rows());
302
304 auto g = matrices.g(x);
305 auto A_e = matrices.A_e(x);
306
308
309 return ExitStatus::SUCCESS;
310 } else if (status == ExitStatus::SUCCESS) {
312 } else {
314 }
315}
316
317/// Finds the iterate that minimizes the constraint violation while not
318/// deviating too far from the starting point. This is a fallback procedure when
319/// the normal interior-point method fails to converge to a feasible point.
320///
321/// @tparam Scalar Scalar type.
322/// @param[in] matrix_callbacks Matrix callbacks.
323/// @param[in] iteration_callbacks The list of callbacks to call at the
324/// beginning of each iteration.
325/// @param[in] options Solver options.
326/// @param[in,out] x The current decision variables from the normal solve.
327/// @param[in,out] s The current inequality constraint slack variables from the
328/// normal solve.
329/// @param[in,out] y The current equality constraint duals from the normal
330/// solve.
331/// @param[in,out] z The current inequality constraint duals from the normal
332/// solve.
333/// @param[in] μ Barrier parameter.
334/// @return The exit status.
335template <typename Scalar>
336ExitStatus feasibility_restoration(
337 const InteriorPointMatrixCallbacks<Scalar>& matrix_callbacks,
338 std::span<std::function<bool(const IterationInfo<Scalar>& info)>>
339 iteration_callbacks,
340 const Options& options, Eigen::Vector<Scalar, Eigen::Dynamic>& x,
341 Eigen::Vector<Scalar, Eigen::Dynamic>& s,
342 Eigen::Vector<Scalar, Eigen::Dynamic>& y,
343 Eigen::Vector<Scalar, Eigen::Dynamic>& z, Scalar μ) {
344 // Feasibility restoration
345 //
346 // min ρ Σ (pₑ + nₑ + pᵢ + nᵢ) + ζ/2 (x - xᵣ)ᵀDᵣ(x - xᵣ)
347 // x
348 // pₑ,nₑ
349 // pᵢ,nᵢ
350 //
351 // s.t. cₑ(x) - pₑ + nₑ = 0
352 // cᵢ(x) - pᵢ + nᵢ ≥ 0
353 // pₑ ≥ 0
354 // nₑ ≥ 0
355 // pᵢ ≥ 0
356 // nᵢ ≥ 0
357 //
358 // where ρ = 1000, ζ = √μ where μ is the barrier parameter, xᵣ is original
359 // iterate before feasibility restoration, and Dᵣ is a scaling matrix defined
360 // by
361 //
362 // Dᵣ = diag(min(1, 1/xᵣ[i]²) for i in x.rows())
363
364 using DenseVector = Eigen::Vector<Scalar, Eigen::Dynamic>;
365 using DiagonalMatrix = Eigen::DiagonalMatrix<Scalar, Eigen::Dynamic>;
366 using SparseMatrix = Eigen::SparseMatrix<Scalar>;
367 using SparseVector = Eigen::SparseVector<Scalar>;
368
369 using std::sqrt;
370
371 const auto& matrices = matrix_callbacks;
372 const auto& num_vars = matrices.num_decision_variables;
373 const auto& num_eq = matrices.num_equality_constraints;
374 const auto& num_ineq = matrices.num_inequality_constraints;
375
376 constexpr Scalar ρ(1e3);
377 const Scalar ζ = sqrt(μ);
378
379 const DenseVector c_e = matrices.c_e(x);
380 const DenseVector c_i = matrices.c_i(x);
381
382 const auto& x_r = x;
383 const auto [p_e_0, n_e_0] = compute_p_n(c_e, ρ, μ);
384 const auto [p_i_0, n_i_0] = compute_p_n((c_i - s).eval(), ρ, μ);
385
386 // Dᵣ = diag(min(1, 1/xᵣ[i]²) for i in x.rows())
387 const DiagonalMatrix D_r =
388 x.cwiseSquare().cwiseInverse().cwiseMin(Scalar(1)).asDiagonal();
389
390 DenseVector fr_x{num_vars + 2 * num_eq + 2 * num_ineq};
391 fr_x << x, p_e_0, n_e_0, p_i_0, n_i_0;
392
393 DenseVector fr_s{s.rows() + 2 * num_eq + 2 * num_ineq};
394 fr_s.segment(0, s.rows()) = s;
395 fr_s.segment(s.rows(), 2 * num_eq + 2 * num_ineq).setOnes();
396
397 DenseVector fr_y = DenseVector::Zero(c_e.rows());
398
399 DenseVector fr_z{c_i.rows() + 2 * num_eq + 2 * num_ineq};
400 fr_z << z.cwiseMin(ρ), μ * p_e_0.cwiseInverse(), μ * n_e_0.cwiseInverse(),
401 μ * p_i_0.cwiseInverse(), μ * n_i_0.cwiseInverse();
402
403 Scalar fr_μ = std::max({μ, c_e.template lpNorm<Eigen::Infinity>(),
404 (c_i - s).template lpNorm<Eigen::Infinity>()});
405
406 InteriorPointMatrixCallbacks<Scalar> fr_matrix_callbacks{
407 static_cast<int>(fr_x.rows()),
408 static_cast<int>(fr_y.rows()),
409 static_cast<int>(fr_z.rows()),
410 [&](const DenseVector& x_p) -> Scalar {
411 auto x = x_p.segment(0, num_vars);
412
413 // Cost function
414 //
415 // ρ Σ (pₑ + nₑ + pᵢ + nᵢ) + ζ/2 (x - xᵣ)ᵀDᵣ(x - xᵣ)
416 auto diff = x - x_r;
417 return ρ * x_p.segment(num_vars, 2 * num_eq + 2 * num_ineq)
418 .array()
419 .sum() +
420 ζ / Scalar(2) * diff.transpose() * D_r * diff;
421 },
422 [&](const DenseVector& x_p) -> SparseVector {
423 auto x = x_p.segment(0, num_vars);
424
425 // Cost function gradient
426 //
427 // [ζDᵣ(x − xᵣ)]
428 // [ ρ ]
429 // [ ρ ]
430 // [ ρ ]
431 // [ ρ ]
432 DenseVector g{x_p.rows()};
433 g.segment(0, num_vars) = ζ * D_r * (x - x_r);
434 g.segment(num_vars, 2 * num_eq + 2 * num_ineq).setConstant(ρ);
435 return g.sparseView();
436 },
437 [&](const DenseVector& x_p, const DenseVector& y_p,
438 const DenseVector& z_p) -> SparseMatrix {
439 auto x = x_p.segment(0, num_vars);
440 const auto& y = y_p;
441 auto z = z_p.segment(0, num_ineq);
442
443 // Cost function Hessian
444 //
445 // [ζDᵣ 0 0 0 0]
446 // [ 0 0 0 0 0]
447 // [ 0 0 0 0 0]
448 // [ 0 0 0 0 0]
449 // [ 0 0 0 0 0]
451 triplets.reserve(x_p.rows());
452 append_as_triplets(triplets, 0, 0, {SparseMatrix{ζ * D_r}});
453 SparseMatrix d2f_dx2{x_p.rows(), x_p.rows()};
454 d2f_dx2.setFromSortedTriplets(triplets.begin(), triplets.end());
455
456 // Constraint part of original problem's Lagrangian Hessian
457 //
458 // −∇ₓₓ²yᵀcₑ(x) − ∇ₓₓ²zᵀcᵢ(x)
459 auto H_c = matrices.H_c(x, y, z);
460 H_c.resize(x_p.rows(), x_p.rows());
461
462 // Lagrangian Hessian
463 //
464 // [ζDᵣ 0 0 0 0]
465 // [ 0 0 0 0 0]
466 // [ 0 0 0 0 0] − ∇ₓₓ²yᵀcₑ(x) − ∇ₓₓ²zᵀcᵢ(x)
467 // [ 0 0 0 0 0]
468 // [ 0 0 0 0 0]
469 return d2f_dx2 + H_c;
470 },
471 [&](const DenseVector& x_p, [[maybe_unused]] const DenseVector& y_p,
472 [[maybe_unused]] const DenseVector& z_p) -> SparseMatrix {
473 return SparseMatrix{x_p.rows(), x_p.rows()};
474 },
475 [&](const DenseVector& x_p) -> DenseVector {
476 auto x = x_p.segment(0, num_vars);
477 auto p_e = x_p.segment(num_vars, num_eq);
478 auto n_e = x_p.segment(num_vars + num_eq, num_eq);
479
480 // Equality constraints
481 //
482 // cₑ(x) - pₑ + nₑ = 0
483 return matrices.c_e(x) - p_e + n_e;
484 },
485 [&](const DenseVector& x_p) -> SparseMatrix {
486 auto x = x_p.segment(0, num_vars);
487
488 // Equality constraint Jacobian
489 //
490 // [Aₑ −I I 0 0]
491
492 SparseMatrix A_e = matrices.A_e(x);
493
495 triplets.reserve(A_e.nonZeros() + 2 * num_eq);
496
497 append_as_triplets(triplets, 0, 0, {A_e});
499 triplets, 0, num_vars,
500 DenseVector::Constant(num_eq, Scalar(-1)).eval());
502 triplets, 0, num_vars + num_eq,
503 DenseVector::Constant(num_eq, Scalar(1)).eval());
504
505 SparseMatrix A_e_p{A_e.rows(), x_p.rows()};
506 A_e_p.setFromSortedTriplets(triplets.begin(), triplets.end());
507 return A_e_p;
508 },
509 [&](const DenseVector& x_p) -> DenseVector {
510 auto x = x_p.segment(0, num_vars);
511 auto p_i = x_p.segment(num_vars + 2 * num_eq, num_ineq);
512 auto n_i = x_p.segment(num_vars + 2 * num_eq + num_ineq, num_ineq);
513
514 // Inequality constraints
515 //
516 // cᵢ(x) - pᵢ + nᵢ ≥ 0
517 // pₑ ≥ 0
518 // nₑ ≥ 0
519 // pᵢ ≥ 0
520 // nᵢ ≥ 0
521 DenseVector c_i_p{c_i.rows() + 2 * num_eq + 2 * num_ineq};
522 c_i_p.segment(0, num_ineq) = matrices.c_i(x) - p_i + n_i;
523 c_i_p.segment(p_i.rows(), 2 * num_eq + 2 * num_ineq) =
524 x_p.segment(num_vars, 2 * num_eq + 2 * num_ineq);
525 return c_i_p;
526 },
527 [&](const DenseVector& x_p) -> SparseMatrix {
528 auto x = x_p.segment(0, num_vars);
529
530 // Inequality constraint Jacobian
531 //
532 // [Aᵢ 0 0 −I I]
533 // [0 I 0 0 0]
534 // [0 0 I 0 0]
535 // [0 0 0 I 0]
536 // [0 0 0 0 I]
537
538 SparseMatrix A_i = matrices.A_i(x);
539
541 triplets.reserve(A_i.nonZeros() + 2 * num_eq + 4 * num_ineq);
542
543 // Column 0
544 append_as_triplets(triplets, 0, 0, {A_i});
545
546 // Columns 1 and 2
548 triplets, num_ineq, num_vars,
549 DenseVector::Constant(2 * num_eq, Scalar(1)).eval());
550
551 SparseMatrix I_ineq{
552 DenseVector::Constant(num_ineq, Scalar(1)).asDiagonal()};
553
554 // Column 3
555 SparseMatrix Z_col3{2 * num_eq, num_ineq};
556 append_as_triplets(triplets, 0, num_vars + 2 * num_eq,
557 {(-I_ineq).eval(), Z_col3, I_ineq});
558
559 // Column 4
560 SparseMatrix Z_col4{2 * num_eq + num_ineq, num_ineq};
561 append_as_triplets(triplets, 0, num_vars + 2 * num_eq + num_ineq,
562 {I_ineq, Z_col4, I_ineq});
563
564 SparseMatrix A_i_p{2 * num_eq + 3 * num_ineq, x_p.rows()};
565 A_i_p.setFromSortedTriplets(triplets.begin(), triplets.end());
566 return A_i_p;
567 }};
568
569 auto status = interior_point<Scalar>(fr_matrix_callbacks, iteration_callbacks,
570 options, true,
571#ifdef SLEIPNIR_ENABLE_BOUND_PROJECTION
572 {},
573#endif
574 fr_x, fr_s, fr_y, fr_z, fr_μ);
575
576 x = fr_x.segment(0, x.rows());
577 s = fr_s.segment(0, s.rows());
578
580 auto g = matrices.g(x);
581 auto A_e = matrices.A_e(x);
582 auto A_i = matrices.A_i(x);
583
584 auto [y_estimate, z_estimate] =
585 lagrange_multiplier_estimate(g, A_e, A_i, s, μ);
586 y = y_estimate;
587 z = z_estimate;
588
589 return ExitStatus::SUCCESS;
590 } else if (status == ExitStatus::SUCCESS) {
592 } else {
594 }
595}
596
597// @endcond
598
599} // namespace slp
600
601#include "sleipnir/optimization/solver/interior_point.hpp"
void append_as_triplets(gch::small_vector< Eigen::Triplet< Scalar > > &triplets, int row_offset, int col_offset, std::initializer_list< Eigen::SparseMatrix< Scalar > > mats)
Appends sparse matrices to list of triplets at the given offset.
Definition append_as_triplets.hpp:22
void append_diagonal_as_triplets(gch::small_vector< Eigen::Triplet< Scalar > > &triplets, int row_offset, int col_offset, const Eigen::Vector< Scalar, Eigen::Dynamic > &diag)
Append diagonal matrix to list of triplets at the given offset.
Definition append_as_triplets.hpp:54
T * p
Definition format.h:758
wpi::util::SmallVector< T > small_vector
Definition small_vector.hpp:10
Definition to_underlying.hpp:7
ExitStatus
Solver exit status. Negative values indicate failure.
Definition exit_status.hpp:14
@ CALLBACK_REQUESTED_STOP
The solver returned its solution so far after the user requested a stop.
Definition exit_status.hpp:18
@ SUCCESS
Solved the problem to the desired tolerance.
Definition exit_status.hpp:16
@ FEASIBILITY_RESTORATION_FAILED
The solver failed to reach the desired tolerance, and feasibility restoration failed to converge.
Definition exit_status.hpp:33
@ LOCALLY_INFEASIBLE
The solver determined the problem to be locally infeasible and gave up.
Definition exit_status.hpp:22
ExitStatus interior_point(const InteriorPointMatrixCallbacks< Scalar > &matrix_callbacks, std::span< std::function< bool(const IterationInfo< Scalar > &info)> > iteration_callbacks, const Options &options, bool in_feasibility_restoration, Eigen::Vector< Scalar, Eigen::Dynamic > &x, Eigen::Vector< Scalar, Eigen::Dynamic > &s, Eigen::Vector< Scalar, Eigen::Dynamic > &y, Eigen::Vector< Scalar, Eigen::Dynamic > &z, Scalar &μ)
std::tuple< Eigen::Vector< Scalar, Eigen::Dynamic >, Eigen::Vector< Scalar, Eigen::Dynamic > > compute_p_n(const Eigen::Vector< Scalar, Eigen::Dynamic > &c, Scalar ρ, Scalar μ)
Computes initial values for p and n in feasibility restoration.
Definition feasibility_restoration.hpp:49
Eigen::Vector< Scalar, Eigen::Dynamic > lagrange_multiplier_estimate(const Eigen::SparseVector< Scalar > &g, const Eigen::SparseMatrix< Scalar > &A_e)
Estimates Lagrange multipliers for SQP.
Definition lagrange_multiplier_estimate.hpp:34
Variable< Scalar > sqrt(const Variable< Scalar > &x)
sqrt() for Variables.
Definition variable.hpp:671
Matrix callbacks for the interior-point method solver.
Definition interior_point_matrix_callbacks.hpp:16
Solver iteration information exposed to an iteration callback.
Definition iteration_info.hpp:14
Solver options.
Definition options.hpp:13