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
21template <typename Scalar>
23 public:
33 : m_num_decision_variables{num_decision_variables},
34 m_num_equality_constraints{num_equality_constraints} {}
35
41 Eigen::ComputationInfo info() const { return m_info; }
42
49 RegularizedLDLT& compute(const Eigen::SparseMatrix<Scalar>& lhs) {
50 // The regularization procedure is based on algorithm B.1 of [1]
51
52 // Max density is 50% due to the caller only providing the lower triangle.
53 // We consider less than 25% to be sparse.
54 m_is_sparse = lhs.nonZeros() < 0.25 * lhs.size();
55
56 m_info = m_is_sparse ? compute_sparse(lhs).info()
57 : m_dense_solver.compute(lhs).info();
58
60
61 if (m_info == Eigen::Success) {
62 inertia = m_is_sparse ? Inertia{m_sparse_solver.vectorD()}
63 : Inertia{m_dense_solver.vectorD()};
64
65 // If the inertia is ideal, don't regularize the system
66 if (inertia == ideal_inertia) {
67 m_prev_δ = Scalar(0);
68 return *this;
69 }
70 }
71
72 // Also regularize the Hessian. If the Hessian wasn't regularized in a
73 // previous run of compute(), start at small values of δ and γ. Otherwise,
74 // attempt a δ and γ half as big as the previous run so δ and γ can trend
75 // downwards over time.
76 Scalar δ = m_prev_δ == Scalar(0) ? Scalar(1e-4) : m_prev_δ / Scalar(2);
77 Scalar γ(1e-10);
78
79 while (true) {
80 // Regularize lhs by adding a multiple of the identity matrix
81 //
82 // lhs = [H + AᵢᵀΣAᵢ + δI Aₑᵀ]
83 // [ Aₑ −γI]
84 if (m_is_sparse) {
85 m_info = compute_sparse(lhs + regularization(δ, γ)).info();
86 if (m_info == Eigen::Success) {
87 inertia = Inertia{m_sparse_solver.vectorD()};
88 }
89 } else {
90 m_info = m_dense_solver.compute(lhs + regularization(δ, γ)).info();
91 if (m_info == Eigen::Success) {
92 inertia = Inertia{m_dense_solver.vectorD()};
93 }
94 }
95
96 if (m_info == Eigen::Success) {
97 if (inertia == ideal_inertia) {
98 // If the inertia is ideal, store δ and return
99 m_prev_δ = δ;
100 return *this;
101 } else if (inertia.zero > 0) {
102 // If there's zero eigenvalues, increase δ and γ by an order of
103 // magnitude and try again
104 δ *= Scalar(10);
105 γ *= Scalar(10);
106 } else if (inertia.negative > ideal_inertia.negative) {
107 // If there's too many negative eigenvalues, increase δ by an order of
108 // magnitude and try again
109 δ *= Scalar(10);
110 } else if (inertia.positive > ideal_inertia.positive) {
111 // If there's too many positive eigenvalues, increase γ by an order of
112 // magnitude and try again
113 γ *= Scalar(10);
114 }
115 } else {
116 // If the decomposition failed, increase δ and γ by an order of
117 // magnitude and try again
118 δ *= Scalar(10);
119 γ *= Scalar(10);
120 }
121
122 // If the Hessian perturbation is too high, report failure. This can be
123 // caused by ill-conditioning.
124 if (δ > Scalar(1e20) || γ > Scalar(1e20)) {
125 m_info = Eigen::NumericalIssue;
126 return *this;
127 }
128 }
129 }
130
137 template <typename Rhs>
138 Eigen::Vector<Scalar, Eigen::Dynamic> solve(
139 const Eigen::MatrixBase<Rhs>& rhs) {
140 if (m_is_sparse) {
141 return m_sparse_solver.solve(rhs);
142 } else {
143 return m_dense_solver.solve(rhs);
144 }
145 }
146
153 template <typename Rhs>
154 Eigen::Vector<Scalar, Eigen::Dynamic> solve(
155 const Eigen::SparseMatrixBase<Rhs>& rhs) {
156 if (m_is_sparse) {
157 return m_sparse_solver.solve(rhs);
158 } else {
159 return m_dense_solver.solve(rhs.toDense());
160 }
161 }
162
168 Scalar hessian_regularization() const { return m_prev_δ; }
169
170 private:
171 using SparseSolver = Eigen::SimplicialLDLT<Eigen::SparseMatrix<Scalar>>;
172 using DenseSolver =
173 Eigen::LDLT<Eigen::Matrix<Scalar, Eigen::Dynamic, Eigen::Dynamic>>;
174
175 SparseSolver m_sparse_solver;
176 DenseSolver m_dense_solver;
177 bool m_is_sparse = true;
178
179 Eigen::ComputationInfo m_info = Eigen::Success;
180
182 int m_num_decision_variables = 0;
183
185 int m_num_equality_constraints = 0;
186
188 Inertia ideal_inertia{m_num_decision_variables, m_num_equality_constraints,
189 0};
190
192 Scalar m_prev_δ{0};
193
194 // Number of non-zeros in LHS.
195 int m_non_zeros = -1;
196
203 SparseSolver& compute_sparse(const Eigen::SparseMatrix<Scalar>& lhs) {
204 // Reanalize lhs's sparsity pattern if it changed
205 int non_zeros = lhs.nonZeros();
206 if (m_non_zeros != non_zeros) {
207 m_sparse_solver.analyzePattern(lhs);
208 m_non_zeros = non_zeros;
209 }
210
211 m_sparse_solver.factorize(lhs);
212
213 return m_sparse_solver;
214 }
215
223 Eigen::SparseMatrix<Scalar> regularization(Scalar δ, Scalar γ) {
224 Eigen::Vector<Scalar, Eigen::Dynamic> vec{m_num_decision_variables +
225 m_num_equality_constraints};
226 vec.segment(0, m_num_decision_variables).setConstant(δ);
227 vec.segment(m_num_decision_variables, m_num_equality_constraints)
228 .setConstant(-γ);
229
230 return Eigen::SparseMatrix<Scalar>{vec.asDiagonal()};
231 }
232};
233
234} // namespace slp
Definition inertia.hpp:13
int positive
The number of positive eigenvalues.
Definition inertia.hpp:16
int negative
The number of negative eigenvalues.
Definition inertia.hpp:18
Definition intrusive_shared_ptr.hpp:29
Definition regularized_ldlt.hpp:22
Eigen::Vector< Scalar, Eigen::Dynamic > solve(const Eigen::SparseMatrixBase< Rhs > &rhs)
Definition regularized_ldlt.hpp:154
Scalar hessian_regularization() const
Definition regularized_ldlt.hpp:168
Eigen::ComputationInfo info() const
Definition regularized_ldlt.hpp:41
RegularizedLDLT & compute(const Eigen::SparseMatrix< Scalar > &lhs)
Definition regularized_ldlt.hpp:49
Eigen::Vector< Scalar, Eigen::Dynamic > solve(const Eigen::MatrixBase< Rhs > &rhs)
Definition regularized_ldlt.hpp:138
RegularizedLDLT(int num_decision_variables, int num_equality_constraints)
Definition regularized_ldlt.hpp:32