eigenpy 3.12.0
Bindings between Numpy and Eigen using Boost.Python
Loading...
Searching...
No Matches
HouseholderQR.hpp
1/*
2 * Copyright 2024 INRIA
3 */
4
5#ifndef __eigenpy_decompositions_houselholder_qr_hpp__
6#define __eigenpy_decompositions_houselholder_qr_hpp__
7
8#include "eigenpy/eigenpy.hpp"
9#include "eigenpy/utils/scalar-name.hpp"
10
11#include <Eigen/QR>
12
13namespace eigenpy {
14
15template <typename _MatrixType>
17 : public boost::python::def_visitor<
18 HouseholderQRSolverVisitor<_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>
23 VectorXs;
24 typedef Eigen::Matrix<Scalar, Eigen::Dynamic, Eigen::Dynamic,
25 MatrixType::Options>
26 MatrixXs;
27 typedef Eigen::HouseholderQR<MatrixType> Solver;
28 typedef Solver Self;
29
30 template <class PyClass>
31 void visit(PyClass &cl) const {
32 cl.def(bp::init<>(bp::arg("self"),
33 "Default constructor.\n"
34 "The default constructor is useful in cases in which the "
35 "user intends to perform decompositions via "
36 "HouseholderQR.compute(matrix)"))
37 .def(bp::init<Eigen::DenseIndex, Eigen::DenseIndex>(
38 bp::args("self", "rows", "cols"),
39 "Default constructor with memory preallocation.\n"
40 "Like the default constructor but with preallocation of the "
41 "internal data according to the specified problem size. "))
42 .def(bp::init<MatrixType>(
43 bp::args("self", "matrix"),
44 "Constructs a QR factorization from a given matrix.\n"
45 "This constructor computes the QR factorization of the matrix "
46 "matrix by calling the method compute()."))
47
48 .def("absDeterminant", &Self::absDeterminant, bp::arg("self"),
49 "Returns the absolute value of the determinant of the matrix of "
50 "which *this is the QR decomposition.\n"
51 "It has only linear complexity (that is, O(n) where n is the "
52 "dimension of the square matrix) as the QR decomposition has "
53 "already been computed.\n"
54 "Note: This is only for square matrices.")
55 .def("logAbsDeterminant", &Self::logAbsDeterminant, bp::arg("self"),
56 "Returns the natural log of the absolute value of the determinant "
57 "of the matrix of which *this is the QR decomposition.\n"
58 "It has only linear complexity (that is, O(n) where n is the "
59 "dimension of the square matrix) as the QR decomposition has "
60 "already been computed.\n"
61 "Note: This is only for square matrices. This method is useful to "
62 "work around the risk of overflow/underflow that's inherent to "
63 "determinant computation.")
64
65 .def("matrixQR", &Self::matrixQR, bp::arg("self"),
66 "Returns the matrix where the Householder QR decomposition is "
67 "stored in a LAPACK-compatible way.",
68 bp::return_value_policy<bp::copy_const_reference>())
69
70 .def(
71 "compute",
72 (Solver & (Solver::*)(const Eigen::EigenBase<MatrixType> &matrix)) &
73 Solver::compute,
74 bp::args("self", "matrix"),
75 "Computes the QR factorization of given matrix.",
76 bp::return_self<>())
77
78 .def("solve", &solve<MatrixXs>, bp::args("self", "B"),
79 "Returns the solution X of A X = B using the current "
80 "decomposition of A where B is a right hand side matrix.");
81 }
82
83 static void expose() {
84 static const std::string classname =
85 "HouseholderQR" + scalar_name<Scalar>::shortname();
86 expose(classname);
87 }
88
89 static void expose(const std::string &name) {
90 bp::class_<Solver>(
91 name.c_str(),
92 "This class performs a QR decomposition of a matrix A into matrices Q "
93 "and R such that A=QR by using Householder transformations.\n"
94 "Here, Q a unitary matrix and R an upper triangular matrix. The result "
95 "is stored in a compact way compatible with LAPACK.\n"
96 "\n"
97 "Note that no pivoting is performed. This is not a rank-revealing "
98 "decomposition. If you want that feature, use FullPivHouseholderQR or "
99 "ColPivHouseholderQR instead.\n"
100 "\n"
101 "This Householder QR decomposition is faster, but less numerically "
102 "stable and less feature-full than FullPivHouseholderQR or "
103 "ColPivHouseholderQR.",
104 bp::no_init)
106 .def(IdVisitor<Solver>());
107 }
108
109 private:
110 template <typename MatrixOrVector>
111 static MatrixOrVector solve(const Solver &self, const MatrixOrVector &vec) {
112 return self.solve(vec);
113 }
114};
115
116} // namespace eigenpy
117
118#endif // ifndef __eigenpy_decompositions_houselholder_qr_hpp__
void expose()
Call the expose function of a given type T.
Definition expose.hpp:23
Add the Python method id to retrieving a unique id for a given object exposed with Boost....
Definition id.hpp:18