LDLT.hpp
Go to the documentation of this file.
1 /*
2  * Copyright 2020-2021 INRIA
3  */
4 
5 #ifndef __eigenpy_decomposition_ldlt_hpp__
6 #define __eigenpy_decomposition_ldlt_hpp__
7 
8 #include <Eigen/Cholesky>
9 #include <Eigen/Core>
10 
11 #include "eigenpy/eigenpy.hpp"
13 
14 namespace eigenpy {
15 
16 template <typename _MatrixType>
18  : public boost::python::def_visitor<LDLTSolverVisitor<_MatrixType> > {
19  typedef _MatrixType MatrixType;
20  typedef typename MatrixType::Scalar Scalar;
21  typedef typename MatrixType::RealScalar RealScalar;
22  typedef Eigen::Matrix<Scalar, Eigen::Dynamic, 1, MatrixType::Options>
24  typedef Eigen::Matrix<Scalar, Eigen::Dynamic, Eigen::Dynamic,
25  MatrixType::Options>
27  typedef Eigen::LDLT<MatrixType> Solver;
28 
29  template <class PyClass>
30  void visit(PyClass &cl) const {
31  cl.def(bp::init<>(bp::arg("self"), "Default constructor"))
32  .def(bp::init<Eigen::DenseIndex>(
33  bp::args("self", "size"),
34  "Default constructor with memory preallocation"))
35  .def(bp::init<MatrixType>(
36  bp::args("self", "matrix"),
37  "Constructs a LDLT factorization from a given matrix."))
38 
39  .def("isNegative", &Solver::isNegative, bp::arg("self"),
40  "Returns true if the matrix is negative (semidefinite).")
41  .def("isPositive", &Solver::isPositive, bp::arg("self"),
42  "Returns true if the matrix is positive (semidefinite).")
43 
44  .def("matrixL", &matrixL, bp::arg("self"),
45  "Returns the lower triangular matrix L.")
46  .def("matrixU", &matrixU, bp::arg("self"),
47  "Returns the upper triangular matrix U.")
48  .def("vectorD", &vectorD, bp::arg("self"),
49  "Returns the coefficients of the diagonal matrix D.")
50  .def("transpositionsP", &transpositionsP, bp::arg("self"),
51  "Returns the permutation matrix P.")
52 
53  .def("matrixLDLT", &Solver::matrixLDLT, bp::arg("self"),
54  "Returns the LDLT decomposition matrix.",
55  bp::return_internal_reference<>())
56 
57  .def("rankUpdate",
58  (Solver & (Solver::*)(const Eigen::MatrixBase<VectorXs> &,
59  const RealScalar &)) &
60  Solver::template rankUpdate<VectorXs>,
61  bp::args("self", "vector", "sigma"), bp::return_self<>())
62 
63 #if EIGEN_VERSION_AT_LEAST(3, 3, 0)
64  .def("adjoint", &Solver::adjoint, bp::arg("self"),
65  "Returns the adjoint, that is, a reference to the decomposition "
66  "itself as if the underlying matrix is self-adjoint.",
67  bp::return_self<>())
68 #endif
69 
70  .def(
71  "compute",
72  (Solver & (Solver::*)(const Eigen::EigenBase<MatrixType> &matrix)) &
73  Solver::compute,
74  bp::args("self", "matrix"), "Computes the LDLT of given matrix.",
75  bp::return_self<>())
76 
77  .def("info", &Solver::info, bp::arg("self"),
78  "NumericalIssue if the input contains INF or NaN values or "
79  "overflow occured. Returns Success otherwise.")
80 #if EIGEN_VERSION_AT_LEAST(3, 3, 0)
81  .def("rcond", &Solver::rcond, bp::arg("self"),
82  "Returns an estimate of the reciprocal condition number of the "
83  "matrix.")
84 #endif
85  .def("reconstructedMatrix", &Solver::reconstructedMatrix,
86  bp::arg("self"),
87  "Returns the matrix represented by the decomposition, i.e., it "
88  "returns the product: L L^*. This function is provided for debug "
89  "purpose.")
90  .def("solve", &solve<VectorXs>, bp::args("self", "b"),
91  "Returns the solution x of A x = b using the current "
92  "decomposition of A.")
93  .def("solve", &solve<MatrixXs>, bp::args("self", "B"),
94  "Returns the solution X of A X = B using the current "
95  "decomposition of A where B is a right hand side matrix.")
96 
97  .def("setZero", &Solver::setZero, bp::arg("self"),
98  "Clear any existing decomposition.");
99  }
100 
101  static void expose() {
102  static const std::string classname =
104  expose(classname);
105  }
106 
107  static void expose(const std::string &name) {
108  bp::class_<Solver>(
109  name.c_str(),
110  "Robust Cholesky decomposition of a matrix with pivoting.\n\n"
111  "Perform a robust Cholesky decomposition of a positive semidefinite or "
112  "negative semidefinite matrix $ A $ such that $ A = P^TLDL^*P $, where "
113  "P is a permutation matrix, L is lower triangular with a unit diagonal "
114  "and D is a diagonal matrix.\n\n"
115  "The decomposition uses pivoting to ensure stability, so that L will "
116  "have zeros in the bottom right rank(A) - n submatrix. Avoiding the "
117  "square root on D also stabilizes the computation.",
118  bp::no_init)
119  .def(LDLTSolverVisitor());
120  }
121 
122  private:
123  static MatrixType matrixL(const Solver &self) { return self.matrixL(); }
124  static MatrixType matrixU(const Solver &self) { return self.matrixU(); }
125  static VectorXs vectorD(const Solver &self) { return self.vectorD(); }
126 
127  static MatrixType transpositionsP(const Solver &self) {
128  return self.transpositionsP() *
129  MatrixType::Identity(self.matrixL().rows(), self.matrixL().rows());
130  }
131 
132  template <typename MatrixOrVector>
133  static MatrixOrVector solve(const Solver &self, const MatrixOrVector &vec) {
134  return self.solve(vec);
135  }
136 };
137 
138 } // namespace eigenpy
139 
140 #endif // ifndef __eigenpy_decomposition_ldlt_hpp__
eigenpy::LDLTSolverVisitor::transpositionsP
static MatrixType transpositionsP(const Solver &self)
Definition: LDLT.hpp:127
eigenpy::LDLTSolverVisitor::expose
static void expose(const std::string &name)
Definition: LDLT.hpp:107
eigenpy::LDLTSolverVisitor::RealScalar
MatrixType::RealScalar RealScalar
Definition: LDLT.hpp:21
eigenpy::LDLTSolverVisitor::vectorD
static VectorXs vectorD(const Solver &self)
Definition: LDLT.hpp:125
eigenpy::LDLTSolverVisitor::Scalar
MatrixType::Scalar Scalar
Definition: LDLT.hpp:20
eigenpy::LDLTSolverVisitor::VectorXs
Eigen::Matrix< Scalar, Eigen::Dynamic, 1, MatrixType::Options > VectorXs
Definition: LDLT.hpp:23
eigenpy::LDLTSolverVisitor::solve
static MatrixOrVector solve(const Solver &self, const MatrixOrVector &vec)
Definition: LDLT.hpp:133
eigenpy::LDLTSolverVisitor::MatrixType
_MatrixType MatrixType
Definition: LDLT.hpp:19
scalar-name.hpp
test_complex.rows
int rows
Definition: test_complex.py:6
test_matrix.vec
vec
Definition: test_matrix.py:180
eigenpy::LDLTSolverVisitor::visit
void visit(PyClass &cl) const
Definition: LDLT.hpp:30
eigenpy::LDLTSolverVisitor::matrixL
static MatrixType matrixL(const Solver &self)
Definition: LDLT.hpp:123
eigenpy::LDLTSolverVisitor::Solver
Eigen::LDLT< MatrixType > Solver
Definition: LDLT.hpp:27
eigenpy
Definition: alignment.hpp:14
eigenpy::LDLTSolverVisitor
Definition: LDLT.hpp:17
eigenpy::scalar_name::shortname
static std::string shortname()
eigenpy::LDLTSolverVisitor::matrixU
static MatrixType matrixU(const Solver &self)
Definition: LDLT.hpp:124
setup.name
name
Definition: setup.in.py:179
eigenpy.hpp
setZero
void setZero(std::vector< MatType, Eigen::aligned_allocator< MatType > > &Ms)
Definition: std_vector.cpp:27
eigenpy::LDLTSolverVisitor::expose
static void expose()
Definition: LDLT.hpp:101
eigenpy::LDLTSolverVisitor::MatrixXs
Eigen::Matrix< Scalar, Eigen::Dynamic, Eigen::Dynamic, MatrixType::Options > MatrixXs
Definition: LDLT.hpp:26


eigenpy
Author(s): Justin Carpentier, Nicolas Mansard
autogenerated on Tue Jan 23 2024 03:15:01