Sleipnir C++ API
Loading...
Searching...
No Matches
sparse_regularized_ldlt.hpp
1// Copyright (c) Sleipnir contributors
2
3#pragma once
4
5#include <Eigen/Core>
6#include <Eigen/SparseCholesky>
7#include <Eigen/SparseCore>
8
9#include "sleipnir/optimization/solver/util/inertia.hpp"
10
11// See docs/algorithms.md#Works_cited for citation definitions
12
13namespace slp {
14
19template <typename Scalar>
21 public:
23 using DenseVector = Eigen::Vector<Scalar, Eigen::Dynamic>;
25 using SparseMatrix = Eigen::SparseMatrix<Scalar>;
26
33 SparseRegularizedLDLT(int num_decision_variables,
34 int num_equality_constraints)
35 : m_num_decision_variables{num_decision_variables},
36 m_num_equality_constraints{num_equality_constraints} {}
37
45 SparseRegularizedLDLT(int num_decision_variables,
46 int num_equality_constraints, Scalar γ_min)
47 : m_num_decision_variables{num_decision_variables},
48 m_num_equality_constraints{num_equality_constraints},
49 m_γ_min{γ_min} {}
50
54 Eigen::ComputationInfo info() const { return m_info; }
55
64 // The regularization procedure is based on algorithm B.1 of [1]
65
66 // Regularization with zeros ensures the pattern analysis in the sparse
67 // solver is reused by all factorizations
68 SparseMatrix unregularized_lhs = lhs + regularization(Scalar(0), Scalar(0));
69
70 if (!m_analyzed_pattern) {
71 m_solver.analyzePattern(unregularized_lhs);
72 m_analyzed_pattern = true;
73 }
74
75 m_solver.factorize(unregularized_lhs);
76 m_info = m_solver.info();
77
78 if (m_info == Eigen::Success) {
79 auto D = m_solver.vectorD();
80
81 // If the inertia is ideal and D from LDLT is sufficiently far from zero,
82 // don't regularize the system
83 if (Inertia{D} == ideal_inertia &&
84 (D.cwiseAbs().array() >= Scalar(1e-4)).all()) {
85 m_prev_δ = Scalar(0);
86 m_prev_γ = Scalar(0);
87 return *this;
88 }
89 }
90
91 // Also regularize the Hessian. If the Hessian wasn't regularized in a
92 // previous run of compute(), start at small values of δ and γ. Otherwise,
93 // attempt a δ and γ half as big as the previous run so δ and γ can trend
94 // downwards over time.
95 Scalar δ = m_prev_δ == Scalar(0) ? Scalar(1e-4) : m_prev_δ / Scalar(2);
96 Scalar γ = m_γ_min;
97
98 while (true) {
99 m_solver.factorize(lhs + regularization(δ, γ));
100 m_info = m_solver.info();
101
102 if (m_info == Eigen::Success) {
103 Inertia inertia{m_solver.vectorD()};
104
105 if (inertia == ideal_inertia) {
106 // If the inertia is ideal, report success
107 m_prev_δ = δ;
108 m_prev_γ = γ;
109 return *this;
110 } else if (inertia.zero > 0) {
111 if (γ == Scalar(0)) {
112 // If there's zero eigenvalues and γ = 0, increase γ
113 γ = Scalar(1e-10);
114 } else {
115 // If there's zero eigenvalues and γ > 0, increase δ and γ
116 δ *= Scalar(10);
117 γ *= Scalar(10);
118 }
119 } else if (inertia.negative > ideal_inertia.negative) {
120 // If there's too many negative eigenvalues, increase δ
121 δ *= Scalar(10);
122 } else if (inertia.positive > ideal_inertia.positive) {
123 // If there's too many positive eigenvalues, increase γ
124 γ = γ == Scalar(0) ? Scalar(1e-10) : γ * Scalar(10);
125 }
126 } else {
127 // If the decomposition failed, increase δ and γ
128 δ *= Scalar(10);
129 γ = γ == Scalar(0) ? Scalar(1e-10) : γ * Scalar(10);
130 }
131
132 // If the Hessian perturbation is too high, report failure. This can be
133 // caused by ill-conditioning.
134 if (δ > Scalar(1e20) || γ > Scalar(1e20)) {
135 m_info = Eigen::NumericalIssue;
136 m_prev_δ = δ;
137 m_prev_γ = γ;
138 return *this;
139 }
140 }
141 }
142
147 template <typename Rhs>
148 DenseVector solve(const Eigen::MatrixBase<Rhs>& rhs) const {
149 return m_solver.solve(rhs);
150 }
151
156 template <typename Rhs>
157 DenseVector solve(const Eigen::SparseMatrixBase<Rhs>& rhs) const {
158 return m_solver.solve(rhs);
159 }
160
164 Scalar hessian_regularization() const { return m_prev_δ; }
165
169 Scalar constraint_jacobian_regularization() const { return m_prev_γ; }
170
171 private:
172 using Solver = Eigen::SimplicialLDLT<SparseMatrix>;
173
174 Solver m_solver;
175 bool m_analyzed_pattern = false;
176
177 Eigen::ComputationInfo m_info = Eigen::Success;
178
180 int m_num_decision_variables = 0;
181
183 int m_num_equality_constraints = 0;
184
186 Scalar m_γ_min{1e-10};
187
189 Inertia ideal_inertia{m_num_decision_variables, m_num_equality_constraints,
190 0};
191
193 Scalar m_prev_δ{0};
194
196 Scalar m_prev_γ{0};
197
206 SparseMatrix regularization(Scalar δ, Scalar γ) const {
207 DenseVector vec{m_num_decision_variables + m_num_equality_constraints};
208 vec.segment(0, m_num_decision_variables).setConstant(δ);
209 vec.segment(m_num_decision_variables, m_num_equality_constraints)
210 .setConstant(-γ);
211
212 return SparseMatrix{vec.asDiagonal()};
213 }
214};
215
216} // 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:20
Eigen::Vector< Scalar, Eigen::Dynamic > DenseVector
Type alias for dense vector.
Definition sparse_regularized_ldlt.hpp:23
Eigen::ComputationInfo info() const
Definition sparse_regularized_ldlt.hpp:54
Scalar constraint_jacobian_regularization() const
Definition sparse_regularized_ldlt.hpp:169
SparseRegularizedLDLT & compute(const SparseMatrix &lhs)
Definition sparse_regularized_ldlt.hpp:63
DenseVector solve(const Eigen::MatrixBase< Rhs > &rhs) const
Definition sparse_regularized_ldlt.hpp:148
DenseVector solve(const Eigen::SparseMatrixBase< Rhs > &rhs) const
Definition sparse_regularized_ldlt.hpp:157
SparseRegularizedLDLT(int num_decision_variables, int num_equality_constraints, Scalar γ_min)
Definition sparse_regularized_ldlt.hpp:45
SparseRegularizedLDLT(int num_decision_variables, int num_equality_constraints)
Definition sparse_regularized_ldlt.hpp:33
Scalar hessian_regularization() const
Definition sparse_regularized_ldlt.hpp:164
Eigen::SparseMatrix< Scalar > SparseMatrix
Type alias for sparse matrix.
Definition sparse_regularized_ldlt.hpp:25