5#include <Eigen/Cholesky>
8#include "sleipnir/optimization/solver/util/inertia.hpp"
18template <
typename Scalar>
22 using DenseMatrix = Eigen::Matrix<Scalar, Eigen::Dynamic, Eigen::Dynamic>;
33 : m_num_decision_variables{num_decision_variables},
34 m_num_equality_constraints{num_equality_constraints} {}
45 : m_num_decision_variables{num_decision_variables},
46 m_num_equality_constraints{num_equality_constraints},
52 Eigen::ComputationInfo
info()
const {
return m_info; }
61 m_info = m_solver.compute(
lhs).info();
63 if (m_info == Eigen::Success) {
64 auto D = m_solver.vectorD();
69 (
D.cwiseAbs().array() >= Scalar(1
e-4)).all()) {
80 Scalar
δ = m_prev_δ == Scalar(0) ? Scalar(1
e-4) : m_prev_δ / Scalar(2);
84 m_info = m_solver.compute(
lhs + regularization(
δ,
γ)).info();
86 if (m_info == Eigen::Success) {
108 γ =
γ == Scalar(0) ? Scalar(1
e-10) :
γ * Scalar(10);
113 γ =
γ == Scalar(0) ? Scalar(1
e-10) :
γ * Scalar(10);
118 if (
δ > Scalar(1
e20) ||
γ > Scalar(1
e20)) {
119 m_info = Eigen::NumericalIssue;
131 template <
typename Rhs>
133 return m_solver.solve(
rhs);
140 template <
typename Rhs>
142 return m_solver.solve(
rhs.toDense());
156 using Solver = Eigen::LDLT<DenseMatrix>;
160 Eigen::ComputationInfo m_info = Eigen::Success;
163 int m_num_decision_variables = 0;
166 int m_num_equality_constraints = 0;
169 Scalar m_γ_min{1
e-10};
172 Inertia ideal_inertia{m_num_decision_variables, m_num_equality_constraints,
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)
195 return vec.asDiagonal().toDenseMatrix();
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