Sleipnir C++ API
Loading...
Searching...
No Matches
regularized_ldlt.hpp
1// Copyright (c) Sleipnir contributors
2
3#pragma once
4
5#include <Eigen/Core>
6#include <Eigen/SparseCore>
7
8#include "sleipnir/optimization/solver/util/dense_regularized_ldlt.hpp"
9#include "sleipnir/optimization/solver/util/sparse_regularized_ldlt.hpp"
10
11namespace slp {
12
16template <typename Scalar>
18 public:
20 using DenseVector = Eigen::Vector<Scalar, Eigen::Dynamic>;
22 using SparseMatrix = Eigen::SparseMatrix<Scalar>;
23
31 RegularizedLDLT(bool use_sparse_solver, int num_decision_variables,
32 int num_equality_constraints)
33 : m_use_sparse_solver{use_sparse_solver},
34 m_sparse_solver{num_decision_variables, num_equality_constraints},
35 m_dense_solver{num_decision_variables, num_equality_constraints} {}
36
45 RegularizedLDLT(bool use_sparse_solver, int num_decision_variables,
46 int num_equality_constraints, Scalar γ_min)
47 : m_use_sparse_solver{use_sparse_solver},
48 m_sparse_solver{num_decision_variables, num_equality_constraints,
49 γ_min},
50 m_dense_solver{num_decision_variables, num_equality_constraints,
51 γ_min} {}
52
56 Eigen::ComputationInfo info() const {
57 if (m_use_sparse_solver) {
58 return m_sparse_solver.info();
59 } else {
60 return m_dense_solver.info();
61 }
62 }
63
73 if (m_use_sparse_solver) {
74 m_sparse_solver.compute(lhs);
75 } else {
76 m_dense_solver.compute(lhs);
77 }
78
79 return *this;
80 }
81
86 template <typename Rhs>
87 DenseVector solve(const Eigen::MatrixBase<Rhs>& rhs) const {
88 if (m_use_sparse_solver) {
89 return m_sparse_solver.solve(rhs);
90 } else {
91 return m_dense_solver.solve(rhs);
92 }
93 }
94
99 template <typename Rhs>
100 DenseVector solve(const Eigen::SparseMatrixBase<Rhs>& rhs) const {
101 if (m_use_sparse_solver) {
102 return m_sparse_solver.solve(rhs);
103 } else {
104 return m_dense_solver.solve(rhs);
105 }
106 }
107
111 Scalar hessian_regularization() const {
112 if (m_use_sparse_solver) {
113 return m_sparse_solver.hessian_regularization();
114 } else {
115 return m_dense_solver.hessian_regularization();
116 }
117 }
118
123 if (m_use_sparse_solver) {
124 return m_sparse_solver.constraint_jacobian_regularization();
125 } else {
126 return m_dense_solver.constraint_jacobian_regularization();
127 }
128 }
129
130 private:
131 bool m_use_sparse_solver;
132 SparseRegularizedLDLT<Scalar> m_sparse_solver;
133 DenseRegularizedLDLT<Scalar> m_dense_solver;
134};
135
136} // namespace slp
Definition intrusive_shared_ptr.hpp:27
Definition regularized_ldlt.hpp:17
Eigen::SparseMatrix< Scalar > SparseMatrix
Type alias for sparse matrix.
Definition regularized_ldlt.hpp:22
RegularizedLDLT(bool use_sparse_solver, int num_decision_variables, int num_equality_constraints)
Definition regularized_ldlt.hpp:31
Eigen::Vector< Scalar, Eigen::Dynamic > DenseVector
Type alias for dense vector.
Definition regularized_ldlt.hpp:20
Scalar constraint_jacobian_regularization() const
Definition regularized_ldlt.hpp:122
RegularizedLDLT & compute(const SparseMatrix &lhs)
Definition regularized_ldlt.hpp:72
Scalar hessian_regularization() const
Definition regularized_ldlt.hpp:111
Eigen::ComputationInfo info() const
Definition regularized_ldlt.hpp:56
DenseVector solve(const Eigen::MatrixBase< Rhs > &rhs) const
Definition regularized_ldlt.hpp:87
RegularizedLDLT(bool use_sparse_solver, int num_decision_variables, int num_equality_constraints, Scalar γ_min)
Definition regularized_ldlt.hpp:45
DenseVector solve(const Eigen::SparseMatrixBase< Rhs > &rhs) const
Definition regularized_ldlt.hpp:100