8#include <Eigen/Cholesky>
11#include "sleipnir/optimization/solver/util/inertia.hpp"
19template <
typename Scalar>
23 using DenseMatrix = Eigen::Matrix<Scalar, Eigen::Dynamic, Eigen::Dynamic>;
34 : m_num_decision_variables{num_decision_variables},
35 m_num_equality_constraints{num_equality_constraints} {}
46 : m_num_decision_variables{num_decision_variables},
47 m_num_equality_constraints{num_equality_constraints},
53 Eigen::ComputationInfo
info()
const {
return m_info; }
60 m_info = m_solver.compute(
lhs).info();
62 if (m_info == Eigen::Success) {
63 auto D = m_solver.vectorD();
68 (
D.cwiseAbs().array() >= Scalar(1
e-4)).all()) {
80 Scalar
δ = m_prev_δ == Scalar(0)
82 : std::max(m_prev_δ / Scalar(2),
83 std::numeric_limits<Scalar>::epsilon());
90 m_info = m_solver.compute(
lhs + regularization(
δ,
γ)).info();
92 if (m_info == Eigen::Success) {
101 if (
γ == Scalar(0)) {
118 γ =
γ == Scalar(0) ? Scalar(1
e-10) :
γ * Scalar(10);
124 γ =
γ == Scalar(0) ? Scalar(1
e-10) :
γ * Scalar(10);
129 if (
δ > Scalar(1
e20) ||
γ > Scalar(1
e20)) {
130 m_info = Eigen::NumericalIssue;
142 template <
typename Rhs>
144 return m_solver.solve(
rhs);
151 template <
typename Rhs>
153 return m_solver.solve(
rhs.toDense());
167 using Solver = Eigen::LDLT<DenseMatrix>;
171 Eigen::ComputationInfo m_info = Eigen::Success;
174 int m_num_decision_variables = 0;
177 int m_num_equality_constraints = 0;
180 Scalar m_γ_min{1
e-10};
183 Inertia ideal_inertia{m_num_decision_variables, m_num_equality_constraints,
200 DenseMatrix regularization(Scalar δ, Scalar γ)
const {
201 DenseVector vec{m_num_decision_variables + m_num_equality_constraints};
202 vec.segment(0, m_num_decision_variables).setConstant(δ);
203 vec.segment(m_num_decision_variables, m_num_equality_constraints)
206 return vec.asDiagonal().toDenseMatrix();
Definition dense_regularized_ldlt.hpp:20
DenseVector solve(const Eigen::MatrixBase< Rhs > &rhs) const
Definition dense_regularized_ldlt.hpp:143
DenseRegularizedLDLT(int num_decision_variables, int num_equality_constraints, Scalar γ_min)
Definition dense_regularized_ldlt.hpp:44
Eigen::Vector< Scalar, Eigen::Dynamic > DenseVector
Type alias for dense vector.
Definition dense_regularized_ldlt.hpp:25
Scalar hessian_regularization() const
Definition dense_regularized_ldlt.hpp:159
Eigen::Matrix< Scalar, Eigen::Dynamic, Eigen::Dynamic > DenseMatrix
Type alias for dense matrix.
Definition dense_regularized_ldlt.hpp:23
DenseRegularizedLDLT(int num_decision_variables, int num_equality_constraints)
Definition dense_regularized_ldlt.hpp:33
Eigen::ComputationInfo info() const
Definition dense_regularized_ldlt.hpp:53
DenseRegularizedLDLT & compute(const DenseMatrix &lhs)
Definition dense_regularized_ldlt.hpp:59
Scalar constraint_jacobian_regularization() const
Definition dense_regularized_ldlt.hpp:164
DenseVector solve(const Eigen::SparseMatrixBase< Rhs > &rhs) const
Definition dense_regularized_ldlt.hpp:152
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