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
32template <typename Scalar>
33Eigen::Vector<Scalar, Eigen::Dynamic> lagrange_multiplier_estimate(
34 const Eigen::SparseVector<Scalar>& g,
35 const Eigen::SparseMatrix<Scalar>& A_e) {
36 // Lagrange multiplier estimates
37 //
38 // ∇f − Aₑᵀy = 0
39 // Aₑᵀy = ∇f
40 // y = (AₑAₑᵀ)⁻¹Aₑ∇f
41 return Eigen::SimplicialLDLT<Eigen::SparseMatrix<Scalar>>{A_e *
42 A_e.transpose()}
43 .solve(A_e * g);
44}
45
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>;
62
63 // Lagrange multiplier estimates
64 //
65 // ∇f − Aₑᵀy − Aᵢᵀz = 0
66 // Sz − μe = 0
67 //
68 // Aₑᵀy + Aᵢᵀz = ∇f
69 // −Sz = −μe
70 //
71 // [Aₑᵀ Aᵢᵀ][y] = [ ∇f]
72 // [ 0 −S ][z] [−μe]
73 //
74 // [Aₑ 0]ᵀ[y] = [ ∇f]
75 // [Aᵢ −S] [z] [−μe]
76 //
77 // Let  = [Aₑ 0]
78 // [Aᵢ −S]
79 //
80 // Âᵀ[y] = [ ∇f]
81 // [z] [−μe]
82 //
83 // [y] = (ÂÂᵀ)⁻¹Â[ ∇f]
84 // [z] [−μe]
85
86 gch::small_vector<Eigen::Triplet<Scalar>> triplets;
87
88 // Â = [Aₑ 0]
89 // [Aᵢ −S]
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());
95
96 // lhs = ÂÂᵀ
97 SparseMatrix lhs = A_hat * A_hat.transpose();
98
99 // rhs = Â[ ∇f]
100 // [−μe]
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;
105
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());
110
111 // A requirement for the convergence proof is that the primal-dual barrier
112 // term Hessian Σₖ₊₁ does not deviate arbitrarily much from the primal barrier
113 // term Hessian μSₖ₊₁⁻².
114 //
115 // Σₖ₊₁ = μSₖ₊₁⁻²
116 // Sₖ₊₁⁻¹Zₖ₊₁ = μSₖ₊₁⁻²
117 // Zₖ₊₁ = μSₖ₊₁⁻¹
118 //
119 // We ensure this by resetting
120 //
121 // zₖ₊₁ = clamp(zₖ₊₁, 1/κ_Σ μ/sₖ₊₁, κ_Σ μ/sₖ₊₁)
122 //
123 // for some fixed κ_Σ ≥ 1 after each step. See equation (16) of [2].
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]);
127 }
128
129 return {std::move(y), std::move(z)};
130}
131
132} // 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