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
36 : m_num_decision_variables{num_decision_variables},
37 m_num_equality_constraints{num_equality_constraints} {}
38
42 Eigen::ComputationInfo info() const { return m_info; }
43
49 // The regularization procedure is based on algorithm B.1 of [1]
50
51 // Max density is 50% due to the caller only providing the lower triangle.
52 // We consider less than 25% to be sparse.
53 m_is_sparse = lhs.nonZeros() < 0.25 * lhs.size();
54
55 m_info = m_is_sparse ? compute_sparse(lhs).info()
56 : m_dense_solver.compute(lhs).info();
57
59
60 if (m_info == Eigen::Success) {
61 inertia = m_is_sparse ? Inertia{m_sparse_solver.vectorD()}
62 : Inertia{m_dense_solver.vectorD()};
63
64 // If the inertia is ideal, don't regularize the system
65 if (inertia == ideal_inertia) {
66 m_prev_δ = Scalar(0);
67 return *this;
68 }
69 }
70
71 // Also regularize the Hessian. If the Hessian wasn't regularized in a
72 // previous run of compute(), start at small values of δ and γ. Otherwise,
73 // attempt a δ and γ half as big as the previous run so δ and γ can trend
74 // downwards over time.
75 Scalar δ = m_prev_δ == Scalar(0) ? Scalar(1e-4) : m_prev_δ / Scalar(2);
76 Scalar γ(1e-10);
77
78 while (true) {
79 // Regularize lhs by adding a multiple of the identity matrix
80 //
81 // lhs = [H + AᵢᵀΣAᵢ + δI Aₑᵀ]
82 // [ Aₑ −γI]
83 if (m_is_sparse) {
84 m_info = compute_sparse(lhs + regularization(δ, γ)).info();
85 if (m_info == Eigen::Success) {
86 inertia = Inertia{m_sparse_solver.vectorD()};
87 }
88 } else {
89 m_info = m_dense_solver.compute(lhs + regularization(δ, γ)).info();
90 if (m_info == Eigen::Success) {
91 inertia = Inertia{m_dense_solver.vectorD()};
92 }
93 }
94
95 if (m_info == Eigen::Success) {
96 if (inertia == ideal_inertia) {
97 // If the inertia is ideal, store δ and return
98 m_prev_δ = δ;
99 return *this;
100 } else if (inertia.zero > 0) {
101 // If there's zero eigenvalues, increase δ and γ by an order of
102 // magnitude and try again
103 δ *= Scalar(10);
104 γ *= Scalar(10);
105 } else if (inertia.negative > ideal_inertia.negative) {
106 // If there's too many negative eigenvalues, increase δ by an order of
107 // magnitude and try again
108 δ *= Scalar(10);
109 } else if (inertia.positive > ideal_inertia.positive) {
110 // If there's too many positive eigenvalues, increase γ by an order of
111 // magnitude and try again
112 γ *= Scalar(10);
113 }
114 } else {
115 // If the decomposition failed, increase δ and γ by an order of
116 // magnitude and try again
117 δ *= Scalar(10);
118 γ *= Scalar(10);
119 }
120
121 // If the Hessian perturbation is too high, report failure. This can be
122 // caused by ill-conditioning.
123 if (δ > Scalar(1e20) || γ > Scalar(1e20)) {
124 m_info = Eigen::NumericalIssue;
125 return *this;
126 }
127 }
128 }
129
134 template <typename Rhs>
135 DenseVector solve(const Eigen::MatrixBase<Rhs>& rhs) {
136 if (m_is_sparse) {
137 return m_sparse_solver.solve(rhs);
138 } else {
139 return m_dense_solver.solve(rhs);
140 }
141 }
142
147 template <typename Rhs>
148 DenseVector solve(const Eigen::SparseMatrixBase<Rhs>& rhs) {
149 if (m_is_sparse) {
150 return m_sparse_solver.solve(rhs);
151 } else {
152 return m_dense_solver.solve(rhs.toDense());
153 }
154 }
155
159 Scalar hessian_regularization() const { return m_prev_δ; }
160
161 private:
162 using SparseSolver = Eigen::SimplicialLDLT<SparseMatrix>;
163 using DenseSolver = Eigen::LDLT<DenseMatrix>;
164
165 SparseSolver m_sparse_solver;
166 DenseSolver m_dense_solver;
167 bool m_is_sparse = true;
168
169 Eigen::ComputationInfo m_info = Eigen::Success;
170
172 int m_num_decision_variables = 0;
173
175 int m_num_equality_constraints = 0;
176
178 Inertia ideal_inertia{m_num_decision_variables, m_num_equality_constraints,
179 0};
180
182 Scalar m_prev_δ{0};
183
184 // Number of non-zeros in LHS.
185 int m_non_zeros = -1;
186
191 SparseSolver& compute_sparse(const SparseMatrix& lhs) {
192 // Reanalize lhs's sparsity pattern if it changed
193 int non_zeros = lhs.nonZeros();
194 if (m_non_zeros != non_zeros) {
195 m_sparse_solver.analyzePattern(lhs);
196 m_non_zeros = non_zeros;
197 }
198
199 m_sparse_solver.factorize(lhs);
200
201 return m_sparse_solver;
202 }
203
209 SparseMatrix regularization(Scalar δ, Scalar γ) {
210 DenseVector vec{m_num_decision_variables + m_num_equality_constraints};
211 vec.segment(0, m_num_decision_variables).setConstant(δ);
212 vec.segment(m_num_decision_variables, m_num_equality_constraints)
213 .setConstant(-γ);
214
215 return SparseMatrix{vec.asDiagonal()};
216 }
217};
218
219} // 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
DenseVector solve(const Eigen::MatrixBase< Rhs > &rhs)
Definition regularized_ldlt.hpp:135
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:48
Scalar hessian_regularization() const
Definition regularized_ldlt.hpp:159
Eigen::ComputationInfo info() const
Definition regularized_ldlt.hpp:42
DenseVector solve(const Eigen::SparseMatrixBase< Rhs > &rhs)
Definition regularized_ldlt.hpp:148
RegularizedLDLT(int num_decision_variables, int num_equality_constraints)
Definition regularized_ldlt.hpp:35