Sleipnir C++ API
Loading...
Searching...
No Matches
dense_regularized_ldlt.hpp
1// Copyright (c) Sleipnir contributors
2
3#pragma once
4
5#include <algorithm>
6#include <limits>
7
8#include <Eigen/Cholesky>
9#include <Eigen/Core>
10
11#include "sleipnir/optimization/solver/util/inertia.hpp"
12
13namespace slp {
14
19template <typename Scalar>
21 public:
23 using DenseMatrix = Eigen::Matrix<Scalar, Eigen::Dynamic, Eigen::Dynamic>;
25 using DenseVector = Eigen::Vector<Scalar, Eigen::Dynamic>;
26
33 DenseRegularizedLDLT(int num_decision_variables, int num_equality_constraints)
34 : m_num_decision_variables{num_decision_variables},
35 m_num_equality_constraints{num_equality_constraints} {}
36
44 DenseRegularizedLDLT(int num_decision_variables, int num_equality_constraints,
45 Scalar γ_min)
46 : m_num_decision_variables{num_decision_variables},
47 m_num_equality_constraints{num_equality_constraints},
48 m_γ_min{γ_min} {}
49
53 Eigen::ComputationInfo info() const { return m_info; }
54
60 m_info = m_solver.compute(lhs).info();
61
62 if (m_info == Eigen::Success) {
63 auto D = m_solver.vectorD();
64
65 // If the inertia is ideal and D from LDLT is sufficiently far from zero,
66 // don't regularize the system
67 if (Inertia{D} == ideal_inertia &&
68 (D.cwiseAbs().array() >= Scalar(1e-4)).all()) {
69 m_prev_δ = Scalar(0);
70 m_prev_γ = Scalar(0);
71 return *this;
72 }
73 }
74
75 // We'll give lhs the correct inertia by adding [δI, 0; 0, −γI] where δ and
76 // γ regularize the Hessian and equality constraint Jacobian respectively.
77
78 // If the previous δ was zero, attempt a small value. Otherwise, attempt a
79 // smaller value than the previous δ so δ trends downward.
80 Scalar δ = m_prev_δ == Scalar(0)
81 ? Scalar(1e-4)
82 : std::max(m_prev_δ / Scalar(2),
83 std::numeric_limits<Scalar>::epsilon());
84
85 // Start γ at the minimum to minimize equality constraint Jacobian
86 // distortion.
87 Scalar γ = m_γ_min;
88
89 while (true) {
90 m_info = m_solver.compute(lhs + regularization(δ, γ)).info();
91
92 if (m_info == Eigen::Success) {
93 Inertia inertia{m_solver.vectorD()};
94
95 if (inertia == ideal_inertia) {
96 // If the inertia is ideal, report success
97 m_prev_δ = δ;
98 m_prev_γ = γ;
99 return *this;
100 } else if (inertia.zero > 0) {
101 if (γ == Scalar(0)) {
102 // If there's zero eigenvalues and γ = 0, increase γ to potentially
103 // compensate for a rank-deficient equality constraint Jacobian
104 γ = Scalar(1e-10);
105 } else {
106 // If there's zero eigenvalues and γ > 0, increase δ and γ to drive
107 // all eigenvalues away from zero
108 δ *= Scalar(10);
109 γ *= Scalar(10);
110 }
111 } else if (inertia.negative > ideal_inertia.negative) {
112 // If there's too many negative eigenvalues, increase δ to add more
113 // positive eigenvalues
114 δ *= Scalar(10);
115 } else if (inertia.positive > ideal_inertia.positive) {
116 // If there's too many positive eigenvalues, increase γ to add more
117 // negative eigenvalues
118 γ = γ == Scalar(0) ? Scalar(1e-10) : γ * Scalar(10);
119 }
120 } else {
121 // If the decomposition failed, increase δ and γ to drive all
122 // eigenvalues away from zero
123 δ *= Scalar(10);
124 γ = γ == Scalar(0) ? Scalar(1e-10) : γ * Scalar(10);
125 }
126
127 // If the lhs perturbation is too high, report failure. This can be caused
128 // by ill-conditioning.
129 if (δ > Scalar(1e20) || γ > Scalar(1e20)) {
130 m_info = Eigen::NumericalIssue;
131 m_prev_δ = δ;
132 m_prev_γ = γ;
133 return *this;
134 }
135 }
136 }
137
142 template <typename Rhs>
143 DenseVector solve(const Eigen::MatrixBase<Rhs>& rhs) const {
144 return m_solver.solve(rhs);
145 }
146
151 template <typename Rhs>
152 DenseVector solve(const Eigen::SparseMatrixBase<Rhs>& rhs) const {
153 return m_solver.solve(rhs.toDense());
154 }
155
159 Scalar hessian_regularization() const { return m_prev_δ; }
160
164 Scalar constraint_jacobian_regularization() const { return m_prev_γ; }
165
166 private:
167 using Solver = Eigen::LDLT<DenseMatrix>;
168
169 Solver m_solver;
170
171 Eigen::ComputationInfo m_info = Eigen::Success;
172
174 int m_num_decision_variables = 0;
175
177 int m_num_equality_constraints = 0;
178
180 Scalar m_γ_min{1e-10};
181
183 Inertia ideal_inertia{m_num_decision_variables, m_num_equality_constraints,
184 0};
185
187 Scalar m_prev_δ{0};
188
190 Scalar m_prev_γ{0};
191
200 DenseMatrix regularization(Scalar δ, Scalar γ) const {
201 DenseVector vec{m_num_decision_variables + m_num_equality_constraints};
202 vec.segment(0, m_num_decision_variables).setConstant(δ);
203 vec.segment(m_num_decision_variables, m_num_equality_constraints)
204 .setConstant(-γ);
205
206 return vec.asDiagonal().toDenseMatrix();
207 }
208};
209
210} // namespace slp
Definition dense_regularized_ldlt.hpp:20
DenseVector solve(const Eigen::MatrixBase< Rhs > &rhs) const
Definition dense_regularized_ldlt.hpp:143
DenseRegularizedLDLT(int num_decision_variables, int num_equality_constraints, Scalar γ_min)
Definition dense_regularized_ldlt.hpp:44
Eigen::Vector< Scalar, Eigen::Dynamic > DenseVector
Type alias for dense vector.
Definition dense_regularized_ldlt.hpp:25
Scalar hessian_regularization() const
Definition dense_regularized_ldlt.hpp:159
Eigen::Matrix< Scalar, Eigen::Dynamic, Eigen::Dynamic > DenseMatrix
Type alias for dense matrix.
Definition dense_regularized_ldlt.hpp:23
DenseRegularizedLDLT(int num_decision_variables, int num_equality_constraints)
Definition dense_regularized_ldlt.hpp:33
Eigen::ComputationInfo info() const
Definition dense_regularized_ldlt.hpp:53
DenseRegularizedLDLT & compute(const DenseMatrix &lhs)
Definition dense_regularized_ldlt.hpp:59
Scalar constraint_jacobian_regularization() const
Definition dense_regularized_ldlt.hpp:164
DenseVector solve(const Eigen::SparseMatrixBase< Rhs > &rhs) const
Definition dense_regularized_ldlt.hpp:152
Definition inertia.hpp:14
int positive
The number of positive eigenvalues.
Definition inertia.hpp:17
int negative
The number of negative eigenvalues.
Definition inertia.hpp:19
Definition intrusive_shared_ptr.hpp:27