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


eigenpy
Author(s): Justin Carpentier, Nicolas Mansard
autogenerated on Fri Jun 14 2024 02:15:58