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::*)(
58  const VectorXs &,
59  const RealScalar &))&Solver::template rankUpdate<VectorXs>,
60  bp::args("self", "vector", "sigma"))
61 #endif
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 LLT 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 
98  static void expose() {
99  static const std::string classname =
101  expose(classname);
102  }
103 
104  static void expose(const std::string &name) {
105  bp::class_<Solver>(
106  name.c_str(),
107  "Standard Cholesky decomposition (LL^T) of a matrix and associated "
108  "features.\n\n"
109  "This class performs a LL^T Cholesky decomposition of a symmetric, "
110  "positive definite matrix A such that A = LL^* = U^*U, where L is "
111  "lower triangular.\n\n"
112  "While the Cholesky decomposition is particularly useful to solve "
113  "selfadjoint problems like D^*D x = b, for that purpose, we recommend "
114  "the Cholesky decomposition without square root which is more stable "
115  "and even faster. Nevertheless, this standard Cholesky decomposition "
116  "remains useful in many other situations like generalised eigen "
117  "problems with hermitian matrices.",
118  bp::no_init)
119  .def(IdVisitor<Solver>())
120  .def(LLTSolverVisitor());
121  }
122 
123  private:
124  static MatrixType matrixL(const Solver &self) { return self.matrixL(); }
125  static MatrixType matrixU(const Solver &self) { return self.matrixU(); }
126 
127  template <typename MatrixOrVector>
128  static MatrixOrVector solve(const Solver &self, const MatrixOrVector &vec) {
129  return self.solve(vec);
130  }
131 };
132 
133 } // namespace eigenpy
134 
135 #endif // ifndef __eigenpy_decompositions_llt_hpp__
eigenpy::LLTSolverVisitor::solve
static MatrixOrVector solve(const Solver &self, const MatrixOrVector &vec)
Definition: LLT.hpp:128
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:104
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:124
setup.name
name
Definition: setup.in.py:179
eigenpy.hpp
eigenpy::LLTSolverVisitor::expose
static void expose()
Definition: LLT.hpp:98
eigenpy::LLTSolverVisitor::Scalar
MatrixType::Scalar Scalar
Definition: LLT.hpp:21
eigenpy::LLTSolverVisitor::matrixU
static MatrixType matrixU(const Solver &self)
Definition: LLT.hpp:125
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 Sat Nov 2 2024 02:14:45