Sleipnir C++ API
Loading...
Searching...
No Matches
lagrange_multiplier_estimate.hpp
1// Copyright (c) Sleipnir contributors
2
3#pragma once
4
5#include <algorithm>
6#include <utility>
7
8#include <Eigen/SparseCholesky>
9#include <Eigen/SparseCore>
10#include <gch/small_vector.hpp>
11
12#include "sleipnir/optimization/solver/util/append_as_triplets.hpp"
13
14namespace slp {
15
19template <typename Scalar>
22 Eigen::Vector<Scalar, Eigen::Dynamic> y;
24 Eigen::Vector<Scalar, Eigen::Dynamic> z;
25};
26
33template <typename Scalar>
34Eigen::Vector<Scalar, Eigen::Dynamic> lagrange_multiplier_estimate(
35 const Eigen::SparseVector<Scalar>& g,
36 const Eigen::SparseMatrix<Scalar>& A_e) {
37 // Lagrange multiplier estimates
38 //
39 // ∇f − Aₑᵀy = 0
40 // Aₑᵀy = ∇f
41 // y = (AₑAₑᵀ)⁻¹Aₑ∇f
42 return Eigen::SimplicialLDLT<Eigen::SparseMatrix<Scalar>>{A_e *
43 A_e.transpose()}
44 .solve(A_e * g);
45}
46
57template <typename Scalar>
58LagrangeMultiplierEstimate<Scalar> lagrange_multiplier_estimate(
59 const Eigen::SparseVector<Scalar>& g,
60 const Eigen::SparseMatrix<Scalar>& A_e,
61 const Eigen::SparseMatrix<Scalar>& A_i,
62 const Eigen::Vector<Scalar, Eigen::Dynamic>& s, Scalar μ) {
63 using DenseVector = Eigen::Vector<Scalar, Eigen::Dynamic>;
64 using SparseMatrix = Eigen::SparseMatrix<Scalar>;
65
66 // Lagrange multiplier estimates
67 //
68 // ∇f − Aₑᵀy − Aᵢᵀz = 0
69 // Sz − μe = 0
70 //
71 // Aₑᵀy + Aᵢᵀz = ∇f
72 // −Sz = −μe
73 //
74 // [Aₑᵀ Aᵢᵀ][y] = [ ∇f]
75 // [ 0 −S ][z] [−μe]
76 //
77 // [Aₑ 0]ᵀ[y] = [ ∇f]
78 // [Aᵢ −S] [z] [−μe]
79 //
80 // Let  = [Aₑ 0]
81 // [Aᵢ −S]
82 //
83 // Âᵀ[y] = [ ∇f]
84 // [z] [−μe]
85 //
86 // [y] = (ÂÂᵀ)⁻¹Â[ ∇f]
87 // [z] [−μe]
88
89 gch::small_vector<Eigen::Triplet<Scalar>> triplets;
90
91 // Â = [Aₑ 0]
92 // [Aᵢ −S]
93 triplets.reserve(A_e.nonZeros() + A_i.nonZeros() + s.rows());
94 append_as_triplets(triplets, 0, 0, {A_e, A_i});
95 append_diagonal_as_triplets(triplets, A_e.rows(), A_i.cols(), (-s).eval());
96 SparseMatrix A_hat{A_e.rows() + A_i.rows(), A_e.cols() + s.rows()};
97 A_hat.setFromSortedTriplets(triplets.begin(), triplets.end());
98
99 // lhs = ÂÂᵀ
100 SparseMatrix lhs = A_hat * A_hat.transpose();
101
102 // rhs = Â[ ∇f]
103 // [−μe]
104 DenseVector rhs_temp{g.rows() + s.rows()};
105 rhs_temp.segment(0, g.rows()) = g;
106 rhs_temp.segment(g.rows(), s.rows()).setConstant(-μ);
107 DenseVector rhs = A_hat * rhs_temp;
108
109 Eigen::SimplicialLDLT<SparseMatrix> yz_estimator{lhs};
110 DenseVector sol = yz_estimator.solve(rhs);
111 DenseVector y = sol.segment(0, A_e.rows());
112 DenseVector z = sol.segment(A_e.rows(), s.rows());
113
114 // A requirement for the convergence proof is that the primal-dual barrier
115 // term Hessian Σₖ₊₁ does not deviate arbitrarily much from the primal barrier
116 // term Hessian μSₖ₊₁⁻².
117 //
118 // Σₖ₊₁ = μSₖ₊₁⁻²
119 // Sₖ₊₁⁻¹Zₖ₊₁ = μSₖ₊₁⁻²
120 // Zₖ₊₁ = μSₖ₊₁⁻¹
121 //
122 // We ensure this by resetting
123 //
124 // zₖ₊₁ = clamp(zₖ₊₁, 1/κ_Σ μ/sₖ₊₁, κ_Σ μ/sₖ₊₁)
125 //
126 // for some fixed κ_Σ ≥ 1 after each step. See equation (16) of [2].
127 for (int row = 0; row < z.rows(); ++row) {
128 constexpr Scalar κ_Σ(1e10);
129 z[row] = std::clamp(z[row], Scalar(1) / κ_Σ * μ / s[row], κ_Σ * μ / s[row]);
130 }
131
132 return {std::move(y), std::move(z)};
133}
134
135} // namespace slp
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