5 #ifndef __eigenpy_decompositions_houselholder_qr_hpp__
6 #define __eigenpy_decompositions_houselholder_qr_hpp__
15 template <
typename _MatrixType>
17 :
public boost::python::def_visitor<
18 HouseholderQRSolverVisitor<_MatrixType> > {
20 typedef typename MatrixType::Scalar
Scalar;
22 typedef Eigen::Matrix<Scalar, Eigen::Dynamic, 1, MatrixType::Options>
24 typedef Eigen::Matrix<
Scalar, Eigen::Dynamic, Eigen::Dynamic,
27 typedef Eigen::HouseholderQR<MatrixType>
Solver;
30 template <
class PyClass>
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()."))
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.")
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>())
72 (
Solver & (
Solver::*)(
const Eigen::EigenBase<MatrixType> &matrix)) &
74 bp::args(
"self",
"matrix"),
75 "Computes the QR factorization of given matrix.",
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.");
84 static const std::string classname =
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"
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"
101 "This Householder QR decomposition is faster, but less numerically "
102 "stable and less feature-full than FullPivHouseholderQR or "
103 "ColPivHouseholderQR.",
110 template <
typename MatrixOrVector>
112 return self.solve(
vec);
118 #endif // ifndef __eigenpy_decompositions_houselholder_qr_hpp__