5#include <Eigen/Cholesky>
7#include <Eigen/SparseCholesky>
8#include <Eigen/SparseCore>
19template <
typename Scalar>
23 using DenseMatrix = Eigen::Matrix<Scalar, Eigen::Dynamic, Eigen::Dynamic>;
36 : m_num_decision_variables{num_decision_variables},
37 m_num_equality_constraints{num_equality_constraints} {}
48 : m_num_decision_variables{num_decision_variables},
49 m_num_equality_constraints{num_equality_constraints},
55 Eigen::ComputationInfo
info()
const {
return m_info; }
66 m_is_sparse = lhs.nonZeros() < 0.25 * lhs.size();
68 m_info = m_is_sparse ? compute_sparse(lhs).
info()
69 : m_dense_solver.compute(lhs).info();
71 if (m_info == Eigen::Success) {
73 m_is_sparse ? m_sparse_solver.vectorD() : m_dense_solver.vectorD();
77 if (
Inertia{D} == ideal_inertia &&
78 (D.cwiseAbs().array() >= Scalar(1e-4)).all()) {
88 Scalar δ = m_prev_δ == Scalar(0) ? Scalar(1e-4) : m_prev_δ / Scalar(2);
99 m_info = compute_sparse(lhs + regularization(δ, γ)).info();
100 if (m_info == Eigen::Success) {
101 inertia =
Inertia{m_sparse_solver.vectorD()};
104 m_info = m_dense_solver.compute(lhs + regularization(δ, γ)).info();
105 if (m_info == Eigen::Success) {
106 inertia =
Inertia{m_dense_solver.vectorD()};
110 if (m_info == Eigen::Success) {
111 if (inertia == ideal_inertia) {
115 }
else if (inertia.
zero > 0) {
120 }
else if (inertia.
negative > ideal_inertia.negative) {
124 }
else if (inertia.
positive > ideal_inertia.positive) {
127 γ = γ == Scalar(0) ? Scalar(1e-10) : γ * Scalar(10);
138 if (δ > Scalar(1e20) || γ > Scalar(1e20)) {
139 m_info = Eigen::NumericalIssue;
149 template <
typename Rhs>
152 return m_sparse_solver.solve(rhs);
154 return m_dense_solver.solve(rhs);
162 template <
typename Rhs>
165 return m_sparse_solver.solve(rhs);
167 return m_dense_solver.solve(rhs.toDense());
177 using SparseSolver = Eigen::SimplicialLDLT<SparseMatrix>;
178 using DenseSolver = Eigen::LDLT<DenseMatrix>;
180 SparseSolver m_sparse_solver;
181 DenseSolver m_dense_solver;
182 bool m_is_sparse =
true;
184 Eigen::ComputationInfo m_info = Eigen::Success;
187 int m_num_decision_variables = 0;
190 int m_num_equality_constraints = 0;
193 Scalar m_γ_min{1e-10};
196 Inertia ideal_inertia{m_num_decision_variables, m_num_equality_constraints,
203 int m_non_zeros = -1;
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;
217 m_sparse_solver.factorize(lhs);
219 return m_sparse_solver;
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)
Represents the inertia of a matrix (the number of positive, negative, and zero eigenvalues).
Definition inertia.hpp:11
int zero
The number of zero eigenvalues.
Definition inertia.hpp:18
int positive
The number of positive eigenvalues.
Definition inertia.hpp:14
int negative
The number of negative eigenvalues.
Definition inertia.hpp:16
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)
Computes the regularized LDLT factorization of a matrix.
Definition regularized_ldlt.hpp:61
Scalar hessian_regularization() const
Returns the Hessian regularization factor.
Definition regularized_ldlt.hpp:174
RegularizedLDLT(int num_decision_variables, int num_equality_constraints, Scalar γ_min)
Constructs a RegularizedLDLT instance.
Definition regularized_ldlt.hpp:46
Eigen::ComputationInfo info() const
Reports whether previous computation was successful.
Definition regularized_ldlt.hpp:55
DenseVector solve(const Eigen::MatrixBase< Rhs > &rhs) const
Solves the system of equations using a regularized LDLT factorization.
Definition regularized_ldlt.hpp:150
DenseVector solve(const Eigen::SparseMatrixBase< Rhs > &rhs) const
Solves the system of equations using a regularized LDLT factorization.
Definition regularized_ldlt.hpp:163
RegularizedLDLT(int num_decision_variables, int num_equality_constraints)
Constructs a RegularizedLDLT instance.
Definition regularized_ldlt.hpp:35
Definition to_underlying.hpp:7