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/all_finite.hpp"
21#include "sleipnir/optimization/solver/util/append_as_triplets.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"
62template <
typename Scalar>
63ExitStatus interior_point(
64 const InteriorPointMatrixCallbacks<Scalar>& matrix_callbacks,
65 std::span<std::function<
bool(
const IterationInfo<Scalar>& info)>>
67 const Options& options,
68#ifdef SLEIPNIR_ENABLE_BOUND_PROJECTION
69 const Eigen::ArrayX<bool>& bound_constraint_mask,
71 Eigen::Vector<Scalar, Eigen::Dynamic>& x) {
72 using DenseVector = Eigen::Vector<Scalar, Eigen::Dynamic>;
75 DenseVector::Ones(matrix_callbacks.num_inequality_constraints);
76 DenseVector y = DenseVector::Zero(matrix_callbacks.num_equality_constraints);
78 DenseVector::Ones(matrix_callbacks.num_inequality_constraints);
81 return interior_point(matrix_callbacks, iteration_callbacks, options,
82#ifdef SLEIPNIR_ENABLE_BOUND_PROJECTION
83 bound_constraint_mask,
118template <
typename Scalar>
119ExitStatus interior_point(
120 const InteriorPointMatrixCallbacks<Scalar>& matrix_callbacks,
121 std::span<std::function<
bool(
const IterationInfo<Scalar>& info)>>
123 const Options& options,
124#ifdef SLEIPNIR_ENABLE_BOUND_PROJECTION
125 const Eigen::ArrayX<bool>& bound_constraint_mask,
127 Eigen::Vector<Scalar, Eigen::Dynamic>& x,
128 Eigen::Vector<Scalar, Eigen::Dynamic>& s,
129 Eigen::Vector<Scalar, Eigen::Dynamic>& y,
130 Eigen::Vector<Scalar, Eigen::Dynamic>& z, Scalar& μ) {
131 using DenseVector = Eigen::Vector<Scalar, Eigen::Dynamic>;
132 using SparseMatrix = Eigen::SparseMatrix<Scalar>;
133 using SparseVector = Eigen::SparseVector<Scalar>;
149 const auto solve_start_time = std::chrono::steady_clock::now();
151 gch::small_vector<SolveProfiler> solve_profilers;
152 solve_profilers.emplace_back(
"solver");
153 solve_profilers.emplace_back(
" ↳ setup");
154 solve_profilers.emplace_back(
" ↳ iteration");
155 solve_profilers.emplace_back(
" ↳ feasibility ✓");
156 solve_profilers.emplace_back(
" ↳ iter callbacks");
157 solve_profilers.emplace_back(
" ↳ KKT matrix build");
158 solve_profilers.emplace_back(
" ↳ KKT matrix decomp");
159 solve_profilers.emplace_back(
" ↳ KKT system solve");
160 solve_profilers.emplace_back(
" ↳ line search");
161 solve_profilers.emplace_back(
" ↳ SOC");
162 solve_profilers.emplace_back(
" ↳ next iter prep");
163 solve_profilers.emplace_back(
" ↳ f(x)");
164 solve_profilers.emplace_back(
" ↳ ∇f(x)");
165 solve_profilers.emplace_back(
" ↳ ∇²ₓₓL");
166 solve_profilers.emplace_back(
" ↳ cₑ(x)");
167 solve_profilers.emplace_back(
" ↳ ∂cₑ/∂x");
168 solve_profilers.emplace_back(
" ↳ cᵢ(x)");
169 solve_profilers.emplace_back(
" ↳ ∂cᵢ/∂x");
171 auto& solver_prof = solve_profilers[0];
172 auto& setup_prof = solve_profilers[1];
173 auto& inner_iter_prof = solve_profilers[2];
174 auto& feasibility_check_prof = solve_profilers[3];
175 auto& iter_callbacks_prof = solve_profilers[4];
176 auto& kkt_matrix_build_prof = solve_profilers[5];
177 auto& kkt_matrix_decomp_prof = solve_profilers[6];
178 auto& kkt_system_solve_prof = solve_profilers[7];
179 auto& line_search_prof = solve_profilers[8];
180 auto& soc_prof = solve_profilers[9];
181 auto& next_iter_prep_prof = solve_profilers[10];
184#ifndef SLEIPNIR_DISABLE_DIAGNOSTICS
185 auto& f_prof = solve_profilers[11];
186 auto& g_prof = solve_profilers[12];
187 auto& H_prof = solve_profilers[13];
188 auto& c_e_prof = solve_profilers[14];
189 auto& A_e_prof = solve_profilers[15];
190 auto& c_i_prof = solve_profilers[16];
191 auto& A_i_prof = solve_profilers[17];
193 InteriorPointMatrixCallbacks<Scalar> matrices{
194 matrix_callbacks.num_decision_variables,
195 matrix_callbacks.num_equality_constraints,
196 matrix_callbacks.num_inequality_constraints,
197 [&](
const DenseVector& x) -> Scalar {
198 ScopedProfiler prof{f_prof};
199 return matrix_callbacks.f(x);
201 [&](
const DenseVector& x) -> SparseVector {
202 ScopedProfiler prof{g_prof};
203 return matrix_callbacks.g(x);
205 [&](
const DenseVector& x,
const DenseVector& y,
206 const DenseVector& z) -> SparseMatrix {
207 ScopedProfiler prof{H_prof};
208 return matrix_callbacks.H(x, y, z);
210 [&](
const DenseVector& x) -> DenseVector {
211 ScopedProfiler prof{c_e_prof};
212 return matrix_callbacks.c_e(x);
214 [&](
const DenseVector& x) -> SparseMatrix {
215 ScopedProfiler prof{A_e_prof};
216 return matrix_callbacks.A_e(x);
218 [&](
const DenseVector& x) -> DenseVector {
219 ScopedProfiler prof{c_i_prof};
220 return matrix_callbacks.c_i(x);
222 [&](
const DenseVector& x) -> SparseMatrix {
223 ScopedProfiler prof{A_i_prof};
224 return matrix_callbacks.A_i(x);
227 const auto& matrices = matrix_callbacks;
233 Scalar f = matrices.f(x);
234 SparseVector g = matrices.g(x);
235 SparseMatrix H = matrices.H(x, y, z);
236 DenseVector c_e = matrices.c_e(x);
237 SparseMatrix A_e = matrices.A_e(x);
238 DenseVector c_i = matrices.c_i(x);
239 SparseMatrix A_i = matrices.A_i(x);
242 slp_assert(g.rows() == matrices.num_decision_variables);
243 slp_assert(H.rows() == matrices.num_decision_variables);
244 slp_assert(H.cols() == matrices.num_decision_variables);
245 slp_assert(c_e.rows() == matrices.num_equality_constraints);
246 slp_assert(A_e.rows() == matrices.num_equality_constraints);
247 slp_assert(A_e.cols() == matrices.num_decision_variables);
248 slp_assert(c_i.rows() == matrices.num_inequality_constraints);
249 slp_assert(A_i.rows() == matrices.num_inequality_constraints);
250 slp_assert(A_i.cols() == matrices.num_decision_variables);
253 if (matrices.num_equality_constraints > matrices.num_decision_variables) {
254 if (options.diagnostics) {
255 print_too_few_dofs_error(c_e);
258 return ExitStatus::TOO_FEW_DOFS;
262 if (!isfinite(f) || !all_finite(g) || !all_finite(H) || !c_e.allFinite() ||
263 !all_finite(A_e) || !c_i.allFinite() || !all_finite(A_i)) {
264 return ExitStatus::NONFINITE_INITIAL_GUESS;
267#ifdef SLEIPNIR_ENABLE_BOUND_PROJECTION
269 s = bound_constraint_mask.select(c_i, s);
275 const Scalar μ_min = Scalar(options.tolerance) / Scalar(10);
278 constexpr Scalar τ_min(0.99);
283 Filter<Scalar> filter;
287 auto update_barrier_parameter_and_reset_filter = [&] {
289 constexpr Scalar κ_μ(0.2);
293 constexpr Scalar θ_μ(1.5);
301 μ = std::max(μ_min, std::min(κ_μ * μ, pow(μ, θ_μ)));
308 τ = std::max(τ_min, Scalar(1) - μ);
315 gch::small_vector<Eigen::Triplet<Scalar>> triplets;
317 RegularizedLDLT<Scalar> solver{matrices.num_decision_variables,
318 matrices.num_equality_constraints};
321 constexpr Scalar α_reduction_factor(0.5);
322 constexpr Scalar α_min(1e-7);
324 int full_step_rejected_counter = 0;
327 Scalar E_0 = std::numeric_limits<Scalar>::infinity();
332 scope_exit exit{[&] {
333 if (options.diagnostics) {
335 if (iterations > 0) {
336 print_bottom_iteration_diagnostics();
338 print_solver_diagnostics(solve_profilers);
342 while (E_0 > Scalar(options.tolerance)) {
343 ScopedProfiler inner_iter_profiler{inner_iter_prof};
344 ScopedProfiler feasibility_check_profiler{feasibility_check_prof};
347 if (is_equality_locally_infeasible(A_e, c_e)) {
348 if (options.diagnostics) {
349 print_c_e_local_infeasibility_error(c_e);
352 return ExitStatus::LOCALLY_INFEASIBLE;
356 if (is_inequality_locally_infeasible(A_i, c_i)) {
357 if (options.diagnostics) {
358 print_c_i_local_infeasibility_error(c_i);
361 return ExitStatus::LOCALLY_INFEASIBLE;
365 if (x.template lpNorm<Eigen::Infinity>() > Scalar(1e10) || !x.allFinite() ||
366 s.template lpNorm<Eigen::Infinity>() > Scalar(1e10) || !s.allFinite()) {
367 return ExitStatus::DIVERGING_ITERATES;
370 feasibility_check_profiler.stop();
371 ScopedProfiler iter_callbacks_profiler{iter_callbacks_prof};
374 for (
const auto& callback : iteration_callbacks) {
375 if (callback({iterations, x, s, y, z, g, H, A_e, A_i})) {
376 return ExitStatus::CALLBACK_REQUESTED_STOP;
380 iter_callbacks_profiler.stop();
381 ScopedProfiler kkt_matrix_build_profiler{kkt_matrix_build_prof};
386 const SparseMatrix Σ{s.cwiseInverse().asDiagonal() * z.asDiagonal()};
392 const SparseMatrix top_left =
393 H + (A_i.transpose() * Σ * A_i).
template triangularView<Eigen::Lower>();
395 triplets.reserve(top_left.nonZeros() + A_e.nonZeros());
396 append_as_triplets(triplets, 0, 0, {top_left, A_e});
398 matrices.num_decision_variables + matrices.num_equality_constraints,
399 matrices.num_decision_variables + matrices.num_equality_constraints);
400 lhs.setFromSortedTriplets(triplets.begin(), triplets.end());
404 DenseVector rhs{x.rows() + y.rows()};
405 rhs.segment(0, x.rows()) =
406 -g + A_e.transpose() * y +
407 A_i.transpose() * (-Σ * c_i + μ * s.cwiseInverse() + z);
408 rhs.segment(x.rows(), y.rows()) = -c_e;
410 kkt_matrix_build_profiler.stop();
411 ScopedProfiler kkt_matrix_decomp_profiler{kkt_matrix_decomp_prof};
422 if (solver.compute(lhs).info() != Eigen::Success) [[unlikely]] {
423 return ExitStatus::FACTORIZATION_FAILED;
426 kkt_matrix_decomp_profiler.stop();
427 ScopedProfiler kkt_system_solve_profiler{kkt_system_solve_prof};
429 auto compute_step = [&](Step& step) {
432 DenseVector p = solver.solve(rhs);
433 step.p_x = p.segment(0, x.rows());
434 step.p_y = -p.segment(x.rows(), y.rows());
438 step.p_s = c_i - s + A_i * step.p_x;
439 step.p_z = -Σ * c_i + μ * s.cwiseInverse() - Σ * A_i * step.p_x;
443 kkt_system_solve_profiler.stop();
444 ScopedProfiler line_search_profiler{line_search_prof};
447 α_max = fraction_to_the_boundary_rule<Scalar>(s, step.p_s, τ);
452 return ExitStatus::LINE_SEARCH_FAILED;
456 α_z = fraction_to_the_boundary_rule<Scalar>(z, step.p_z, τ);
460 DenseVector trial_x = x + α * step.p_x;
461 DenseVector trial_y = y + α_z * step.p_y;
462 DenseVector trial_z = z + α_z * step.p_z;
464 Scalar trial_f = matrices.f(trial_x);
465 DenseVector trial_c_e = matrices.c_e(trial_x);
466 DenseVector trial_c_i = matrices.c_i(trial_x);
470 if (!isfinite(trial_f) || !trial_c_e.allFinite() ||
471 !trial_c_i.allFinite()) {
473 α *= α_reduction_factor;
476 return ExitStatus::LINE_SEARCH_FAILED;
482 if (options.feasible_ipm && c_i.cwiseGreater(Scalar(0)).all()) {
489 trial_s = s + α * step.p_s;
493 if (filter.try_add(FilterEntry{trial_f, trial_s, trial_c_e, trial_c_i, μ},
499 Scalar prev_constraint_violation =
500 c_e.template lpNorm<1>() + (c_i - s).
template lpNorm<1>();
501 Scalar next_constraint_violation =
502 trial_c_e.template lpNorm<1>() +
503 (trial_c_i - trial_s).
template lpNorm<1>();
510 next_constraint_violation >= prev_constraint_violation) {
512 auto soc_step = step;
515 Scalar α_z_soc = α_z;
516 DenseVector c_e_soc = c_e;
518 bool step_acceptable =
false;
519 for (
int soc_iteration = 0; soc_iteration < 5 && !step_acceptable;
521 ScopedProfiler soc_profiler{soc_prof};
523 scope_exit soc_exit{[&] {
526 if (options.diagnostics) {
527 print_iteration_diagnostics(
529 step_acceptable ? IterationType::ACCEPTED_SOC
530 : IterationType::REJECTED_SOC,
531 soc_profiler.current_duration(),
532 kkt_error<Scalar, KKTErrorType::INF_NORM_SCALED>(
533 g, A_e, trial_c_e, A_i, trial_c_i, trial_s, trial_y,
536 trial_c_e.template lpNorm<1>() +
537 (trial_c_i - trial_s).template lpNorm<1>(),
538 trial_s.dot(trial_z), μ, solver.hessian_regularization(),
539 α_soc, Scalar(1), α_reduction_factor, α_z_soc);
549 c_e_soc = α_soc * c_e_soc + trial_c_e;
550 rhs.bottomRows(y.rows()) = -c_e_soc;
553 compute_step(soc_step);
556 α_soc = fraction_to_the_boundary_rule<Scalar>(s, soc_step.p_s, τ);
557 trial_x = x + α_soc * soc_step.p_x;
558 trial_s = s + α_soc * soc_step.p_s;
561 α_z_soc = fraction_to_the_boundary_rule<Scalar>(z, soc_step.p_z, τ);
562 trial_y = y + α_z_soc * soc_step.p_y;
563 trial_z = z + α_z_soc * soc_step.p_z;
565 trial_f = matrices.f(trial_x);
566 trial_c_e = matrices.c_e(trial_x);
567 trial_c_i = matrices.c_i(trial_x);
570 constexpr Scalar κ_soc(0.99);
574 next_constraint_violation =
575 trial_c_e.template lpNorm<1>() +
576 (trial_c_i - trial_s).
template lpNorm<1>();
577 if (next_constraint_violation > κ_soc * prev_constraint_violation) {
583 FilterEntry{trial_f, trial_s, trial_c_e, trial_c_i, μ}, α)) {
587 step_acceptable =
true;
591 if (step_acceptable) {
601 ++full_step_rejected_counter;
608 if (full_step_rejected_counter >= 4 &&
609 filter.max_constraint_violation >
610 filter.back().constraint_violation / Scalar(10)) {
611 filter.max_constraint_violation *= Scalar(0.1);
617 α *= α_reduction_factor;
622 Scalar current_kkt_error = kkt_error<Scalar, KKTErrorType::ONE_NORM>(
623 g, A_e, c_e, A_i, c_i, s, y, z, μ);
625 trial_x = x + α_max * step.p_x;
626 trial_s = s + α_max * step.p_s;
628 trial_y = y + α_z * step.p_y;
629 trial_z = z + α_z * step.p_z;
631 trial_c_e = matrices.c_e(trial_x);
632 trial_c_i = matrices.c_i(trial_x);
634 Scalar next_kkt_error = kkt_error<Scalar, KKTErrorType::ONE_NORM>(
635 matrices.g(trial_x), matrices.A_e(trial_x), matrices.c_e(trial_x),
636 matrices.A_i(trial_x), trial_c_i, trial_s, trial_y, trial_z, μ);
639 if (next_kkt_error <= Scalar(0.999) * current_kkt_error) {
646 return ExitStatus::LINE_SEARCH_FAILED;
650 line_search_profiler.stop();
654 full_step_rejected_counter = 0;
679 for (
int row = 0; row < z.rows(); ++row) {
680 constexpr Scalar κ_Σ(1e10);
682 std::clamp(z[row], Scalar(1) / κ_Σ * μ / s[row], κ_Σ * μ / s[row]);
687 A_e = matrices.A_e(x);
688 A_i = matrices.A_i(x);
690 H = matrices.H(x, y, z);
692 ScopedProfiler next_iter_prep_profiler{next_iter_prep_prof};
694 c_e = matrices.c_e(x);
695 c_i = matrices.c_i(x);
698 E_0 = kkt_error<Scalar, KKTErrorType::INF_NORM_SCALED>(
699 g, A_e, c_e, A_i, c_i, s, y, z, Scalar(0));
702 if (E_0 > Scalar(options.tolerance)) {
704 constexpr Scalar κ_ε(10);
708 Scalar E_μ = kkt_error<Scalar, KKTErrorType::INF_NORM_SCALED>(
709 g, A_e, c_e, A_i, c_i, s, y, z, μ);
710 while (μ > μ_min && E_μ <= κ_ε * μ) {
711 update_barrier_parameter_and_reset_filter();
712 E_μ = kkt_error<Scalar, KKTErrorType::INF_NORM_SCALED>(g, A_e, c_e, A_i,
717 next_iter_prep_profiler.stop();
718 inner_iter_profiler.stop();
720 if (options.diagnostics) {
721 print_iteration_diagnostics(
722 iterations, IterationType::NORMAL,
723 inner_iter_profiler.current_duration(), E_0, f,
724 c_e.template lpNorm<1>() + (c_i - s).template lpNorm<1>(), s.dot(z),
725 μ, solver.hessian_regularization(), α, α_max, α_reduction_factor,
732 if (iterations >= options.max_iterations) {
733 return ExitStatus::MAX_ITERATIONS_EXCEEDED;
737 if (std::chrono::steady_clock::now() - solve_start_time > options.timeout) {
738 return ExitStatus::TIMEOUT;
742 return ExitStatus::SUCCESS;
745extern template SLEIPNIR_DLLEXPORT ExitStatus
746interior_point(
const InteriorPointMatrixCallbacks<double>& matrix_callbacks,
747 std::span<std::function<
bool(
const IterationInfo<double>& info)>>
749 const Options& options,
750#ifdef SLEIPNIR_ENABLE_BOUND_PROJECTION
751 const Eigen::ArrayX<bool>& bound_constraint_mask,
753 Eigen::Vector<double, Eigen::Dynamic>& x);