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