8#include <Eigen/SparseCholesky>
9#include <Eigen/SparseCore>
10#include <gch/small_vector.hpp>
12#include "sleipnir/optimization/solver/util/append_as_triplets.hpp"
19template <
typename Scalar>
22 Eigen::Vector<Scalar, Eigen::Dynamic>
y;
24 Eigen::Vector<Scalar, Eigen::Dynamic>
z;
32template <
typename Scalar>
33Eigen::Vector<Scalar, Eigen::Dynamic> lagrange_multiplier_estimate(
34 const Eigen::SparseVector<Scalar>& g,
35 const Eigen::SparseMatrix<Scalar>& A_e) {
41 return Eigen::SimplicialLDLT<Eigen::SparseMatrix<Scalar>>{A_e *
54template <
typename Scalar>
55LagrangeMultiplierEstimate<Scalar> lagrange_multiplier_estimate(
56 const Eigen::SparseVector<Scalar>& g,
57 const Eigen::SparseMatrix<Scalar>& A_e,
58 const Eigen::SparseMatrix<Scalar>& A_i,
59 const Eigen::Vector<Scalar, Eigen::Dynamic>& s, Scalar μ) {
60 using DenseVector = Eigen::Vector<Scalar, Eigen::Dynamic>;
61 using SparseMatrix = Eigen::SparseMatrix<Scalar>;
86 gch::small_vector<Eigen::Triplet<Scalar>> triplets;
90 triplets.reserve(A_e.nonZeros() + A_i.nonZeros() + s.rows());
91 append_as_triplets(triplets, 0, 0, {A_e, A_i});
92 append_diagonal_as_triplets(triplets, A_e.rows(), A_i.cols(), (-s).eval());
93 SparseMatrix A_hat{A_e.rows() + A_i.rows(), A_e.cols() + s.rows()};
94 A_hat.setFromSortedTriplets(triplets.begin(), triplets.end());
97 SparseMatrix lhs = A_hat * A_hat.transpose();
101 DenseVector rhs_temp{g.rows() + s.rows()};
102 rhs_temp.segment(0, g.rows()) = g;
103 rhs_temp.segment(g.rows(), s.rows()).setConstant(-μ);
104 DenseVector rhs = A_hat * rhs_temp;
106 Eigen::SimplicialLDLT<SparseMatrix> yz_estimator{lhs};
107 DenseVector sol = yz_estimator.solve(rhs);
108 DenseVector y = sol.segment(0, A_e.rows());
109 DenseVector z = sol.segment(A_e.rows(), s.rows());
124 for (
int row = 0; row < z.rows(); ++row) {
125 constexpr Scalar κ_Σ(1e10);
126 z[row] = std::clamp(z[row], Scalar(1) / κ_Σ * μ / s[row], κ_Σ * μ / s[row]);
129 return {std::move(y), std::move(z)};
Definition lagrange_multiplier_estimate.hpp:20
Eigen::Vector< Scalar, Eigen::Dynamic > y
Equality constraint dual estimate.
Definition lagrange_multiplier_estimate.hpp:22
Eigen::Vector< Scalar, Eigen::Dynamic > z
Inequality constraint dual estimate.
Definition lagrange_multiplier_estimate.hpp:24