9#include <Eigen/SparseCholesky>
10#include <Eigen/SparseCore>
12#include "sleipnir/optimization/solver/util/inertia.hpp"
20template <
typename Scalar>
35 int num_equality_constraints)
36 : m_num_decision_variables{num_decision_variables},
37 m_num_equality_constraints{num_equality_constraints} {}
47 int num_equality_constraints, Scalar
γ_min)
48 : m_num_decision_variables{num_decision_variables},
49 m_num_equality_constraints{num_equality_constraints},
55 Eigen::ComputationInfo
info()
const {
return m_info; }
69 if (!m_analyzed_pattern) {
71 m_analyzed_pattern =
true;
75 m_info = m_solver.info();
77 if (m_info == Eigen::Success) {
78 auto D = m_solver.vectorD();
83 (
D.cwiseAbs().array() >= Scalar(1
e-4)).all()) {
95 Scalar
δ = m_prev_δ == Scalar(0)
97 : std::max(m_prev_δ / Scalar(2),
98 std::numeric_limits<Scalar>::epsilon());
105 m_solver.factorize(
lhs + regularization(
δ,
γ));
106 m_info = m_solver.info();
108 if (m_info == Eigen::Success) {
111 if (
inertia == ideal_inertia) {
117 if (
γ == Scalar(0)) {
134 γ =
γ == Scalar(0) ? Scalar(1
e-10) :
γ * Scalar(10);
140 γ =
γ == Scalar(0) ? Scalar(1
e-10) :
γ * Scalar(10);
145 if (
δ > Scalar(1
e20) ||
γ > Scalar(1
e20)) {
146 m_info = Eigen::NumericalIssue;
158 template <
typename Rhs>
160 return m_solver.solve(
rhs);
167 template <
typename Rhs>
169 return m_solver.solve(
rhs);
183 using Solver = Eigen::SimplicialLDLT<SparseMatrix>;
186 bool m_analyzed_pattern =
false;
188 Eigen::ComputationInfo m_info = Eigen::Success;
191 int m_num_decision_variables = 0;
194 int m_num_equality_constraints = 0;
197 Scalar m_γ_min{1
e-10};
200 Inertia ideal_inertia{m_num_decision_variables, m_num_equality_constraints,
218 DenseVector vec{m_num_decision_variables + m_num_equality_constraints};
219 vec.segment(0, m_num_decision_variables).setConstant(δ);
220 vec.segment(m_num_decision_variables, m_num_equality_constraints)
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:21
Eigen::Vector< Scalar, Eigen::Dynamic > DenseVector
Type alias for dense vector.
Definition sparse_regularized_ldlt.hpp:24
Eigen::ComputationInfo info() const
Definition sparse_regularized_ldlt.hpp:55
Scalar constraint_jacobian_regularization() const
Definition sparse_regularized_ldlt.hpp:180
SparseRegularizedLDLT & compute(const SparseMatrix &lhs)
Definition sparse_regularized_ldlt.hpp:64
DenseVector solve(const Eigen::MatrixBase< Rhs > &rhs) const
Definition sparse_regularized_ldlt.hpp:159
DenseVector solve(const Eigen::SparseMatrixBase< Rhs > &rhs) const
Definition sparse_regularized_ldlt.hpp:168
SparseRegularizedLDLT(int num_decision_variables, int num_equality_constraints, Scalar γ_min)
Definition sparse_regularized_ldlt.hpp:46
SparseRegularizedLDLT(int num_decision_variables, int num_equality_constraints)
Definition sparse_regularized_ldlt.hpp:34
Scalar hessian_regularization() const
Definition sparse_regularized_ldlt.hpp:175
Eigen::SparseMatrix< Scalar > SparseMatrix
Type alias for sparse matrix.
Definition sparse_regularized_ldlt.hpp:26