19 :
public boost::python::def_visitor<FullPivLUSolverVisitor<_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,
28 typedef Eigen::FullPivLU<MatrixType> Solver;
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, Eigen::DenseIndex>(
34 bp::args(
"self",
"rows",
"cols"),
35 "Default constructor with memory preallocation"))
36 .def(bp::init<MatrixType>(
37 bp::args(
"self",
"matrix"),
38 "Constructs a LU factorization from a given matrix."))
44 (Solver & (Solver::*)(
const Eigen::EigenBase<MatrixType> &matrix)) &
46 bp::args(
"self",
"matrix"),
47 "Computes the LU decomposition of the given matrix.",
50 .def(
"determinant", &Solver::determinant, bp::arg(
"self"),
51 "Returns the determinant of the matrix of which *this is the LU " 53 .def(
"dimensionOfKernel", &Solver::dimensionOfKernel, bp::arg(
"self"),
54 "Returns the dimension of the kernel of the matrix of which *this " 55 "is the LU decomposition.")
58 +[](Solver &self,
const MatrixType &mat) -> MatrixType {
59 return self.image(mat);
61 bp::args(
"self",
"originalMatrix"),
62 "Returns the image of the matrix, also called its column-space. " 63 "The columns of the returned matrix will form a basis of the " 64 "image (column-space).")
67 +[](Solver &self) -> MatrixType {
return self.inverse(); },
69 "Returns the inverse of the matrix of which *this is the LU " 72 .def(
"isInjective", &Solver::isInjective, bp::arg(
"self"))
73 .def(
"isInvertible", &Solver::isInvertible, bp::arg(
"self"))
74 .def(
"isSurjective", &Solver::isSurjective, bp::arg(
"self"))
77 "kernel", +[](Solver &self) -> MatrixType {
return self.kernel(); },
79 "Returns the kernel of the matrix, also called its null-space. " 80 "The columns of the returned matrix will form a basis of the " 83 .def(
"matrixLU", &Solver::matrixLU, bp::arg(
"self"),
84 "Returns the LU decomposition matrix.",
85 bp::return_internal_reference<>())
87 .def(
"maxPivot", &Solver::maxPivot, bp::arg(
"self"))
88 .def(
"nonzeroPivots", &Solver::nonzeroPivots, bp::arg(
"self"))
90 .def(
"permutationP", &Solver::permutationP, bp::arg(
"self"),
91 "Returns the permutation P.",
92 bp::return_value_policy<bp::copy_const_reference>())
93 .def(
"permutationQ", &Solver::permutationQ, bp::arg(
"self"),
94 "Returns the permutation Q.",
95 bp::return_value_policy<bp::copy_const_reference>())
97 .def(
"rank", &Solver::rank, bp::arg(
"self"))
99 .def(
"rcond", &Solver::rcond, bp::arg(
"self"),
100 "Returns an estimate of the reciprocal condition number of the " 102 .def(
"reconstructedMatrix", &Solver::reconstructedMatrix,
104 "Returns the matrix represented by the decomposition, i.e., it " 105 "returns the product: P-1LUQ-1. This function is provided for " 110 (Solver & (Solver::*)(
const RealScalar &)) & Solver::setThreshold,
111 bp::args(
"self",
"threshold"),
112 "Allows to prescribe a threshold to be used by certain methods, " 113 "such as rank(), who need to determine when pivots are to be " 114 "considered nonzero. This is not used for the LU decomposition " 117 "When it needs to get the threshold value, Eigen calls " 118 "threshold(). By default, this uses a formula to automatically " 119 "determine a reasonable threshold. Once you have called the " 120 "present method setThreshold(const RealScalar&), your value is " 123 "Note: A pivot will be considered nonzero if its absolute value " 124 "is strictly greater than |pivot| ⩽ threshold×|maxpivot| where " 125 "maxpivot is the biggest pivot.",
129 +[](Solver &self) -> Solver & {
130 return self.setThreshold(Eigen::Default);
133 "Allows to come back to the default behavior, letting Eigen use " 134 "its default formula for determining the threshold.",
137 .def(
"solve", &solve<VectorXs>, bp::args(
"self",
"b"),
138 "Returns the solution x of A x = b using the current " 139 "decomposition of A.")
140 .def(
"solve", &solve<MatrixXs>, bp::args(
"self",
"B"),
141 "Returns the solution X of A X = B using the current " 142 "decomposition of A where B is a right hand side matrix.")
144 .def(
"threshold", &Solver::threshold, bp::arg(
"self"),
145 "Returns the threshold that will be used by certain methods such " 150 static const std::string classname =
151 "FullPivLU" + scalar_name<Scalar>::shortname();
155 static void expose(
const std::string &name) {
158 "LU decomposition of a matrix with complete pivoting, and related " 160 "This class represents a LU decomposition of any matrix, with complete " 161 "pivoting: the matrix A is decomposed as A=P−1LUQ−1 where L is " 162 "unit-lower-triangular, U is upper-triangular, and P and Q are " 163 "permutation matrices. This is a rank-revealing LU decomposition. " 164 "The eigenvalues (diagonal coefficients) of U are sorted in such a " 165 "way that any zeros are at the end.\n\n" 166 "This decomposition provides the generic approach to solving systems " 168 "linear equations, computing the rank, invertibility, inverse, kernel, " 169 "and determinant. \n\n" 170 "This LU decomposition is very stable and well tested with large " 172 "However there are use cases where the SVD decomposition is inherently " 174 "stable and/or flexible. For example, when computing the kernel of a " 176 "working with the SVD allows to select the smallest singular values of " 177 "the matrix, something that the LU decomposition doesn't see.",
184 template <
typename MatrixOrVector>
185 static MatrixOrVector solve(
const Solver &self,
const MatrixOrVector &vec) {
186 return self.solve(vec);