Sleipnir C++ API
Loading...
Searching...
No Matches
interior_point.hpp
1// Copyright (c) Sleipnir contributors
2
3#pragma once
4
5#include <algorithm>
6#include <chrono>
7#include <cmath>
8#include <functional>
9#include <limits>
10#include <span>
11
12#include <Eigen/Core>
13#include <Eigen/SparseCore>
14#include <gch/small_vector.hpp>
15
16#include "sleipnir/optimization/solver/exit_status.hpp"
17#include "sleipnir/optimization/solver/interior_point_matrix_callbacks.hpp"
18#include "sleipnir/optimization/solver/iteration_info.hpp"
19#include "sleipnir/optimization/solver/options.hpp"
20#include "sleipnir/optimization/solver/util/error_estimate.hpp"
21#include "sleipnir/optimization/solver/util/filter.hpp"
22#include "sleipnir/optimization/solver/util/fraction_to_the_boundary_rule.hpp"
23#include "sleipnir/optimization/solver/util/is_locally_infeasible.hpp"
24#include "sleipnir/optimization/solver/util/kkt_error.hpp"
25#include "sleipnir/optimization/solver/util/regularized_ldlt.hpp"
26#include "sleipnir/util/assert.hpp"
27#include "sleipnir/util/print_diagnostics.hpp"
28#include "sleipnir/util/profiler.hpp"
29#include "sleipnir/util/scope_exit.hpp"
30#include "sleipnir/util/symbol_exports.hpp"
31
32// See docs/algorithms.md#Works_cited for citation definitions.
33//
34// See docs/algorithms.md#Interior-point_method for a derivation of the
35// interior-point method formulation being used.
36
37namespace slp {
38
61template <typename Scalar>
62ExitStatus interior_point(
63 const InteriorPointMatrixCallbacks<Scalar>& matrix_callbacks,
64 std::span<std::function<bool(const IterationInfo<Scalar>& info)>>
65 iteration_callbacks,
66 const Options& options,
67#ifdef SLEIPNIR_ENABLE_BOUND_PROJECTION
68 const Eigen::ArrayX<bool>& bound_constraint_mask,
69#endif
70 Eigen::Vector<Scalar, Eigen::Dynamic>& x) {
71 using DenseVector = Eigen::Vector<Scalar, Eigen::Dynamic>;
72 using SparseMatrix = Eigen::SparseMatrix<Scalar>;
73 using SparseVector = Eigen::SparseVector<Scalar>;
74
76 struct Step {
78 DenseVector p_x;
80 DenseVector p_y;
82 DenseVector p_s;
84 DenseVector p_z;
85 };
86
87 using std::isfinite;
88
89 const auto solve_start_time = std::chrono::steady_clock::now();
90
91 gch::small_vector<SolveProfiler> solve_profilers;
92 solve_profilers.emplace_back("solver");
93 solve_profilers.emplace_back(" ↳ setup");
94 solve_profilers.emplace_back(" ↳ iteration");
95 solve_profilers.emplace_back(" ↳ feasibility ✓");
96 solve_profilers.emplace_back(" ↳ iter callbacks");
97 solve_profilers.emplace_back(" ↳ KKT matrix build");
98 solve_profilers.emplace_back(" ↳ KKT matrix decomp");
99 solve_profilers.emplace_back(" ↳ KKT system solve");
100 solve_profilers.emplace_back(" ↳ line search");
101 solve_profilers.emplace_back(" ↳ SOC");
102 solve_profilers.emplace_back(" ↳ next iter prep");
103 solve_profilers.emplace_back(" ↳ f(x)");
104 solve_profilers.emplace_back(" ↳ ∇f(x)");
105 solve_profilers.emplace_back(" ↳ ∇²ₓₓL");
106 solve_profilers.emplace_back(" ↳ cₑ(x)");
107 solve_profilers.emplace_back(" ↳ ∂cₑ/∂x");
108 solve_profilers.emplace_back(" ↳ cᵢ(x)");
109 solve_profilers.emplace_back(" ↳ ∂cᵢ/∂x");
110
111 auto& solver_prof = solve_profilers[0];
112 auto& setup_prof = solve_profilers[1];
113 auto& inner_iter_prof = solve_profilers[2];
114 auto& feasibility_check_prof = solve_profilers[3];
115 auto& iter_callbacks_prof = solve_profilers[4];
116 auto& kkt_matrix_build_prof = solve_profilers[5];
117 auto& kkt_matrix_decomp_prof = solve_profilers[6];
118 auto& kkt_system_solve_prof = solve_profilers[7];
119 auto& line_search_prof = solve_profilers[8];
120 auto& soc_prof = solve_profilers[9];
121 auto& next_iter_prep_prof = solve_profilers[10];
122
123 // Set up profiled matrix callbacks
124#ifndef SLEIPNIR_DISABLE_DIAGNOSTICS
125 auto& f_prof = solve_profilers[11];
126 auto& g_prof = solve_profilers[12];
127 auto& H_prof = solve_profilers[13];
128 auto& c_e_prof = solve_profilers[14];
129 auto& A_e_prof = solve_profilers[15];
130 auto& c_i_prof = solve_profilers[16];
131 auto& A_i_prof = solve_profilers[17];
132
133 InteriorPointMatrixCallbacks<Scalar> matrices{
134 [&](const DenseVector& x) -> Scalar {
135 ScopedProfiler prof{f_prof};
136 return matrix_callbacks.f(x);
137 },
138 [&](const DenseVector& x) -> SparseVector {
139 ScopedProfiler prof{g_prof};
140 return matrix_callbacks.g(x);
141 },
142 [&](const DenseVector& x, const DenseVector& y,
143 const DenseVector& z) -> SparseMatrix {
144 ScopedProfiler prof{H_prof};
145 return matrix_callbacks.H(x, y, z);
146 },
147 [&](const DenseVector& x) -> DenseVector {
148 ScopedProfiler prof{c_e_prof};
149 return matrix_callbacks.c_e(x);
150 },
151 [&](const DenseVector& x) -> SparseMatrix {
152 ScopedProfiler prof{A_e_prof};
153 return matrix_callbacks.A_e(x);
154 },
155 [&](const DenseVector& x) -> DenseVector {
156 ScopedProfiler prof{c_i_prof};
157 return matrix_callbacks.c_i(x);
158 },
159 [&](const DenseVector& x) -> SparseMatrix {
160 ScopedProfiler prof{A_i_prof};
161 return matrix_callbacks.A_i(x);
162 }};
163#else
164 const auto& matrices = matrix_callbacks;
165#endif
166
167 solver_prof.start();
168 setup_prof.start();
169
170 Scalar f = matrices.f(x);
171 DenseVector c_e = matrices.c_e(x);
172 DenseVector c_i = matrices.c_i(x);
173
174 int num_decision_variables = x.rows();
175 int num_equality_constraints = c_e.rows();
176 int num_inequality_constraints = c_i.rows();
177
178 // Check for overconstrained problem
179 if (num_equality_constraints > num_decision_variables) {
180 if (options.diagnostics) {
181 print_too_few_dofs_error(c_e);
182 }
183
184 return ExitStatus::TOO_FEW_DOFS;
185 }
186
187 SparseVector g = matrices.g(x);
188 SparseMatrix A_e = matrices.A_e(x);
189 SparseMatrix A_i = matrices.A_i(x);
190
191 DenseVector s = DenseVector::Ones(num_inequality_constraints);
192#ifdef SLEIPNIR_ENABLE_BOUND_PROJECTION
193 // We set sʲ = cᵢʲ(x) for each bound inequality constraint index j
194 s = bound_constraint_mask.select(c_i, s);
195#endif
196 DenseVector y = DenseVector::Zero(num_equality_constraints);
197 DenseVector z = DenseVector::Ones(num_inequality_constraints);
198
199 SparseMatrix H = matrices.H(x, y, z);
200
201 // Ensure matrix callback dimensions are consistent
202 slp_assert(g.rows() == num_decision_variables);
203 slp_assert(A_e.rows() == num_equality_constraints);
204 slp_assert(A_e.cols() == num_decision_variables);
205 slp_assert(A_i.rows() == num_inequality_constraints);
206 slp_assert(A_i.cols() == num_decision_variables);
207 slp_assert(H.rows() == num_decision_variables);
208 slp_assert(H.cols() == num_decision_variables);
209
210 // Check whether initial guess has finite f(xₖ), cₑ(xₖ), and cᵢ(xₖ)
211 if (!isfinite(f) || !c_e.allFinite() || !c_i.allFinite()) {
212 return ExitStatus::NONFINITE_INITIAL_COST_OR_CONSTRAINTS;
213 }
214
215 int iterations = 0;
216
217 // Barrier parameter minimum
218 const Scalar μ_min = Scalar(options.tolerance) / Scalar(10);
219
220 // Barrier parameter μ
221 Scalar μ(0.1);
222
223 // Fraction-to-the-boundary rule scale factor minimum
224 constexpr Scalar τ_min(0.99);
225
226 // Fraction-to-the-boundary rule scale factor τ
227 Scalar τ = τ_min;
228
229 Filter<Scalar> filter;
230
231 // This should be run when the error estimate is below a desired threshold for
232 // the current barrier parameter
233 auto update_barrier_parameter_and_reset_filter = [&] {
234 // Barrier parameter linear decrease power in "κ_μ μ". Range of (0, 1).
235 constexpr Scalar κ_μ(0.2);
236
237 // Barrier parameter superlinear decrease power in "μ^(θ_μ)". Range of (1,
238 // 2).
239 constexpr Scalar θ_μ(1.5);
240
241 // Update the barrier parameter.
242 //
243 // μⱼ₊₁ = max(εₜₒₗ/10, min(κ_μ μⱼ, μⱼ^θ_μ))
244 //
245 // See equation (7) of [2].
246 using std::pow;
247 μ = std::max(μ_min, std::min(κ_μ * μ, pow(μ, θ_μ)));
248
249 // Update the fraction-to-the-boundary rule scaling factor.
250 //
251 // τⱼ = max(τₘᵢₙ, 1 − μⱼ)
252 //
253 // See equation (8) of [2].
254 τ = std::max(τ_min, Scalar(1) - μ);
255
256 // Reset the filter when the barrier parameter is updated
257 filter.reset();
258 };
259
260 // Kept outside the loop so its storage can be reused
261 gch::small_vector<Eigen::Triplet<Scalar>> triplets;
262
263 RegularizedLDLT<Scalar> solver{num_decision_variables,
264 num_equality_constraints};
265
266 // Variables for determining when a step is acceptable
267 constexpr Scalar α_reduction_factor(0.5);
268 constexpr Scalar α_min(1e-7);
269
270 int full_step_rejected_counter = 0;
271
272 // Error estimate
273 Scalar E_0 = std::numeric_limits<Scalar>::infinity();
274
275 setup_prof.stop();
276
277 // Prints final solver diagnostics when the solver exits
278 scope_exit exit{[&] {
279 if (options.diagnostics) {
280 solver_prof.stop();
281 if (iterations > 0) {
282 print_bottom_iteration_diagnostics();
283 }
284 print_solver_diagnostics(solve_profilers);
285 }
286 }};
287
288 while (E_0 > Scalar(options.tolerance)) {
289 ScopedProfiler inner_iter_profiler{inner_iter_prof};
290 ScopedProfiler feasibility_check_profiler{feasibility_check_prof};
291
292 // Check for local equality constraint infeasibility
293 if (is_equality_locally_infeasible(A_e, c_e)) {
294 if (options.diagnostics) {
295 print_c_e_local_infeasibility_error(c_e);
296 }
297
298 return ExitStatus::LOCALLY_INFEASIBLE;
299 }
300
301 // Check for local inequality constraint infeasibility
302 if (is_inequality_locally_infeasible(A_i, c_i)) {
303 if (options.diagnostics) {
304 print_c_i_local_infeasibility_error(c_i);
305 }
306
307 return ExitStatus::LOCALLY_INFEASIBLE;
308 }
309
310 // Check for diverging iterates
311 if (x.template lpNorm<Eigen::Infinity>() > Scalar(1e10) || !x.allFinite() ||
312 s.template lpNorm<Eigen::Infinity>() > Scalar(1e10) || !s.allFinite()) {
313 return ExitStatus::DIVERGING_ITERATES;
314 }
315
316 feasibility_check_profiler.stop();
317 ScopedProfiler iter_callbacks_profiler{iter_callbacks_prof};
318
319 // Call iteration callbacks
320 for (const auto& callback : iteration_callbacks) {
321 if (callback({iterations, x, g, H, A_e, A_i})) {
322 return ExitStatus::CALLBACK_REQUESTED_STOP;
323 }
324 }
325
326 iter_callbacks_profiler.stop();
327 ScopedProfiler kkt_matrix_build_profiler{kkt_matrix_build_prof};
328
329 // S = diag(s)
330 // Z = diag(z)
331 // Σ = S⁻¹Z
332 const SparseMatrix Σ{s.cwiseInverse().asDiagonal() * z.asDiagonal()};
333
334 // lhs = [H + AᵢᵀΣAᵢ Aₑᵀ]
335 // [ Aₑ 0 ]
336 //
337 // Don't assign upper triangle because solver only uses lower triangle.
338 const SparseMatrix top_left =
339 H + (A_i.transpose() * Σ * A_i).template triangularView<Eigen::Lower>();
340 triplets.clear();
341 triplets.reserve(top_left.nonZeros() + A_e.nonZeros());
342 for (int col = 0; col < H.cols(); ++col) {
343 // Append column of H + AᵢᵀΣAᵢ lower triangle in top-left quadrant
344 for (typename SparseMatrix::InnerIterator it{top_left, col}; it; ++it) {
345 triplets.emplace_back(it.row(), it.col(), it.value());
346 }
347 // Append column of Aₑ in bottom-left quadrant
348 for (typename SparseMatrix::InnerIterator it{A_e, col}; it; ++it) {
349 triplets.emplace_back(H.rows() + it.row(), it.col(), it.value());
350 }
351 }
352 SparseMatrix lhs(num_decision_variables + num_equality_constraints,
353 num_decision_variables + num_equality_constraints);
354 lhs.setFromSortedTriplets(triplets.begin(), triplets.end());
355
356 // rhs = −[∇f − Aₑᵀy − Aᵢᵀ(−Σcᵢ + μS⁻¹e + z)]
357 // [ cₑ ]
358 DenseVector rhs{x.rows() + y.rows()};
359 rhs.segment(0, x.rows()) =
360 -g + A_e.transpose() * y +
361 A_i.transpose() * (-Σ * c_i + μ * s.cwiseInverse() + z);
362 rhs.segment(x.rows(), y.rows()) = -c_e;
363
364 kkt_matrix_build_profiler.stop();
365 ScopedProfiler kkt_matrix_decomp_profiler{kkt_matrix_decomp_prof};
366
367 Step step;
368 Scalar α_max(1);
369 Scalar α(1);
370 Scalar α_z(1);
371
372 // Solve the Newton-KKT system
373 //
374 // [H + AᵢᵀΣAᵢ Aₑᵀ][ pˣ] = −[∇f − Aₑᵀy − Aᵢᵀ(−Σcᵢ + μS⁻¹e + z)]
375 // [ Aₑ 0 ][−pʸ] [ cₑ ]
376 if (solver.compute(lhs).info() != Eigen::Success) [[unlikely]] {
377 return ExitStatus::FACTORIZATION_FAILED;
378 }
379
380 kkt_matrix_decomp_profiler.stop();
381 ScopedProfiler kkt_system_solve_profiler{kkt_system_solve_prof};
382
383 auto compute_step = [&](Step& step) {
384 // p = [ pˣ]
385 // [−pʸ]
386 DenseVector p = solver.solve(rhs);
387 step.p_x = p.segment(0, x.rows());
388 step.p_y = -p.segment(x.rows(), y.rows());
389
390 // pˢ = cᵢ − s + Aᵢpˣ
391 // pᶻ = −Σcᵢ + μS⁻¹e − ΣAᵢpˣ
392 step.p_s = c_i - s + A_i * step.p_x;
393 step.p_z = -Σ * c_i + μ * s.cwiseInverse() - Σ * A_i * step.p_x;
394 };
395 compute_step(step);
396
397 kkt_system_solve_profiler.stop();
398 ScopedProfiler line_search_profiler{line_search_prof};
399
400 // αᵐᵃˣ = max(α ∈ (0, 1] : sₖ + αpₖˢ ≥ (1−τⱼ)sₖ)
401 α_max = fraction_to_the_boundary_rule<Scalar>(s, step.p_s, τ);
402 α = α_max;
403
404 // If maximum step size is below minimum, report line search failure
405 if (α < α_min) {
406 return ExitStatus::LINE_SEARCH_FAILED;
407 }
408
409 // αₖᶻ = max(α ∈ (0, 1] : zₖ + αpₖᶻ ≥ (1−τⱼ)zₖ)
410 α_z = fraction_to_the_boundary_rule<Scalar>(z, step.p_z, τ);
411
412 // Loop until a step is accepted
413 while (1) {
414 DenseVector trial_x = x + α * step.p_x;
415 DenseVector trial_y = y + α_z * step.p_y;
416 DenseVector trial_z = z + α_z * step.p_z;
417
418 Scalar trial_f = matrices.f(trial_x);
419 DenseVector trial_c_e = matrices.c_e(trial_x);
420 DenseVector trial_c_i = matrices.c_i(trial_x);
421
422 // If f(xₖ + αpₖˣ), cₑ(xₖ + αpₖˣ), or cᵢ(xₖ + αpₖˣ) aren't finite, reduce
423 // step size immediately
424 if (!isfinite(trial_f) || !trial_c_e.allFinite() ||
425 !trial_c_i.allFinite()) {
426 // Reduce step size
427 α *= α_reduction_factor;
428
429 if (α < α_min) {
430 return ExitStatus::LINE_SEARCH_FAILED;
431 }
432 continue;
433 }
434
435 DenseVector trial_s;
436 if (options.feasible_ipm && c_i.cwiseGreater(Scalar(0)).all()) {
437 // If the inequality constraints are all feasible, prevent them from
438 // becoming infeasible again.
439 //
440 // See equation (19.30) in [1].
441 trial_s = trial_c_i;
442 } else {
443 trial_s = s + α * step.p_s;
444 }
445
446 // Check whether filter accepts trial iterate
447 if (filter.try_add(FilterEntry{trial_f, trial_s, trial_c_e, trial_c_i, μ},
448 α)) {
449 // Accept step
450 break;
451 }
452
453 Scalar prev_constraint_violation =
454 c_e.template lpNorm<1>() + (c_i - s).template lpNorm<1>();
455 Scalar next_constraint_violation =
456 trial_c_e.template lpNorm<1>() +
457 (trial_c_i - trial_s).template lpNorm<1>();
458
459 // Second-order corrections
460 //
461 // If first trial point was rejected and constraint violation stayed the
462 // same or went up, apply second-order corrections
463 if (α == α_max &&
464 next_constraint_violation >= prev_constraint_violation) {
465 // Apply second-order corrections. See section 2.4 of [2].
466 auto soc_step = step;
467
468 Scalar α_soc = α;
469 Scalar α_z_soc = α_z;
470 DenseVector c_e_soc = c_e;
471
472 bool step_acceptable = false;
473 for (int soc_iteration = 0; soc_iteration < 5 && !step_acceptable;
474 ++soc_iteration) {
475 ScopedProfiler soc_profiler{soc_prof};
476
477 scope_exit soc_exit{[&] {
478 soc_profiler.stop();
479
480 if (options.diagnostics) {
481 print_iteration_diagnostics(
482 iterations,
483 step_acceptable ? IterationType::ACCEPTED_SOC
484 : IterationType::REJECTED_SOC,
485 soc_profiler.current_duration(),
486 error_estimate<Scalar>(g, A_e, trial_c_e, A_i, trial_c_i,
487 trial_s, trial_y, trial_z, Scalar(0)),
488 trial_f,
489 trial_c_e.template lpNorm<1>() +
490 (trial_c_i - trial_s).template lpNorm<1>(),
491 trial_s.dot(trial_z), μ, solver.hessian_regularization(),
492 α_soc, Scalar(1), α_reduction_factor, α_z_soc);
493 }
494 }};
495
496 // Rebuild Newton-KKT rhs with updated constraint values.
497 //
498 // rhs = −[∇f − Aₑᵀy − Aᵢᵀ(−Σcᵢ + μS⁻¹e + z)]
499 // [ cₑˢᵒᶜ ]
500 //
501 // where cₑˢᵒᶜ = αc(xₖ) + c(xₖ + αpₖˣ)
502 c_e_soc = α_soc * c_e_soc + trial_c_e;
503 rhs.bottomRows(y.rows()) = -c_e_soc;
504
505 // Solve the Newton-KKT system
506 compute_step(soc_step);
507
508 // αˢᵒᶜ = max(α ∈ (0, 1] : sₖ + αpₖˢ ≥ (1−τⱼ)sₖ)
509 α_soc = fraction_to_the_boundary_rule<Scalar>(s, soc_step.p_s, τ);
510 trial_x = x + α_soc * soc_step.p_x;
511 trial_s = s + α_soc * soc_step.p_s;
512
513 // αₖᶻ = max(α ∈ (0, 1] : zₖ + αpₖᶻ ≥ (1−τⱼ)zₖ)
514 α_z_soc = fraction_to_the_boundary_rule<Scalar>(z, soc_step.p_z, τ);
515 trial_y = y + α_z_soc * soc_step.p_y;
516 trial_z = z + α_z_soc * soc_step.p_z;
517
518 trial_f = matrices.f(trial_x);
519 trial_c_e = matrices.c_e(trial_x);
520 trial_c_i = matrices.c_i(trial_x);
521
522 // Constraint violation scale factor for second-order corrections
523 constexpr Scalar κ_soc(0.99);
524
525 // If constraint violation hasn't been sufficiently reduced, stop
526 // making second-order corrections
527 next_constraint_violation =
528 trial_c_e.template lpNorm<1>() +
529 (trial_c_i - trial_s).template lpNorm<1>();
530 if (next_constraint_violation > κ_soc * prev_constraint_violation) {
531 break;
532 }
533
534 // Check whether filter accepts trial iterate
535 if (filter.try_add(
536 FilterEntry{trial_f, trial_s, trial_c_e, trial_c_i, μ}, α)) {
537 step = soc_step;
538 α = α_soc;
539 α_z = α_z_soc;
540 step_acceptable = true;
541 }
542 }
543
544 if (step_acceptable) {
545 // Accept step
546 break;
547 }
548 }
549
550 // If we got here and α is the full step, the full step was rejected.
551 // Increment the full-step rejected counter to keep track of how many full
552 // steps have been rejected in a row.
553 if (α == α_max) {
554 ++full_step_rejected_counter;
555 }
556
557 // If the full step was rejected enough times in a row, reset the filter
558 // because it may be impeding progress.
559 //
560 // See section 3.2 case I of [2].
561 if (full_step_rejected_counter >= 4 &&
562 filter.max_constraint_violation >
563 filter.back().constraint_violation / Scalar(10)) {
564 filter.max_constraint_violation *= Scalar(0.1);
565 filter.reset();
566 continue;
567 }
568
569 // Reduce step size
570 α *= α_reduction_factor;
571
572 // If step size hit a minimum, check if the KKT error was reduced. If it
573 // wasn't, report line search failure.
574 if (α < α_min) {
575 Scalar current_kkt_error =
576 kkt_error<Scalar>(g, A_e, c_e, A_i, c_i, s, y, z, μ);
577
578 trial_x = x + α_max * step.p_x;
579 trial_s = s + α_max * step.p_s;
580
581 trial_y = y + α_z * step.p_y;
582 trial_z = z + α_z * step.p_z;
583
584 trial_c_e = matrices.c_e(trial_x);
585 trial_c_i = matrices.c_i(trial_x);
586
587 Scalar next_kkt_error = kkt_error<Scalar>(
588 matrices.g(trial_x), matrices.A_e(trial_x), matrices.c_e(trial_x),
589 matrices.A_i(trial_x), trial_c_i, trial_s, trial_y, trial_z, μ);
590
591 // If the step using αᵐᵃˣ reduced the KKT error, accept it anyway
592 if (next_kkt_error <= Scalar(0.999) * current_kkt_error) {
593 α = α_max;
594
595 // Accept step
596 break;
597 }
598
599 return ExitStatus::LINE_SEARCH_FAILED;
600 }
601 }
602
603 line_search_profiler.stop();
604
605 // If full step was accepted, reset full-step rejected counter
606 if (α == α_max) {
607 full_step_rejected_counter = 0;
608 }
609
610 // xₖ₊₁ = xₖ + αₖpₖˣ
611 // sₖ₊₁ = sₖ + αₖpₖˢ
612 // yₖ₊₁ = yₖ + αₖᶻpₖʸ
613 // zₖ₊₁ = zₖ + αₖᶻpₖᶻ
614 x += α * step.p_x;
615 s += α * step.p_s;
616 y += α_z * step.p_y;
617 z += α_z * step.p_z;
618
619 // A requirement for the convergence proof is that the primal-dual barrier
620 // term Hessian Σₖ₊₁ does not deviate arbitrarily much from the primal
621 // barrier term Hessian μSₖ₊₁⁻².
622 //
623 // Σₖ₊₁ = μSₖ₊₁⁻²
624 // Sₖ₊₁⁻¹Zₖ₊₁ = μSₖ₊₁⁻²
625 // Zₖ₊₁ = μSₖ₊₁⁻¹
626 //
627 // We ensure this by resetting
628 //
629 // zₖ₊₁ = clamp(zₖ₊₁, 1/κ_Σ μ/sₖ₊₁, κ_Σ μ/sₖ₊₁)
630 //
631 // for some fixed κ_Σ ≥ 1 after each step. See equation (16) of [2].
632 for (int row = 0; row < z.rows(); ++row) {
633 constexpr Scalar κ_Σ(1e10);
634 z[row] =
635 std::clamp(z[row], Scalar(1) / κ_Σ * μ / s[row], κ_Σ * μ / s[row]);
636 }
637
638 // Update autodiff for Jacobians and Hessian
639 f = matrices.f(x);
640 A_e = matrices.A_e(x);
641 A_i = matrices.A_i(x);
642 g = matrices.g(x);
643 H = matrices.H(x, y, z);
644
645 ScopedProfiler next_iter_prep_profiler{next_iter_prep_prof};
646
647 c_e = matrices.c_e(x);
648 c_i = matrices.c_i(x);
649
650 // Update the error estimate
651 E_0 = error_estimate<Scalar>(g, A_e, c_e, A_i, c_i, s, y, z, Scalar(0));
652
653 // Update the barrier parameter if necessary
654 if (E_0 > Scalar(options.tolerance)) {
655 // Barrier parameter scale factor for tolerance checks
656 constexpr Scalar κ_ε(10);
657
658 // While the error estimate is below the desired threshold for this
659 // barrier parameter value, decrease the barrier parameter further
660 Scalar E_μ = error_estimate<Scalar>(g, A_e, c_e, A_i, c_i, s, y, z, μ);
661 while (μ > μ_min && E_μ <= κ_ε * μ) {
662 update_barrier_parameter_and_reset_filter();
663 E_μ = error_estimate<Scalar>(g, A_e, c_e, A_i, c_i, s, y, z, μ);
664 }
665 }
666
667 next_iter_prep_profiler.stop();
668 inner_iter_profiler.stop();
669
670 if (options.diagnostics) {
671 print_iteration_diagnostics(
672 iterations, IterationType::NORMAL,
673 inner_iter_profiler.current_duration(), E_0, f,
674 c_e.template lpNorm<1>() + (c_i - s).template lpNorm<1>(), s.dot(z),
675 μ, solver.hessian_regularization(), α, α_max, α_reduction_factor,
676 α_z);
677 }
678
679 ++iterations;
680
681 // Check for max iterations
682 if (iterations >= options.max_iterations) {
683 return ExitStatus::MAX_ITERATIONS_EXCEEDED;
684 }
685
686 // Check for max wall clock time
687 if (std::chrono::steady_clock::now() - solve_start_time > options.timeout) {
688 return ExitStatus::TIMEOUT;
689 }
690 }
691
692 return ExitStatus::SUCCESS;
693}
694
695extern template SLEIPNIR_DLLEXPORT ExitStatus
696interior_point(const InteriorPointMatrixCallbacks<double>& matrix_callbacks,
697 std::span<std::function<bool(const IterationInfo<double>& info)>>
698 iteration_callbacks,
699 const Options& options,
700#ifdef SLEIPNIR_ENABLE_BOUND_PROJECTION
701 const Eigen::ArrayX<bool>& bound_constraint_mask,
702#endif
703 Eigen::Vector<double, Eigen::Dynamic>& x);
704
705} // namespace slp