12#include <Eigen/SparseCore>
13#include <gch/small_vector.hpp>
15#include "sleipnir/optimization/solver/exit_status.hpp"
16#include "sleipnir/optimization/solver/iteration_info.hpp"
17#include "sleipnir/optimization/solver/options.hpp"
18#include "sleipnir/optimization/solver/sqp_matrix_callbacks.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/filter.hpp"
22#include "sleipnir/optimization/solver/util/is_locally_infeasible.hpp"
23#include "sleipnir/optimization/solver/util/kkt_error.hpp"
24#include "sleipnir/optimization/solver/util/regularized_ldlt.hpp"
25#include "sleipnir/util/assert.hpp"
26#include "sleipnir/util/print_diagnostics.hpp"
27#include "sleipnir/util/profiler.hpp"
28#include "sleipnir/util/scope_exit.hpp"
29#include "sleipnir/util/symbol_exports.hpp"
55template <
typename Scalar>
56ExitStatus sqp(
const SQPMatrixCallbacks<Scalar>& matrix_callbacks,
57 std::span<std::function<
bool(
const IterationInfo<Scalar>& info)>>
59 const Options& options,
60 Eigen::Vector<Scalar, Eigen::Dynamic>& x) {
61 using DenseVector = Eigen::Vector<Scalar, Eigen::Dynamic>;
63 DenseVector y = DenseVector::Zero(matrix_callbacks.num_equality_constraints);
65 return sqp(matrix_callbacks, iteration_callbacks, options, x, y);
90template <
typename Scalar>
91ExitStatus sqp(
const SQPMatrixCallbacks<Scalar>& matrix_callbacks,
92 std::span<std::function<
bool(
const IterationInfo<Scalar>& info)>>
94 const Options& options, Eigen::Vector<Scalar, Eigen::Dynamic>& x,
95 Eigen::Vector<Scalar, Eigen::Dynamic>& y) {
96 using DenseVector = Eigen::Vector<Scalar, Eigen::Dynamic>;
97 using SparseMatrix = Eigen::SparseMatrix<Scalar>;
98 using SparseVector = Eigen::SparseVector<Scalar>;
110 const auto solve_start_time = std::chrono::steady_clock::now();
112 gch::small_vector<SolveProfiler> solve_profilers;
113 solve_profilers.emplace_back(
"solver");
114 solve_profilers.emplace_back(
" ↳ setup");
115 solve_profilers.emplace_back(
" ↳ iteration");
116 solve_profilers.emplace_back(
" ↳ feasibility ✓");
117 solve_profilers.emplace_back(
" ↳ iter callbacks");
118 solve_profilers.emplace_back(
" ↳ KKT matrix build");
119 solve_profilers.emplace_back(
" ↳ KKT matrix decomp");
120 solve_profilers.emplace_back(
" ↳ KKT system solve");
121 solve_profilers.emplace_back(
" ↳ line search");
122 solve_profilers.emplace_back(
" ↳ SOC");
123 solve_profilers.emplace_back(
" ↳ next iter prep");
124 solve_profilers.emplace_back(
" ↳ f(x)");
125 solve_profilers.emplace_back(
" ↳ ∇f(x)");
126 solve_profilers.emplace_back(
" ↳ ∇²ₓₓL");
127 solve_profilers.emplace_back(
" ↳ cₑ(x)");
128 solve_profilers.emplace_back(
" ↳ ∂cₑ/∂x");
130 auto& solver_prof = solve_profilers[0];
131 auto& setup_prof = solve_profilers[1];
132 auto& inner_iter_prof = solve_profilers[2];
133 auto& feasibility_check_prof = solve_profilers[3];
134 auto& iter_callbacks_prof = solve_profilers[4];
135 auto& kkt_matrix_build_prof = solve_profilers[5];
136 auto& kkt_matrix_decomp_prof = solve_profilers[6];
137 auto& kkt_system_solve_prof = solve_profilers[7];
138 auto& line_search_prof = solve_profilers[8];
139 auto& soc_prof = solve_profilers[9];
140 auto& next_iter_prep_prof = solve_profilers[10];
143#ifndef SLEIPNIR_DISABLE_DIAGNOSTICS
144 auto& f_prof = solve_profilers[11];
145 auto& g_prof = solve_profilers[12];
146 auto& H_prof = solve_profilers[13];
147 auto& c_e_prof = solve_profilers[14];
148 auto& A_e_prof = solve_profilers[15];
150 SQPMatrixCallbacks<Scalar> matrices{
151 matrix_callbacks.num_decision_variables,
152 matrix_callbacks.num_equality_constraints,
153 [&](
const DenseVector& x) -> Scalar {
154 ScopedProfiler prof{f_prof};
155 return matrix_callbacks.f(x);
157 [&](
const DenseVector& x) -> SparseVector {
158 ScopedProfiler prof{g_prof};
159 return matrix_callbacks.g(x);
161 [&](
const DenseVector& x,
const DenseVector& y) -> SparseMatrix {
162 ScopedProfiler prof{H_prof};
163 return matrix_callbacks.H(x, y);
165 [&](
const DenseVector& x) -> DenseVector {
166 ScopedProfiler prof{c_e_prof};
167 return matrix_callbacks.c_e(x);
169 [&](
const DenseVector& x) -> SparseMatrix {
170 ScopedProfiler prof{A_e_prof};
171 return matrix_callbacks.A_e(x);
174 const auto& matrices = matrix_callbacks;
180 Scalar f = matrices.f(x);
181 SparseVector g = matrices.g(x);
182 SparseMatrix H = matrices.H(x, y);
183 DenseVector c_e = matrices.c_e(x);
184 SparseMatrix A_e = matrices.A_e(x);
187 slp_assert(g.rows() == matrices.num_decision_variables);
188 slp_assert(H.rows() == matrices.num_decision_variables);
189 slp_assert(H.cols() == matrices.num_decision_variables);
190 slp_assert(c_e.rows() == matrices.num_equality_constraints);
191 slp_assert(A_e.rows() == matrices.num_equality_constraints);
192 slp_assert(A_e.cols() == matrices.num_decision_variables);
195 if (matrices.num_equality_constraints > matrices.num_decision_variables) {
196 if (options.diagnostics) {
197 print_too_few_dofs_error(c_e);
200 return ExitStatus::TOO_FEW_DOFS;
204 if (!isfinite(f) || !all_finite(g) || !all_finite(H) || !c_e.allFinite() ||
206 return ExitStatus::NONFINITE_INITIAL_GUESS;
211 Filter<Scalar> filter;
214 gch::small_vector<Eigen::Triplet<Scalar>> triplets;
216 RegularizedLDLT<Scalar> solver{matrices.num_decision_variables,
217 matrices.num_equality_constraints};
220 constexpr Scalar α_reduction_factor(0.5);
221 constexpr Scalar α_min(1e-7);
223 int full_step_rejected_counter = 0;
226 Scalar E_0 = std::numeric_limits<Scalar>::infinity();
231 scope_exit exit{[&] {
232 if (options.diagnostics) {
234 if (iterations > 0) {
235 print_bottom_iteration_diagnostics();
237 print_solver_diagnostics(solve_profilers);
241 while (E_0 > Scalar(options.tolerance)) {
242 ScopedProfiler inner_iter_profiler{inner_iter_prof};
243 ScopedProfiler feasibility_check_profiler{feasibility_check_prof};
246 if (is_equality_locally_infeasible(A_e, c_e)) {
247 if (options.diagnostics) {
248 print_c_e_local_infeasibility_error(c_e);
251 return ExitStatus::LOCALLY_INFEASIBLE;
255 if (x.template lpNorm<Eigen::Infinity>() > Scalar(1e10) || !x.allFinite()) {
256 return ExitStatus::DIVERGING_ITERATES;
259 feasibility_check_profiler.stop();
260 ScopedProfiler iter_callbacks_profiler{iter_callbacks_prof};
263 for (
const auto& callback : iteration_callbacks) {
264 if (callback({iterations, x, {}, y, {}, g, H, A_e, {}})) {
265 return ExitStatus::CALLBACK_REQUESTED_STOP;
269 iter_callbacks_profiler.stop();
270 ScopedProfiler kkt_matrix_build_profiler{kkt_matrix_build_prof};
277 triplets.reserve(H.nonZeros() + A_e.nonZeros());
278 append_as_triplets(triplets, 0, 0, {H, A_e});
280 matrices.num_decision_variables + matrices.num_equality_constraints,
281 matrices.num_decision_variables + matrices.num_equality_constraints);
282 lhs.setFromSortedTriplets(triplets.begin(), triplets.end());
286 DenseVector rhs{x.rows() + y.rows()};
287 rhs.segment(0, x.rows()) = -g + A_e.transpose() * y;
288 rhs.segment(x.rows(), y.rows()) = -c_e;
290 kkt_matrix_build_profiler.stop();
291 ScopedProfiler kkt_matrix_decomp_profiler{kkt_matrix_decomp_prof};
294 constexpr Scalar α_max(1);
301 if (solver.compute(lhs).info() != Eigen::Success) [[unlikely]] {
302 return ExitStatus::FACTORIZATION_FAILED;
305 kkt_matrix_decomp_profiler.stop();
306 ScopedProfiler kkt_system_solve_profiler{kkt_system_solve_prof};
308 auto compute_step = [&](Step& step) {
311 DenseVector p = solver.solve(rhs);
312 step.p_x = p.segment(0, x.rows());
313 step.p_y = -p.segment(x.rows(), y.rows());
317 kkt_system_solve_profiler.stop();
318 ScopedProfiler line_search_profiler{line_search_prof};
324 DenseVector trial_x = x + α * step.p_x;
325 DenseVector trial_y = y + α * step.p_y;
327 Scalar trial_f = matrices.f(trial_x);
328 DenseVector trial_c_e = matrices.c_e(trial_x);
332 if (!isfinite(trial_f) || !trial_c_e.allFinite()) {
334 α *= α_reduction_factor;
337 return ExitStatus::LINE_SEARCH_FAILED;
343 if (filter.try_add(FilterEntry{trial_f, trial_c_e}, α)) {
348 Scalar prev_constraint_violation = c_e.template lpNorm<1>();
349 Scalar next_constraint_violation = trial_c_e.template lpNorm<1>();
356 next_constraint_violation >= prev_constraint_violation) {
358 auto soc_step = step;
361 DenseVector c_e_soc = c_e;
363 bool step_acceptable =
false;
364 for (
int soc_iteration = 0; soc_iteration < 5 && !step_acceptable;
366 ScopedProfiler soc_profiler{soc_prof};
368 scope_exit soc_exit{[&] {
371 if (options.diagnostics) {
372 print_iteration_diagnostics(
374 step_acceptable ? IterationType::ACCEPTED_SOC
375 : IterationType::REJECTED_SOC,
376 soc_profiler.current_duration(),
377 kkt_error<Scalar, KKTErrorType::INF_NORM_SCALED>(
378 g, A_e, trial_c_e, trial_y),
379 trial_f, trial_c_e.template lpNorm<1>(), Scalar(0), Scalar(0),
380 solver.hessian_regularization(), α_soc, Scalar(1),
381 α_reduction_factor, Scalar(1));
391 c_e_soc = α_soc * c_e_soc + trial_c_e;
392 rhs.bottomRows(y.rows()) = -c_e_soc;
395 compute_step(soc_step);
397 trial_x = x + α_soc * soc_step.p_x;
398 trial_y = y + α_soc * soc_step.p_y;
400 trial_f = matrices.f(trial_x);
401 trial_c_e = matrices.c_e(trial_x);
404 constexpr Scalar κ_soc(0.99);
408 next_constraint_violation = trial_c_e.template lpNorm<1>();
409 if (next_constraint_violation > κ_soc * prev_constraint_violation) {
414 if (filter.try_add(FilterEntry{trial_f, trial_c_e}, α)) {
417 step_acceptable =
true;
421 if (step_acceptable) {
431 ++full_step_rejected_counter;
438 if (full_step_rejected_counter >= 4 &&
439 filter.max_constraint_violation >
440 filter.back().constraint_violation / Scalar(10)) {
441 filter.max_constraint_violation *= Scalar(0.1);
447 α *= α_reduction_factor;
452 Scalar current_kkt_error =
453 kkt_error<Scalar, KKTErrorType::ONE_NORM>(g, A_e, c_e, y);
455 trial_x = x + α_max * step.p_x;
456 trial_y = y + α_max * step.p_y;
458 trial_c_e = matrices.c_e(trial_x);
460 Scalar next_kkt_error = kkt_error<Scalar, KKTErrorType::ONE_NORM>(
461 matrices.g(trial_x), matrices.A_e(trial_x), trial_c_e, trial_y);
464 if (next_kkt_error <= Scalar(0.999) * current_kkt_error) {
471 return ExitStatus::LINE_SEARCH_FAILED;
475 line_search_profiler.stop();
479 full_step_rejected_counter = 0;
489 A_e = matrices.A_e(x);
491 H = matrices.H(x, y);
493 ScopedProfiler next_iter_prep_profiler{next_iter_prep_prof};
495 c_e = matrices.c_e(x);
498 E_0 = kkt_error<Scalar, KKTErrorType::INF_NORM_SCALED>(g, A_e, c_e, y);
500 next_iter_prep_profiler.stop();
501 inner_iter_profiler.stop();
503 if (options.diagnostics) {
504 print_iteration_diagnostics(iterations, IterationType::NORMAL,
505 inner_iter_profiler.current_duration(), E_0,
506 f, c_e.template lpNorm<1>(), Scalar(0),
507 Scalar(0), solver.hessian_regularization(), α,
508 α_max, α_reduction_factor, α);
514 if (iterations >= options.max_iterations) {
515 return ExitStatus::MAX_ITERATIONS_EXCEEDED;
519 if (std::chrono::steady_clock::now() - solve_start_time > options.timeout) {
520 return ExitStatus::TIMEOUT;
524 return ExitStatus::SUCCESS;
527extern template SLEIPNIR_DLLEXPORT ExitStatus
528sqp(
const SQPMatrixCallbacks<double>& matrix_callbacks,
529 std::span<std::function<
bool(
const IterationInfo<double>& info)>>
531 const Options& options, Eigen::Vector<double, Eigen::Dynamic>& x);