Sleipnir C++ API
Loading...
Searching...
No Matches
regularized_ldlt.hpp
1// Copyright (c) Sleipnir contributors
2
3#pragma once
4
5#include <Eigen/Cholesky>
6#include <Eigen/Core>
7#include <Eigen/SparseCholesky>
8#include <Eigen/SparseCore>
9
10#include "sleipnir/optimization/solver/util/inertia.hpp"
11
12// See docs/algorithms.md#Works_cited for citation definitions
13
14namespace slp {
15
19template <typename Scalar>
21 public:
23 using DenseMatrix = Eigen::Matrix<Scalar, Eigen::Dynamic, Eigen::Dynamic>;
25 using DenseVector = Eigen::Vector<Scalar, Eigen::Dynamic>;
27 using SparseMatrix = Eigen::SparseMatrix<Scalar>;
28
35 RegularizedLDLT(int num_decision_variables, int num_equality_constraints)
36 : m_num_decision_variables{num_decision_variables},
37 m_num_equality_constraints{num_equality_constraints} {}
38
46 RegularizedLDLT(int num_decision_variables, int num_equality_constraints,
47 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
62 // The regularization procedure is based on algorithm B.1 of [1]
63
64 // Max density is 50% due to the caller only providing the lower triangle.
65 // We consider less than 25% to be sparse.
66 m_is_sparse = lhs.nonZeros() < 0.25 * lhs.size();
67
68 m_info = m_is_sparse ? compute_sparse(lhs).info()
69 : m_dense_solver.compute(lhs).info();
70
71 if (m_info == Eigen::Success) {
72 auto D =
73 m_is_sparse ? m_sparse_solver.vectorD() : m_dense_solver.vectorD();
74
75 // If the inertia is ideal and D from LDLT is sufficiently far from zero,
76 // don't regularize the system
77 if (Inertia{D} == ideal_inertia &&
78 (D.cwiseAbs().array() >= Scalar(1e-4)).all()) {
79 m_prev_δ = Scalar(0);
80 return *this;
81 }
82 }
83
84 // Also regularize the Hessian. If the Hessian wasn't regularized in a
85 // previous run of compute(), start at small values of δ and γ. Otherwise,
86 // attempt a δ and γ half as big as the previous run so δ and γ can trend
87 // downwards over time.
88 Scalar δ = m_prev_δ == Scalar(0) ? Scalar(1e-4) : m_prev_δ / Scalar(2);
89 Scalar γ = m_γ_min;
90
91 while (true) {
93
94 // Regularize lhs by adding a multiple of the identity matrix
95 //
96 // lhs = [H + AᵢᵀΣAᵢ + δI Aₑᵀ]
97 // [ Aₑ −γI]
98 if (m_is_sparse) {
99 m_info = compute_sparse(lhs + regularization(δ, γ)).info();
100 if (m_info == Eigen::Success) {
101 inertia = Inertia{m_sparse_solver.vectorD()};
102 }
103 } else {
104 m_info = m_dense_solver.compute(lhs + regularization(δ, γ)).info();
105 if (m_info == Eigen::Success) {
106 inertia = Inertia{m_dense_solver.vectorD()};
107 }
108 }
109
110 if (m_info == Eigen::Success) {
111 if (inertia == ideal_inertia) {
112 // If the inertia is ideal, store δ and return
113 m_prev_δ = δ;
114 return *this;
115 } else if (inertia.zero > 0) {
116 // If there's zero eigenvalues, increase δ and γ by an order of
117 // magnitude and try again
118 δ *= Scalar(10);
119 γ *= Scalar(10);
120 } else if (inertia.negative > ideal_inertia.negative) {
121 // If there's too many negative eigenvalues, increase δ by an order of
122 // magnitude and try again
123 δ *= Scalar(10);
124 } else if (inertia.positive > ideal_inertia.positive) {
125 // If there's too many positive eigenvalues, increase γ by an order of
126 // magnitude and try again
127 γ = γ == Scalar(0) ? Scalar(1e-10) : γ * Scalar(10);
128 }
129 } else {
130 // If the decomposition failed, increase δ and γ by an order of
131 // magnitude and try again
132 δ *= Scalar(10);
133 γ *= Scalar(10);
134 }
135
136 // If the Hessian perturbation is too high, report failure. This can be
137 // caused by ill-conditioning.
138 if (δ > Scalar(1e20) || γ > Scalar(1e20)) {
139 m_info = Eigen::NumericalIssue;
140 return *this;
141 }
142 }
143 }
144
149 template <typename Rhs>
150 DenseVector solve(const Eigen::MatrixBase<Rhs>& rhs) const {
151 if (m_is_sparse) {
152 return m_sparse_solver.solve(rhs);
153 } else {
154 return m_dense_solver.solve(rhs);
155 }
156 }
157
162 template <typename Rhs>
163 DenseVector solve(const Eigen::SparseMatrixBase<Rhs>& rhs) const {
164 if (m_is_sparse) {
165 return m_sparse_solver.solve(rhs);
166 } else {
167 return m_dense_solver.solve(rhs.toDense());
168 }
169 }
170
174 Scalar hessian_regularization() const { return m_prev_δ; }
175
176 private:
177 using SparseSolver = Eigen::SimplicialLDLT<SparseMatrix>;
178 using DenseSolver = Eigen::LDLT<DenseMatrix>;
179
180 SparseSolver m_sparse_solver;
181 DenseSolver m_dense_solver;
182 bool m_is_sparse = true;
183
184 Eigen::ComputationInfo m_info = Eigen::Success;
185
187 int m_num_decision_variables = 0;
188
190 int m_num_equality_constraints = 0;
191
193 Scalar m_γ_min{1e-10};
194
196 Inertia ideal_inertia{m_num_decision_variables, m_num_equality_constraints,
197 0};
198
200 Scalar m_prev_δ{0};
201
202 // Number of non-zeros in LHS.
203 int m_non_zeros = -1;
204
209 SparseSolver& compute_sparse(const SparseMatrix& lhs) {
210 // Reanalize lhs's sparsity pattern if it changed
211 int non_zeros = lhs.nonZeros();
212 if (m_non_zeros != non_zeros) {
213 m_sparse_solver.analyzePattern(lhs);
214 m_non_zeros = non_zeros;
215 }
216
217 m_sparse_solver.factorize(lhs);
218
219 return m_sparse_solver;
220 }
221
227 SparseMatrix regularization(Scalar δ, Scalar γ) const {
228 DenseVector vec{m_num_decision_variables + m_num_equality_constraints};
229 vec.segment(0, m_num_decision_variables).setConstant(δ);
230 vec.segment(m_num_decision_variables, m_num_equality_constraints)
231 .setConstant(-γ);
232
233 return SparseMatrix{vec.asDiagonal()};
234 }
235};
236
237} // namespace slp
Definition inertia.hpp:11
int positive
The number of positive eigenvalues.
Definition inertia.hpp:14
int negative
The number of negative eigenvalues.
Definition inertia.hpp:16
Definition intrusive_shared_ptr.hpp:27
Definition regularized_ldlt.hpp:20
Eigen::Matrix< Scalar, Eigen::Dynamic, Eigen::Dynamic > DenseMatrix
Type alias for dense matrix.
Definition regularized_ldlt.hpp:23
Eigen::SparseMatrix< Scalar > SparseMatrix
Type alias for sparse matrix.
Definition regularized_ldlt.hpp:27
Eigen::Vector< Scalar, Eigen::Dynamic > DenseVector
Type alias for dense vector.
Definition regularized_ldlt.hpp:25
RegularizedLDLT & compute(const SparseMatrix &lhs)
Definition regularized_ldlt.hpp:61
Scalar hessian_regularization() const
Definition regularized_ldlt.hpp:174
RegularizedLDLT(int num_decision_variables, int num_equality_constraints, Scalar γ_min)
Definition regularized_ldlt.hpp:46
Eigen::ComputationInfo info() const
Definition regularized_ldlt.hpp:55
DenseVector solve(const Eigen::MatrixBase< Rhs > &rhs) const
Definition regularized_ldlt.hpp:150
DenseVector solve(const Eigen::SparseMatrixBase< Rhs > &rhs) const
Definition regularized_ldlt.hpp:163
RegularizedLDLT(int num_decision_variables, int num_equality_constraints)
Definition regularized_ldlt.hpp:35