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