Sleipnir C++ API
Loading...
Searching...
No Matches
sparse_regularized_ldlt.hpp
1// Copyright (c) Sleipnir contributors
2
3#pragma once
4
5#include <algorithm>
6#include <limits>
7
8#include <Eigen/Core>
9#include <Eigen/SparseCholesky>
10#include <Eigen/SparseCore>
11
12#include "sleipnir/optimization/solver/util/inertia.hpp"
13
14namespace slp {
15
20template <typename Scalar>
22 public:
24 using DenseVector = Eigen::Vector<Scalar, Eigen::Dynamic>;
26 using SparseMatrix = Eigen::SparseMatrix<Scalar>;
27
34 SparseRegularizedLDLT(int num_decision_variables,
35 int num_equality_constraints)
36 : m_num_decision_variables{num_decision_variables},
37 m_num_equality_constraints{num_equality_constraints} {}
38
46 SparseRegularizedLDLT(int num_decision_variables,
47 int num_equality_constraints, Scalar γ_min)
48 : m_num_decision_variables{num_decision_variables},
49 m_num_equality_constraints{num_equality_constraints},
50 m_γ_min{γ_min} {}
51
55 Eigen::ComputationInfo info() const { return m_info; }
56
65 // Regularization with zeros ensures the pattern analysis in the sparse
66 // solver is reused by all factorizations
67 SparseMatrix unregularized_lhs = lhs + regularization(Scalar(0), Scalar(0));
68
69 if (!m_analyzed_pattern) {
70 m_solver.analyzePattern(unregularized_lhs);
71 m_analyzed_pattern = true;
72 }
73
74 m_solver.factorize(unregularized_lhs);
75 m_info = m_solver.info();
76
77 if (m_info == Eigen::Success) {
78 auto D = m_solver.vectorD();
79
80 // If the inertia is ideal and D from LDLT is sufficiently far from zero,
81 // don't regularize the system
82 if (Inertia{D} == ideal_inertia &&
83 (D.cwiseAbs().array() >= Scalar(1e-4)).all()) {
84 m_prev_δ = Scalar(0);
85 m_prev_γ = Scalar(0);
86 return *this;
87 }
88 }
89
90 // We'll give lhs the correct inertia by adding [δI, 0; 0, −γI] where δ and
91 // γ regularize the Hessian and equality constraint Jacobian respectively.
92
93 // If the previous δ was zero, attempt a small value. Otherwise, attempt a
94 // smaller value than the previous δ so δ trends downward.
95 Scalar δ = m_prev_δ == Scalar(0)
96 ? Scalar(1e-4)
97 : std::max(m_prev_δ / Scalar(2),
98 std::numeric_limits<Scalar>::epsilon());
99
100 // Start γ at the minimum to minimize equality constraint Jacobian
101 // distortion.
102 Scalar γ = m_γ_min;
103
104 while (true) {
105 m_solver.factorize(lhs + regularization(δ, γ));
106 m_info = m_solver.info();
107
108 if (m_info == Eigen::Success) {
109 Inertia inertia{m_solver.vectorD()};
110
111 if (inertia == ideal_inertia) {
112 // If the inertia is ideal, report success
113 m_prev_δ = δ;
114 m_prev_γ = γ;
115 return *this;
116 } else if (inertia.zero > 0) {
117 if (γ == Scalar(0)) {
118 // If there's zero eigenvalues and γ = 0, increase γ to potentially
119 // compensate for a rank-deficient equality constraint Jacobian
120 γ = Scalar(1e-10);
121 } else {
122 // If there's zero eigenvalues and γ > 0, increase δ and γ to drive
123 // all eigenvalues away from zero
124 δ *= Scalar(10);
125 γ *= Scalar(10);
126 }
127 } else if (inertia.negative > ideal_inertia.negative) {
128 // If there's too many negative eigenvalues, increase δ to add more
129 // positive eigenvalues
130 δ *= Scalar(10);
131 } else if (inertia.positive > ideal_inertia.positive) {
132 // If there's too many positive eigenvalues, increase γ to add more
133 // negative eigenvalues
134 γ = γ == Scalar(0) ? Scalar(1e-10) : γ * Scalar(10);
135 }
136 } else {
137 // If the decomposition failed, increase δ and γ to drive all
138 // eigenvalues away from zero
139 δ *= Scalar(10);
140 γ = γ == Scalar(0) ? Scalar(1e-10) : γ * Scalar(10);
141 }
142
143 // If the lhs perturbation is too high, report failure. This can be caused
144 // by ill-conditioning.
145 if (δ > Scalar(1e20) || γ > Scalar(1e20)) {
146 m_info = Eigen::NumericalIssue;
147 m_prev_δ = δ;
148 m_prev_γ = γ;
149 return *this;
150 }
151 }
152 }
153
158 template <typename Rhs>
159 DenseVector solve(const Eigen::MatrixBase<Rhs>& rhs) const {
160 return m_solver.solve(rhs);
161 }
162
167 template <typename Rhs>
168 DenseVector solve(const Eigen::SparseMatrixBase<Rhs>& rhs) const {
169 return m_solver.solve(rhs);
170 }
171
175 Scalar hessian_regularization() const { return m_prev_δ; }
176
180 Scalar constraint_jacobian_regularization() const { return m_prev_γ; }
181
182 private:
183 using Solver = Eigen::SimplicialLDLT<SparseMatrix>;
184
185 Solver m_solver;
186 bool m_analyzed_pattern = false;
187
188 Eigen::ComputationInfo m_info = Eigen::Success;
189
191 int m_num_decision_variables = 0;
192
194 int m_num_equality_constraints = 0;
195
197 Scalar m_γ_min{1e-10};
198
200 Inertia ideal_inertia{m_num_decision_variables, m_num_equality_constraints,
201 0};
202
204 Scalar m_prev_δ{0};
205
207 Scalar m_prev_γ{0};
208
217 SparseMatrix regularization(Scalar δ, Scalar γ) const {
218 DenseVector vec{m_num_decision_variables + m_num_equality_constraints};
219 vec.segment(0, m_num_decision_variables).setConstant(δ);
220 vec.segment(m_num_decision_variables, m_num_equality_constraints)
221 .setConstant(-γ);
222
223 return SparseMatrix{vec.asDiagonal()};
224 }
225};
226
227} // namespace slp
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
Definition sparse_regularized_ldlt.hpp:21
Eigen::Vector< Scalar, Eigen::Dynamic > DenseVector
Type alias for dense vector.
Definition sparse_regularized_ldlt.hpp:24
Eigen::ComputationInfo info() const
Definition sparse_regularized_ldlt.hpp:55
Scalar constraint_jacobian_regularization() const
Definition sparse_regularized_ldlt.hpp:180
SparseRegularizedLDLT & compute(const SparseMatrix &lhs)
Definition sparse_regularized_ldlt.hpp:64
DenseVector solve(const Eigen::MatrixBase< Rhs > &rhs) const
Definition sparse_regularized_ldlt.hpp:159
DenseVector solve(const Eigen::SparseMatrixBase< Rhs > &rhs) const
Definition sparse_regularized_ldlt.hpp:168
SparseRegularizedLDLT(int num_decision_variables, int num_equality_constraints, Scalar γ_min)
Definition sparse_regularized_ldlt.hpp:46
SparseRegularizedLDLT(int num_decision_variables, int num_equality_constraints)
Definition sparse_regularized_ldlt.hpp:34
Scalar hessian_regularization() const
Definition sparse_regularized_ldlt.hpp:175
Eigen::SparseMatrix< Scalar > SparseMatrix
Type alias for sparse matrix.
Definition sparse_regularized_ldlt.hpp:26