eigenpy 3.12.0
Bindings between Numpy and Eigen using Boost.Python
Loading...
Searching...
No Matches
FullPivLU.hpp
1/*
2 * Copyright 2025 INRIA
3 */
4
5#ifndef __eigenpy_decompositions_fullpivlu_hpp__
6#define __eigenpy_decompositions_fullpivlu_hpp__
7
8#include <Eigen/LU>
9#include <Eigen/Core>
10
11#include "eigenpy/eigenpy.hpp"
12#include "eigenpy/utils/scalar-name.hpp"
13#include "eigenpy/eigen/EigenBase.hpp"
14
15namespace eigenpy {
16
17template <typename _MatrixType>
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>
24 VectorXs;
25 typedef Eigen::Matrix<Scalar, Eigen::Dynamic, Eigen::Dynamic,
26 MatrixType::Options>
27 MatrixXs;
28 typedef Eigen::FullPivLU<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, 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."))
39
41
42 .def(
43 "compute",
44 (Solver & (Solver::*)(const Eigen::EigenBase<MatrixType> &matrix)) &
45 Solver::compute,
46 bp::args("self", "matrix"),
47 "Computes the LU decomposition of the given matrix.",
48 bp::return_self<>())
49
50 .def("determinant", &Solver::determinant, bp::arg("self"),
51 "Returns the determinant of the matrix of which *this is the LU "
52 "decomposition.")
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.")
56 .def(
57 "image",
58 +[](Solver &self, const MatrixType &mat) -> MatrixType {
59 return self.image(mat);
60 },
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).")
65 .def(
66 "inverse",
67 +[](Solver &self) -> MatrixType { return self.inverse(); },
68 bp::arg("self"),
69 "Returns the inverse of the matrix of which *this is the LU "
70 "decomposition.")
71
72 .def("isInjective", &Solver::isInjective, bp::arg("self"))
73 .def("isInvertible", &Solver::isInvertible, bp::arg("self"))
74 .def("isSurjective", &Solver::isSurjective, bp::arg("self"))
75
76 .def(
77 "kernel", +[](Solver &self) -> MatrixType { return self.kernel(); },
78 bp::arg("self"),
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 "
81 "kernel.")
82
83 .def("matrixLU", &Solver::matrixLU, bp::arg("self"),
84 "Returns the LU decomposition matrix.",
85 bp::return_internal_reference<>())
86
87 .def("maxPivot", &Solver::maxPivot, bp::arg("self"))
88 .def("nonzeroPivots", &Solver::nonzeroPivots, bp::arg("self"))
89
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>())
96
97 .def("rank", &Solver::rank, bp::arg("self"))
98
99 .def("rcond", &Solver::rcond, bp::arg("self"),
100 "Returns an estimate of the reciprocal condition number of the "
101 "matrix.")
102 .def("reconstructedMatrix", &Solver::reconstructedMatrix,
103 bp::arg("self"),
104 "Returns the matrix represented by the decomposition, i.e., it "
105 "returns the product: P-1LUQ-1. This function is provided for "
106 "debug "
107 "purpose.")
108
109 .def("setThreshold",
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 "
115 "itself.\n"
116 "\n"
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 "
121 "used instead.\n"
122 "\n"
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.",
126 bp::return_self<>())
127 .def(
128 "setThreshold",
129 +[](Solver &self) -> Solver & {
130 return self.setThreshold(Eigen::Default);
131 },
132 bp::arg("self"),
133 "Allows to come back to the default behavior, letting Eigen use "
134 "its default formula for determining the threshold.",
135 bp::return_self<>())
136
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.")
143
144 .def("threshold", &Solver::threshold, bp::arg("self"),
145 "Returns the threshold that will be used by certain methods such "
146 "as rank().");
147 }
148
149 static void expose() {
150 static const std::string classname =
151 "FullPivLU" + scalar_name<Scalar>::shortname();
152 expose(classname);
153 }
154
155 static void expose(const std::string &name) {
156 bp::class_<Solver>(
157 name.c_str(),
158 "LU decomposition of a matrix with complete pivoting, and related "
159 "features.\n\n"
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 "
167 "of "
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 "
171 "matrices. "
172 "However there are use cases where the SVD decomposition is inherently "
173 "more "
174 "stable and/or flexible. For example, when computing the kernel of a "
175 "matrix, "
176 "working with the SVD allows to select the smallest singular values of "
177 "the matrix, something that the LU decomposition doesn't see.",
178 bp::no_init)
179 .def(IdVisitor<Solver>())
181 }
182
183 private:
184 template <typename MatrixOrVector>
185 static MatrixOrVector solve(const Solver &self, const MatrixOrVector &vec) {
186 return self.solve(vec);
187 }
188};
189
190} // namespace eigenpy
191
192#endif // ifndef __eigenpy_decompositions_fullpivlu_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