13#include <Eigen/SparseCore>
14#include <gch/small_vector.hpp>
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"
61template <
typename Scalar>
62ExitStatus interior_point(
63 const InteriorPointMatrixCallbacks<Scalar>& matrix_callbacks,
64 std::span<std::function<
bool(
const IterationInfo<Scalar>& info)>>
66 const Options& options,
67#ifdef SLEIPNIR_ENABLE_BOUND_PROJECTION
68 const Eigen::ArrayX<bool>& bound_constraint_mask,
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>;
89 const auto solve_start_time = std::chrono::steady_clock::now();
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");
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];
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];
133 InteriorPointMatrixCallbacks<Scalar> matrices{
134 [&](
const DenseVector& x) -> Scalar {
135 ScopedProfiler prof{f_prof};
136 return matrix_callbacks.f(x);
138 [&](
const DenseVector& x) -> SparseVector {
139 ScopedProfiler prof{g_prof};
140 return matrix_callbacks.g(x);
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);
147 [&](
const DenseVector& x) -> DenseVector {
148 ScopedProfiler prof{c_e_prof};
149 return matrix_callbacks.c_e(x);
151 [&](
const DenseVector& x) -> SparseMatrix {
152 ScopedProfiler prof{A_e_prof};
153 return matrix_callbacks.A_e(x);
155 [&](
const DenseVector& x) -> DenseVector {
156 ScopedProfiler prof{c_i_prof};
157 return matrix_callbacks.c_i(x);
159 [&](
const DenseVector& x) -> SparseMatrix {
160 ScopedProfiler prof{A_i_prof};
161 return matrix_callbacks.A_i(x);
164 const auto& matrices = matrix_callbacks;
170 Scalar f = matrices.f(x);
171 DenseVector c_e = matrices.c_e(x);
172 DenseVector c_i = matrices.c_i(x);
174 int num_decision_variables = x.rows();
175 int num_equality_constraints = c_e.rows();
176 int num_inequality_constraints = c_i.rows();
179 if (num_equality_constraints > num_decision_variables) {
180 if (options.diagnostics) {
181 print_too_few_dofs_error(c_e);
184 return ExitStatus::TOO_FEW_DOFS;
187 SparseVector g = matrices.g(x);
188 SparseMatrix A_e = matrices.A_e(x);
189 SparseMatrix A_i = matrices.A_i(x);
191 DenseVector s = DenseVector::Ones(num_inequality_constraints);
192#ifdef SLEIPNIR_ENABLE_BOUND_PROJECTION
194 s = bound_constraint_mask.select(c_i, s);
196 DenseVector y = DenseVector::Zero(num_equality_constraints);
197 DenseVector z = DenseVector::Ones(num_inequality_constraints);
199 SparseMatrix H = matrices.H(x, y, z);
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);
211 if (!isfinite(f) || !c_e.allFinite() || !c_i.allFinite()) {
212 return ExitStatus::NONFINITE_INITIAL_COST_OR_CONSTRAINTS;
218 const Scalar μ_min = Scalar(options.tolerance) / Scalar(10);
224 constexpr Scalar τ_min(0.99);
229 Filter<Scalar> filter;
233 auto update_barrier_parameter_and_reset_filter = [&] {
235 constexpr Scalar κ_μ(0.2);
239 constexpr Scalar θ_μ(1.5);
247 μ = std::max(μ_min, std::min(κ_μ * μ, pow(μ, θ_μ)));
254 τ = std::max(τ_min, Scalar(1) - μ);
261 gch::small_vector<Eigen::Triplet<Scalar>> triplets;
263 RegularizedLDLT<Scalar> solver{num_decision_variables,
264 num_equality_constraints};
267 constexpr Scalar α_reduction_factor(0.5);
268 constexpr Scalar α_min(1e-7);
270 int full_step_rejected_counter = 0;
273 Scalar E_0 = std::numeric_limits<Scalar>::infinity();
278 scope_exit exit{[&] {
279 if (options.diagnostics) {
281 if (iterations > 0) {
282 print_bottom_iteration_diagnostics();
284 print_solver_diagnostics(solve_profilers);
288 while (E_0 > Scalar(options.tolerance)) {
289 ScopedProfiler inner_iter_profiler{inner_iter_prof};
290 ScopedProfiler feasibility_check_profiler{feasibility_check_prof};
293 if (is_equality_locally_infeasible(A_e, c_e)) {
294 if (options.diagnostics) {
295 print_c_e_local_infeasibility_error(c_e);
298 return ExitStatus::LOCALLY_INFEASIBLE;
302 if (is_inequality_locally_infeasible(A_i, c_i)) {
303 if (options.diagnostics) {
304 print_c_i_local_infeasibility_error(c_i);
307 return ExitStatus::LOCALLY_INFEASIBLE;
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;
316 feasibility_check_profiler.stop();
317 ScopedProfiler iter_callbacks_profiler{iter_callbacks_prof};
320 for (
const auto& callback : iteration_callbacks) {
321 if (callback({iterations, x, g, H, A_e, A_i})) {
322 return ExitStatus::CALLBACK_REQUESTED_STOP;
326 iter_callbacks_profiler.stop();
327 ScopedProfiler kkt_matrix_build_profiler{kkt_matrix_build_prof};
332 const SparseMatrix Σ{s.cwiseInverse().asDiagonal() * z.asDiagonal()};
338 const SparseMatrix top_left =
339 H + (A_i.transpose() * Σ * A_i).
template triangularView<Eigen::Lower>();
341 triplets.reserve(top_left.nonZeros() + A_e.nonZeros());
342 for (
int col = 0; col < H.cols(); ++col) {
344 for (
typename SparseMatrix::InnerIterator it{top_left, col}; it; ++it) {
345 triplets.emplace_back(it.row(), it.col(), it.value());
348 for (
typename SparseMatrix::InnerIterator it{A_e, col}; it; ++it) {
349 triplets.emplace_back(H.rows() + it.row(), it.col(), it.value());
352 SparseMatrix lhs(num_decision_variables + num_equality_constraints,
353 num_decision_variables + num_equality_constraints);
354 lhs.setFromSortedTriplets(triplets.begin(), triplets.end());
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;
364 kkt_matrix_build_profiler.stop();
365 ScopedProfiler kkt_matrix_decomp_profiler{kkt_matrix_decomp_prof};
376 if (solver.compute(lhs).info() != Eigen::Success) [[unlikely]] {
377 return ExitStatus::FACTORIZATION_FAILED;
380 kkt_matrix_decomp_profiler.stop();
381 ScopedProfiler kkt_system_solve_profiler{kkt_system_solve_prof};
383 auto compute_step = [&](Step& step) {
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());
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;
397 kkt_system_solve_profiler.stop();
398 ScopedProfiler line_search_profiler{line_search_prof};
401 α_max = fraction_to_the_boundary_rule<Scalar>(s, step.p_s, τ);
406 return ExitStatus::LINE_SEARCH_FAILED;
410 α_z = fraction_to_the_boundary_rule<Scalar>(z, step.p_z, τ);
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;
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);
424 if (!isfinite(trial_f) || !trial_c_e.allFinite() ||
425 !trial_c_i.allFinite()) {
427 α *= α_reduction_factor;
430 return ExitStatus::LINE_SEARCH_FAILED;
436 if (options.feasible_ipm && c_i.cwiseGreater(Scalar(0)).all()) {
443 trial_s = s + α * step.p_s;
447 if (filter.try_add(FilterEntry{trial_f, trial_s, trial_c_e, trial_c_i, μ},
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>();
464 next_constraint_violation >= prev_constraint_violation) {
466 auto soc_step = step;
469 Scalar α_z_soc = α_z;
470 DenseVector c_e_soc = c_e;
472 bool step_acceptable =
false;
473 for (
int soc_iteration = 0; soc_iteration < 5 && !step_acceptable;
475 ScopedProfiler soc_profiler{soc_prof};
477 scope_exit soc_exit{[&] {
480 if (options.diagnostics) {
481 print_iteration_diagnostics(
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)),
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);
502 c_e_soc = α_soc * c_e_soc + trial_c_e;
503 rhs.bottomRows(y.rows()) = -c_e_soc;
506 compute_step(soc_step);
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;
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;
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);
523 constexpr Scalar κ_soc(0.99);
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) {
536 FilterEntry{trial_f, trial_s, trial_c_e, trial_c_i, μ}, α)) {
540 step_acceptable =
true;
544 if (step_acceptable) {
554 ++full_step_rejected_counter;
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);
570 α *= α_reduction_factor;
575 Scalar current_kkt_error =
576 kkt_error<Scalar>(g, A_e, c_e, A_i, c_i, s, y, z, μ);
578 trial_x = x + α_max * step.p_x;
579 trial_s = s + α_max * step.p_s;
581 trial_y = y + α_z * step.p_y;
582 trial_z = z + α_z * step.p_z;
584 trial_c_e = matrices.c_e(trial_x);
585 trial_c_i = matrices.c_i(trial_x);
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, μ);
592 if (next_kkt_error <= Scalar(0.999) * current_kkt_error) {
599 return ExitStatus::LINE_SEARCH_FAILED;
603 line_search_profiler.stop();
607 full_step_rejected_counter = 0;
632 for (
int row = 0; row < z.rows(); ++row) {
633 constexpr Scalar κ_Σ(1e10);
635 std::clamp(z[row], Scalar(1) / κ_Σ * μ / s[row], κ_Σ * μ / s[row]);
640 A_e = matrices.A_e(x);
641 A_i = matrices.A_i(x);
643 H = matrices.H(x, y, z);
645 ScopedProfiler next_iter_prep_profiler{next_iter_prep_prof};
647 c_e = matrices.c_e(x);
648 c_i = matrices.c_i(x);
651 E_0 = error_estimate<Scalar>(g, A_e, c_e, A_i, c_i, s, y, z, Scalar(0));
654 if (E_0 > Scalar(options.tolerance)) {
656 constexpr Scalar κ_ε(10);
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, μ);
667 next_iter_prep_profiler.stop();
668 inner_iter_profiler.stop();
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,
682 if (iterations >= options.max_iterations) {
683 return ExitStatus::MAX_ITERATIONS_EXCEEDED;
687 if (std::chrono::steady_clock::now() - solve_start_time > options.timeout) {
688 return ExitStatus::TIMEOUT;
692 return ExitStatus::SUCCESS;
695extern template SLEIPNIR_DLLEXPORT ExitStatus
696interior_point(
const InteriorPointMatrixCallbacks<double>& matrix_callbacks,
697 std::span<std::function<
bool(
const IterationInfo<double>& info)>>
699 const Options& options,
700#ifdef SLEIPNIR_ENABLE_BOUND_PROJECTION
701 const Eigen::ArrayX<bool>& bound_constraint_mask,
703 Eigen::Vector<double, Eigen::Dynamic>& x);