LLT.hpp
Go to the documentation of this file.
1 /*
2  * Copyright 2020-2024 INRIA
3  */
4 
5 #ifndef __eigenpy_decompositions_llt_hpp__
6 #define __eigenpy_decompositions_llt_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<LLTSolverVisitor<_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::LLT<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 LLT factorization from a given matrix."))
39 
41 
42  .def("matrixL", &matrixL, bp::arg("self"),
43  "Returns the lower triangular matrix L.")
44  .def("matrixU", &matrixU, bp::arg("self"),
45  "Returns the upper triangular matrix U.")
46  .def("matrixLLT", &Solver::matrixLLT, bp::arg("self"),
47  "Returns the LLT decomposition matrix.",
48  bp::return_internal_reference<>())
49 
50 #if EIGEN_VERSION_AT_LEAST(3, 3, 90)
51  .def("rankUpdate",
52  (Solver & (Solver::*)(const VectorXs &, const RealScalar &)) &
53  Solver::template rankUpdate<VectorXs>,
54  bp::args("self", "vector", "sigma"), bp::return_self<>())
55 #else
56  .def("rankUpdate",
57  (Solver(Solver::*)(const VectorXs &, const RealScalar &)) &
58  Solver::template rankUpdate<VectorXs>,
59  bp::args("self", "vector", "sigma"))
60 #endif
61 
62 #if EIGEN_VERSION_AT_LEAST(3, 3, 0)
63  .def("adjoint", &Solver::adjoint, bp::arg("self"),
64  "Returns the adjoint, that is, a reference to the decomposition "
65  "itself as if the underlying matrix is self-adjoint.",
66  bp::return_self<>())
67 #endif
68 
69  .def(
70  "compute",
71  (Solver & (Solver::*)(const Eigen::EigenBase<MatrixType> &matrix)) &
72  Solver::compute,
73  bp::args("self", "matrix"), "Computes the LLT of given matrix.",
74  bp::return_self<>())
75 
76  .def("info", &Solver::info, bp::arg("self"),
77  "NumericalIssue if the input contains INF or NaN values or "
78  "overflow occured. Returns Success otherwise.")
79 #if EIGEN_VERSION_AT_LEAST(3, 3, 0)
80  .def("rcond", &Solver::rcond, bp::arg("self"),
81  "Returns an estimate of the reciprocal condition number of the "
82  "matrix.")
83 #endif
84  .def("reconstructedMatrix", &Solver::reconstructedMatrix,
85  bp::arg("self"),
86  "Returns the matrix represented by the decomposition, i.e., it "
87  "returns the product: L L^*. This function is provided for debug "
88  "purpose.")
89  .def("solve", &solve<VectorXs>, bp::args("self", "b"),
90  "Returns the solution x of A x = b using the current "
91  "decomposition of A.")
92  .def("solve", &solve<MatrixXs>, bp::args("self", "B"),
93  "Returns the solution X of A X = B using the current "
94  "decomposition of A where B is a right hand side matrix.");
95  }
96 
97  static void expose() {
98  static const std::string classname =
100  expose(classname);
101  }
102 
103  static void expose(const std::string &name) {
104  bp::class_<Solver>(
105  name.c_str(),
106  "Standard Cholesky decomposition (LL^T) of a matrix and associated "
107  "features.\n\n"
108  "This class performs a LL^T Cholesky decomposition of a symmetric, "
109  "positive definite matrix A such that A = LL^* = U^*U, where L is "
110  "lower triangular.\n\n"
111  "While the Cholesky decomposition is particularly useful to solve "
112  "selfadjoint problems like D^*D x = b, for that purpose, we recommend "
113  "the Cholesky decomposition without square root which is more stable "
114  "and even faster. Nevertheless, this standard Cholesky decomposition "
115  "remains useful in many other situations like generalised eigen "
116  "problems with hermitian matrices.",
117  bp::no_init)
118  .def(IdVisitor<Solver>())
119  .def(LLTSolverVisitor());
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 
126  template <typename MatrixOrVector>
127  static MatrixOrVector solve(const Solver &self, const MatrixOrVector &vec) {
128  return self.solve(vec);
129  }
130 };
131 
132 } // namespace eigenpy
133 
134 #endif // ifndef __eigenpy_decompositions_llt_hpp__
eigenpy::LLTSolverVisitor::solve
static MatrixOrVector solve(const Solver &self, const MatrixOrVector &vec)
Definition: LLT.hpp:127
scalar-name.hpp
eigenpy::LLTSolverVisitor
Definition: LLT.hpp:18
test_matrix.vec
vec
Definition: test_matrix.py:180
eigenpy::LLTSolverVisitor::visit
void visit(PyClass &cl) const
Definition: LLT.hpp:31
eigenpy::LLTSolverVisitor::expose
static void expose(const std::string &name)
Definition: LLT.hpp:103
eigenpy::LLTSolverVisitor::RealScalar
MatrixType::RealScalar RealScalar
Definition: LLT.hpp:22
eigenpy
Definition: alignment.hpp:14
eigenpy::LLTSolverVisitor::MatrixType
_MatrixType MatrixType
Definition: LLT.hpp:20
eigenpy::scalar_name::shortname
static std::string shortname()
eigenpy::EigenBaseVisitor
Definition: EigenBase.hpp:13
eigenpy::LLTSolverVisitor::MatrixXs
Eigen::Matrix< Scalar, Eigen::Dynamic, Eigen::Dynamic, MatrixType::Options > MatrixXs
Definition: LLT.hpp:27
eigenpy::IdVisitor
Add the Python method id to retrieving a unique id for a given object exposed with Boost....
Definition: id.hpp:18
eigenpy::LLTSolverVisitor::Solver
Eigen::LLT< MatrixType > Solver
Definition: LLT.hpp:28
eigenpy::LLTSolverVisitor::matrixL
static MatrixType matrixL(const Solver &self)
Definition: LLT.hpp:123
setup.name
name
Definition: setup.in.py:179
eigenpy.hpp
eigenpy::LLTSolverVisitor::expose
static void expose()
Definition: LLT.hpp:97
eigenpy::LLTSolverVisitor::Scalar
MatrixType::Scalar Scalar
Definition: LLT.hpp:21
eigenpy::LLTSolverVisitor::matrixU
static MatrixType matrixU(const Solver &self)
Definition: LLT.hpp:124
EigenBase.hpp
eigenpy::LLTSolverVisitor::VectorXs
Eigen::Matrix< Scalar, Eigen::Dynamic, 1, MatrixType::Options > VectorXs
Definition: LLT.hpp:24


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