Sleipnir C++ API
Loading...
Searching...
No Matches
dense_regularized_ldlt.hpp
1// Copyright (c) Sleipnir contributors
2
3#pragma once
4
5#include <Eigen/Cholesky>
6#include <Eigen/Core>
7
8#include "sleipnir/optimization/solver/util/inertia.hpp"
9
10// See docs/algorithms.md#Works_cited for citation definitions
11
12namespace slp {
13
18template <typename Scalar>
20 public:
22 using DenseMatrix = Eigen::Matrix<Scalar, Eigen::Dynamic, Eigen::Dynamic>;
24 using DenseVector = Eigen::Vector<Scalar, Eigen::Dynamic>;
25
32 DenseRegularizedLDLT(int num_decision_variables, int num_equality_constraints)
33 : m_num_decision_variables{num_decision_variables},
34 m_num_equality_constraints{num_equality_constraints} {}
35
43 DenseRegularizedLDLT(int num_decision_variables, int num_equality_constraints,
44 Scalar γ_min)
45 : m_num_decision_variables{num_decision_variables},
46 m_num_equality_constraints{num_equality_constraints},
47 m_γ_min{γ_min} {}
48
52 Eigen::ComputationInfo info() const { return m_info; }
53
59 // The regularization procedure is based on algorithm B.1 of [1]
60
61 m_info = m_solver.compute(lhs).info();
62
63 if (m_info == Eigen::Success) {
64 auto D = m_solver.vectorD();
65
66 // If the inertia is ideal and D from LDLT is sufficiently far from zero,
67 // don't regularize the system
68 if (Inertia{D} == ideal_inertia &&
69 (D.cwiseAbs().array() >= Scalar(1e-4)).all()) {
70 m_prev_δ = Scalar(0);
71 m_prev_γ = Scalar(0);
72 return *this;
73 }
74 }
75
76 // Also regularize the Hessian. If the Hessian wasn't regularized in a
77 // previous run of compute(), start at small values of δ and γ. Otherwise,
78 // attempt a δ and γ half as big as the previous run so δ and γ can trend
79 // downwards over time.
80 Scalar δ = m_prev_δ == Scalar(0) ? Scalar(1e-4) : m_prev_δ / Scalar(2);
81 Scalar γ = m_γ_min;
82
83 while (true) {
84 m_info = m_solver.compute(lhs + regularization(δ, γ)).info();
85
86 if (m_info == Eigen::Success) {
87 Inertia inertia{m_solver.vectorD()};
88
89 if (inertia == ideal_inertia) {
90 // If the inertia is ideal, report success
91 m_prev_δ = δ;
92 m_prev_γ = γ;
93 return *this;
94 } else if (inertia.zero > 0) {
95 if (γ == Scalar(0)) {
96 // If there's zero eigenvalues and γ = 0, increase γ
97 γ = Scalar(1e-10);
98 } else {
99 // If there's zero eigenvalues and γ > 0, increase δ and γ
100 δ *= Scalar(10);
101 γ *= Scalar(10);
102 }
103 } else if (inertia.negative > ideal_inertia.negative) {
104 // If there's too many negative eigenvalues, increase δ
105 δ *= Scalar(10);
106 } else if (inertia.positive > ideal_inertia.positive) {
107 // If there's too many positive eigenvalues, increase γ
108 γ = γ == Scalar(0) ? Scalar(1e-10) : γ * Scalar(10);
109 }
110 } else {
111 // If the decomposition failed, increase δ and γ
112 δ *= Scalar(10);
113 γ = γ == Scalar(0) ? Scalar(1e-10) : γ * Scalar(10);
114 }
115
116 // If the Hessian perturbation is too high, report failure. This can be
117 // caused by ill-conditioning.
118 if (δ > Scalar(1e20) || γ > Scalar(1e20)) {
119 m_info = Eigen::NumericalIssue;
120 m_prev_δ = δ;
121 m_prev_γ = γ;
122 return *this;
123 }
124 }
125 }
126
131 template <typename Rhs>
132 DenseVector solve(const Eigen::MatrixBase<Rhs>& rhs) const {
133 return m_solver.solve(rhs);
134 }
135
140 template <typename Rhs>
141 DenseVector solve(const Eigen::SparseMatrixBase<Rhs>& rhs) const {
142 return m_solver.solve(rhs.toDense());
143 }
144
148 Scalar hessian_regularization() const { return m_prev_δ; }
149
153 Scalar constraint_jacobian_regularization() const { return m_prev_γ; }
154
155 private:
156 using Solver = Eigen::LDLT<DenseMatrix>;
157
158 Solver m_solver;
159
160 Eigen::ComputationInfo m_info = Eigen::Success;
161
163 int m_num_decision_variables = 0;
164
166 int m_num_equality_constraints = 0;
167
169 Scalar m_γ_min{1e-10};
170
172 Inertia ideal_inertia{m_num_decision_variables, m_num_equality_constraints,
173 0};
174
176 Scalar m_prev_δ{0};
177
179 Scalar m_prev_γ{0};
180
189 DenseMatrix regularization(Scalar δ, Scalar γ) const {
190 DenseVector vec{m_num_decision_variables + m_num_equality_constraints};
191 vec.segment(0, m_num_decision_variables).setConstant(δ);
192 vec.segment(m_num_decision_variables, m_num_equality_constraints)
193 .setConstant(-γ);
194
195 return vec.asDiagonal().toDenseMatrix();
196 }
197};
198
199} // namespace slp
Definition dense_regularized_ldlt.hpp:19
DenseVector solve(const Eigen::MatrixBase< Rhs > &rhs) const
Definition dense_regularized_ldlt.hpp:132
DenseRegularizedLDLT(int num_decision_variables, int num_equality_constraints, Scalar γ_min)
Definition dense_regularized_ldlt.hpp:43
Eigen::Vector< Scalar, Eigen::Dynamic > DenseVector
Type alias for dense vector.
Definition dense_regularized_ldlt.hpp:24
Scalar hessian_regularization() const
Definition dense_regularized_ldlt.hpp:148
Eigen::Matrix< Scalar, Eigen::Dynamic, Eigen::Dynamic > DenseMatrix
Type alias for dense matrix.
Definition dense_regularized_ldlt.hpp:22
DenseRegularizedLDLT(int num_decision_variables, int num_equality_constraints)
Definition dense_regularized_ldlt.hpp:32
Eigen::ComputationInfo info() const
Definition dense_regularized_ldlt.hpp:52
DenseRegularizedLDLT & compute(const DenseMatrix &lhs)
Definition dense_regularized_ldlt.hpp:58
Scalar constraint_jacobian_regularization() const
Definition dense_regularized_ldlt.hpp:153
DenseVector solve(const Eigen::SparseMatrixBase< Rhs > &rhs) const
Definition dense_regularized_ldlt.hpp:141
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