24template <
typename Scalar>
29 const Options& options,
bool in_feasibility_restoration,
30#ifdef SLEIPNIR_ENABLE_BOUND_PROJECTION
31 const Eigen::ArrayX<bool>& bound_constraint_mask,
33 Eigen::Vector<Scalar, Eigen::Dynamic>& x,
34 Eigen::Vector<Scalar, Eigen::Dynamic>& s,
35 Eigen::Vector<Scalar, Eigen::Dynamic>& y,
36 Eigen::Vector<Scalar, Eigen::Dynamic>& z, Scalar& μ);
46template <
typename Scalar>
47std::tuple<Eigen::Vector<Scalar, Eigen::Dynamic>,
48 Eigen::Vector<Scalar, Eigen::Dynamic>>
49compute_p_n(
const Eigen::Vector<Scalar, Eigen::Dynamic>& c, Scalar ρ,
81 using DenseVector = Eigen::Vector<Scalar, Eigen::Dynamic>;
85 DenseVector
p{c.rows()};
86 DenseVector n{c.rows()};
87 for (
int row = 0; row <
p.rows(); ++row) {
89 Scalar _b = ρ * c[row] - μ;
90 Scalar _c = -μ * c[row] / Scalar(2);
92 n[row] = (-_b +
sqrt(_b * _b - Scalar(4) * _a * _c)) / (Scalar(2) * _a);
93 p[row] = c[row] + n[row];
96 return {std::move(
p), std::move(n)};
115template <
typename Scalar>
117 const SQPMatrixCallbacks<Scalar>& matrix_callbacks,
118 std::span<std::function<
bool(
const IterationInfo<Scalar>& info)>>
120 const Options& options, Eigen::Vector<Scalar, Eigen::Dynamic>& x,
121 Eigen::Vector<Scalar, Eigen::Dynamic>& y) {
138 using DenseVector = Eigen::Vector<Scalar, Eigen::Dynamic>;
139 using DiagonalMatrix = Eigen::DiagonalMatrix<Scalar, Eigen::Dynamic>;
140 using SparseMatrix = Eigen::SparseMatrix<Scalar>;
141 using SparseVector = Eigen::SparseVector<Scalar>;
145 const auto& matrices = matrix_callbacks;
146 const auto& num_vars = matrices.num_decision_variables;
147 const auto& num_eq = matrices.num_equality_constraints;
149 constexpr Scalar ρ(1e3);
150 const Scalar μ(options.tolerance / 10.0);
151 const Scalar ζ =
sqrt(μ);
153 const DenseVector c_e = matrices.c_e(x);
156 const auto [p_e_0, n_e_0] =
compute_p_n(c_e, ρ, μ);
159 const DiagonalMatrix D_r =
160 x.cwiseSquare().cwiseInverse().cwiseMin(Scalar(1)).asDiagonal();
162 DenseVector fr_x{num_vars + 2 * num_eq};
163 fr_x << x, p_e_0, n_e_0;
165 DenseVector fr_s = DenseVector::Ones(2 * num_eq);
167 DenseVector fr_y = DenseVector::Zero(num_eq);
169 DenseVector fr_z{2 * num_eq};
170 fr_z << μ * p_e_0.cwiseInverse(), μ * n_e_0.cwiseInverse();
172 Scalar fr_μ = std::max(μ, c_e.template lpNorm<Eigen::Infinity>());
175 static_cast<int>(fr_x.rows()),
176 static_cast<int>(fr_y.rows()),
177 static_cast<int>(fr_z.rows()),
178 [&](
const DenseVector& x_p) -> Scalar {
179 auto x = x_p.segment(0, num_vars);
186 return ρ * x_p.segment(num_vars, 2 * num_eq).array().sum() +
187 ζ / Scalar(2) * diff.transpose() * D_r * diff;
189 [&](
const DenseVector& x_p) -> SparseVector {
190 auto x = x_p.segment(0, num_vars);
197 DenseVector g{x_p.rows()};
198 g.segment(0, num_vars) = ζ * D_r * (x - x_r);
199 g.segment(num_vars, 2 * num_eq).setConstant(ρ);
200 return g.sparseView();
202 [&](
const DenseVector& x_p,
const DenseVector& y_p,
203 [[maybe_unused]]
const DenseVector& z_p) -> SparseMatrix {
204 auto x = x_p.segment(0, num_vars);
213 triplets.reserve(x_p.rows());
215 SparseMatrix d2f_dx2{x_p.rows(), x_p.rows()};
216 d2f_dx2.setFromSortedTriplets(triplets.begin(), triplets.end());
221 auto H_c = matrices.H_c(x, y);
222 H_c.resize(x_p.rows(), x_p.rows());
229 return d2f_dx2 + H_c;
231 [&](
const DenseVector& x_p, [[maybe_unused]]
const DenseVector& y_p,
232 [[maybe_unused]]
const DenseVector& z_p) -> SparseMatrix {
233 return SparseMatrix{x_p.rows(), x_p.rows()};
235 [&](
const DenseVector& x_p) -> DenseVector {
236 auto x = x_p.segment(0, num_vars);
237 auto p_e = x_p.segment(num_vars, num_eq);
238 auto n_e = x_p.segment(num_vars + num_eq, num_eq);
243 return matrices.c_e(x) - p_e + n_e;
245 [&](
const DenseVector& x_p) -> SparseMatrix {
246 auto x = x_p.segment(0, num_vars);
252 SparseMatrix A_e = matrices.A_e(x);
255 triplets.reserve(A_e.nonZeros() + 2 * num_eq);
259 triplets, 0, num_vars,
260 DenseVector::Constant(num_eq, Scalar(-1)).eval());
262 triplets, 0, num_vars + num_eq,
263 DenseVector::Constant(num_eq, Scalar(1)).eval());
265 SparseMatrix A_e_p{A_e.rows(), x_p.rows()};
266 A_e_p.setFromSortedTriplets(triplets.begin(), triplets.end());
269 [&](
const DenseVector& x_p) -> DenseVector {
274 return x_p.segment(num_vars, 2 * num_eq);
276 [&](
const DenseVector& x_p) -> SparseMatrix {
283 triplets.reserve(2 * num_eq);
286 triplets, 0, num_vars,
287 DenseVector::Constant(2 * num_eq, Scalar(1)).eval());
289 SparseMatrix A_i_p{2 * num_eq, x_p.rows()};
290 A_i_p.setFromSortedTriplets(triplets.begin(), triplets.end());
296#ifdef SLEIPNIR_ENABLE_BOUND_PROJECTION
299 fr_x, fr_s, fr_y, fr_z, fr_μ);
301 x = fr_x.segment(0, x.rows());
304 auto g = matrices.g(x);
305 auto A_e = matrices.A_e(x);
335template <
typename Scalar>
340 const Options& options, Eigen::Vector<Scalar, Eigen::Dynamic>& x,
341 Eigen::Vector<Scalar, Eigen::Dynamic>& s,
342 Eigen::Vector<Scalar, Eigen::Dynamic>& y,
343 Eigen::Vector<Scalar, Eigen::Dynamic>& z, Scalar μ) {
364 using DenseVector = Eigen::Vector<Scalar, Eigen::Dynamic>;
365 using DiagonalMatrix = Eigen::DiagonalMatrix<Scalar, Eigen::Dynamic>;
366 using SparseMatrix = Eigen::SparseMatrix<Scalar>;
367 using SparseVector = Eigen::SparseVector<Scalar>;
371 const auto& matrices = matrix_callbacks;
372 const auto& num_vars = matrices.num_decision_variables;
373 const auto& num_eq = matrices.num_equality_constraints;
374 const auto& num_ineq = matrices.num_inequality_constraints;
376 constexpr Scalar ρ(1e3);
377 const Scalar ζ =
sqrt(μ);
379 const DenseVector c_e = matrices.c_e(x);
380 const DenseVector c_i = matrices.c_i(x);
383 const auto [p_e_0, n_e_0] =
compute_p_n(c_e, ρ, μ);
384 const auto [p_i_0, n_i_0] =
compute_p_n((c_i - s).eval(), ρ, μ);
387 const DiagonalMatrix D_r =
388 x.cwiseSquare().cwiseInverse().cwiseMin(Scalar(1)).asDiagonal();
390 DenseVector fr_x{num_vars + 2 * num_eq + 2 * num_ineq};
391 fr_x << x, p_e_0, n_e_0, p_i_0, n_i_0;
393 DenseVector fr_s{s.rows() + 2 * num_eq + 2 * num_ineq};
394 fr_s.segment(0, s.rows()) = s;
395 fr_s.segment(s.rows(), 2 * num_eq + 2 * num_ineq).setOnes();
397 DenseVector fr_y = DenseVector::Zero(c_e.rows());
399 DenseVector fr_z{c_i.rows() + 2 * num_eq + 2 * num_ineq};
400 fr_z << z.cwiseMin(ρ), μ * p_e_0.cwiseInverse(), μ * n_e_0.cwiseInverse(),
401 μ * p_i_0.cwiseInverse(), μ * n_i_0.cwiseInverse();
403 Scalar fr_μ = std::max({μ, c_e.template lpNorm<Eigen::Infinity>(),
404 (c_i - s).
template lpNorm<Eigen::Infinity>()});
407 static_cast<int>(fr_x.rows()),
408 static_cast<int>(fr_y.rows()),
409 static_cast<int>(fr_z.rows()),
410 [&](
const DenseVector& x_p) -> Scalar {
411 auto x = x_p.segment(0, num_vars);
417 return ρ * x_p.segment(num_vars, 2 * num_eq + 2 * num_ineq)
420 ζ / Scalar(2) * diff.transpose() * D_r * diff;
422 [&](
const DenseVector& x_p) -> SparseVector {
423 auto x = x_p.segment(0, num_vars);
432 DenseVector g{x_p.rows()};
433 g.segment(0, num_vars) = ζ * D_r * (x - x_r);
434 g.segment(num_vars, 2 * num_eq + 2 * num_ineq).setConstant(ρ);
435 return g.sparseView();
437 [&](
const DenseVector& x_p,
const DenseVector& y_p,
438 const DenseVector& z_p) -> SparseMatrix {
439 auto x = x_p.segment(0, num_vars);
441 auto z = z_p.segment(0, num_ineq);
451 triplets.reserve(x_p.rows());
453 SparseMatrix d2f_dx2{x_p.rows(), x_p.rows()};
454 d2f_dx2.setFromSortedTriplets(triplets.begin(), triplets.end());
459 auto H_c = matrices.H_c(x, y, z);
460 H_c.resize(x_p.rows(), x_p.rows());
469 return d2f_dx2 + H_c;
471 [&](
const DenseVector& x_p, [[maybe_unused]]
const DenseVector& y_p,
472 [[maybe_unused]]
const DenseVector& z_p) -> SparseMatrix {
473 return SparseMatrix{x_p.rows(), x_p.rows()};
475 [&](
const DenseVector& x_p) -> DenseVector {
476 auto x = x_p.segment(0, num_vars);
477 auto p_e = x_p.segment(num_vars, num_eq);
478 auto n_e = x_p.segment(num_vars + num_eq, num_eq);
483 return matrices.c_e(x) - p_e + n_e;
485 [&](
const DenseVector& x_p) -> SparseMatrix {
486 auto x = x_p.segment(0, num_vars);
492 SparseMatrix A_e = matrices.A_e(x);
495 triplets.reserve(A_e.nonZeros() + 2 * num_eq);
499 triplets, 0, num_vars,
500 DenseVector::Constant(num_eq, Scalar(-1)).eval());
502 triplets, 0, num_vars + num_eq,
503 DenseVector::Constant(num_eq, Scalar(1)).eval());
505 SparseMatrix A_e_p{A_e.rows(), x_p.rows()};
506 A_e_p.setFromSortedTriplets(triplets.begin(), triplets.end());
509 [&](
const DenseVector& x_p) -> DenseVector {
510 auto x = x_p.segment(0, num_vars);
511 auto p_i = x_p.segment(num_vars + 2 * num_eq, num_ineq);
512 auto n_i = x_p.segment(num_vars + 2 * num_eq + num_ineq, num_ineq);
521 DenseVector c_i_p{c_i.rows() + 2 * num_eq + 2 * num_ineq};
522 c_i_p.segment(0, num_ineq) = matrices.c_i(x) - p_i + n_i;
523 c_i_p.segment(p_i.rows(), 2 * num_eq + 2 * num_ineq) =
524 x_p.segment(num_vars, 2 * num_eq + 2 * num_ineq);
527 [&](
const DenseVector& x_p) -> SparseMatrix {
528 auto x = x_p.segment(0, num_vars);
538 SparseMatrix A_i = matrices.A_i(x);
541 triplets.reserve(A_i.nonZeros() + 2 * num_eq + 4 * num_ineq);
548 triplets, num_ineq, num_vars,
549 DenseVector::Constant(2 * num_eq, Scalar(1)).eval());
552 DenseVector::Constant(num_ineq, Scalar(1)).asDiagonal()};
555 SparseMatrix Z_col3{2 * num_eq, num_ineq};
557 {(-I_ineq).eval(), Z_col3, I_ineq});
560 SparseMatrix Z_col4{2 * num_eq + num_ineq, num_ineq};
562 {I_ineq, Z_col4, I_ineq});
564 SparseMatrix A_i_p{2 * num_eq + 3 * num_ineq, x_p.rows()};
565 A_i_p.setFromSortedTriplets(triplets.begin(), triplets.end());
571#ifdef SLEIPNIR_ENABLE_BOUND_PROJECTION
574 fr_x, fr_s, fr_y, fr_z, fr_μ);
576 x = fr_x.segment(0, x.rows());
577 s = fr_s.segment(0, s.rows());
580 auto g = matrices.g(x);
581 auto A_e = matrices.A_e(x);
582 auto A_i = matrices.A_i(x);
584 auto [y_estimate, z_estimate] =
601#include "sleipnir/optimization/solver/interior_point.hpp"
void append_as_triplets(gch::small_vector< Eigen::Triplet< Scalar > > &triplets, int row_offset, int col_offset, std::initializer_list< Eigen::SparseMatrix< Scalar > > mats)
Appends sparse matrices to list of triplets at the given offset.
Definition append_as_triplets.hpp:22
void append_diagonal_as_triplets(gch::small_vector< Eigen::Triplet< Scalar > > &triplets, int row_offset, int col_offset, const Eigen::Vector< Scalar, Eigen::Dynamic > &diag)
Append diagonal matrix to list of triplets at the given offset.
Definition append_as_triplets.hpp:54
wpi::util::SmallVector< T > small_vector
Definition small_vector.hpp:10
Definition to_underlying.hpp:7
ExitStatus
Solver exit status. Negative values indicate failure.
Definition exit_status.hpp:14
@ CALLBACK_REQUESTED_STOP
The solver returned its solution so far after the user requested a stop.
Definition exit_status.hpp:18
@ SUCCESS
Solved the problem to the desired tolerance.
Definition exit_status.hpp:16
@ FEASIBILITY_RESTORATION_FAILED
The solver failed to reach the desired tolerance, and feasibility restoration failed to converge.
Definition exit_status.hpp:33
@ LOCALLY_INFEASIBLE
The solver determined the problem to be locally infeasible and gave up.
Definition exit_status.hpp:22
ExitStatus interior_point(const InteriorPointMatrixCallbacks< Scalar > &matrix_callbacks, std::span< std::function< bool(const IterationInfo< Scalar > &info)> > iteration_callbacks, const Options &options, bool in_feasibility_restoration, Eigen::Vector< Scalar, Eigen::Dynamic > &x, Eigen::Vector< Scalar, Eigen::Dynamic > &s, Eigen::Vector< Scalar, Eigen::Dynamic > &y, Eigen::Vector< Scalar, Eigen::Dynamic > &z, Scalar &μ)
std::tuple< Eigen::Vector< Scalar, Eigen::Dynamic >, Eigen::Vector< Scalar, Eigen::Dynamic > > compute_p_n(const Eigen::Vector< Scalar, Eigen::Dynamic > &c, Scalar ρ, Scalar μ)
Computes initial values for p and n in feasibility restoration.
Definition feasibility_restoration.hpp:49
Eigen::Vector< Scalar, Eigen::Dynamic > lagrange_multiplier_estimate(const Eigen::SparseVector< Scalar > &g, const Eigen::SparseMatrix< Scalar > &A_e)
Estimates Lagrange multipliers for SQP.
Definition lagrange_multiplier_estimate.hpp:34
Variable< Scalar > sqrt(const Variable< Scalar > &x)
sqrt() for Variables.
Definition variable.hpp:671
Matrix callbacks for the interior-point method solver.
Definition interior_point_matrix_callbacks.hpp:16
Solver iteration information exposed to an iteration callback.
Definition iteration_info.hpp:14
Solver options.
Definition options.hpp:13