5#include <Eigen/Cholesky>
7#include <Eigen/SparseCholesky>
8#include <Eigen/SparseCore>
10#include "sleipnir/optimization/solver/util/inertia.hpp"
19template <
typename Scalar>
23 using DenseMatrix = Eigen::Matrix<Scalar, Eigen::Dynamic, Eigen::Dynamic>;
42 Eigen::ComputationInfo
info()
const {
return m_info; }
53 m_is_sparse =
lhs.nonZeros() < 0.25 *
lhs.size();
55 m_info = m_is_sparse ? compute_sparse(
lhs).info()
56 : m_dense_solver.compute(
lhs).info();
60 if (m_info == Eigen::Success) {
62 :
Inertia{m_dense_solver.vectorD()};
75 Scalar
δ = m_prev_δ == Scalar(0) ? Scalar(1
e-4) : m_prev_δ / Scalar(2);
84 m_info = compute_sparse(
lhs + regularization(
δ,
γ)).info();
85 if (m_info == Eigen::Success) {
89 m_info = m_dense_solver.compute(
lhs + regularization(
δ,
γ)).info();
90 if (m_info == Eigen::Success) {
95 if (m_info == Eigen::Success) {
123 if (
δ > Scalar(1
e20) ||
γ > Scalar(1
e20)) {
124 m_info = Eigen::NumericalIssue;
134 template <
typename Rhs>
137 return m_sparse_solver.solve(
rhs);
139 return m_dense_solver.solve(
rhs);
147 template <
typename Rhs>
150 return m_sparse_solver.solve(
rhs);
152 return m_dense_solver.solve(
rhs.toDense());
162 using SparseSolver = Eigen::SimplicialLDLT<SparseMatrix>;
163 using DenseSolver = Eigen::LDLT<DenseMatrix>;
165 SparseSolver m_sparse_solver;
166 DenseSolver m_dense_solver;
167 bool m_is_sparse =
true;
169 Eigen::ComputationInfo m_info = Eigen::Success;
172 int m_num_decision_variables = 0;
175 int m_num_equality_constraints = 0;
178 Inertia ideal_inertia{m_num_decision_variables, m_num_equality_constraints,
185 int m_non_zeros = -1;
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;
199 m_sparse_solver.factorize(lhs);
201 return m_sparse_solver;
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)
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