eigenpy 3.12.0
Bindings between Numpy and Eigen using Boost.Python
Loading...
Searching...
No Matches
SVDBase.hpp
1/*
2 * Copyright 2025 INRIA
3 */
4
5#ifndef __eigenpy_decompositions_svdbase_hpp__
6#define __eigenpy_decompositions_svdbase_hpp__
7
8#include <Eigen/SVD>
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 Derived>
19 : public boost::python::def_visitor<SVDBaseVisitor<Derived>> {
20 typedef Derived Solver;
21
22 typedef typename Derived::MatrixType MatrixType;
23 typedef typename MatrixType::Scalar Scalar;
24 typedef typename MatrixType::RealScalar RealScalar;
25
26 typedef Eigen::Matrix<Scalar, Eigen::Dynamic, 1, MatrixType::Options>
27 VectorXs;
28 typedef Eigen::Matrix<Scalar, Eigen::Dynamic, Eigen::Dynamic,
29 MatrixType::Options>
30 MatrixXs;
31
32 template <class PyClass>
33 void visit(PyClass &cl) const {
34 cl.def(bp::init<>(bp::arg("self"), "Default constructor"))
35
36 .def("computeU", &Solver::computeU, bp::arg("self"),
37 "Returns true if U (full or thin) is asked for in this "
38 "SVD decomposition ")
39 .def("computeV", &Solver::computeV, bp::arg("self"),
40 "Returns true if V (full or thin) is asked for in this "
41 "SVD decomposition ")
42
43 .def("info", &Solver::info, bp::arg("self"),
44 "Reports whether previous computation was successful. ")
45
46 .def("matrixU", &matrixU, bp::arg("self"), "Returns the matrix U.")
47 .def("matrixV", &matrixV, bp::arg("self"), "Returns the matrix V.")
48
49 .def("nonzeroSingularValues", &Solver::nonzeroSingularValues,
50 bp::arg("self"),
51 "Returns the number of singular values that are not exactly 0 ")
52 .def("rank", &Solver::rank, bp::arg("self"),
53 "the rank of the matrix of which *this is the SVD. ")
54
55 .def("setThreshold",
56 (Solver & (Solver::*)(const RealScalar &)) & Solver::setThreshold,
57 bp::args("self", "threshold"),
58 "Allows to prescribe a threshold to be used by certain methods, "
59 "such as "
60 "rank() and solve(), which need to determine when singular values "
61 "are "
62 "to be considered nonzero. This is not used for the SVD "
63 "decomposition "
64 "itself.\n"
65 "\n"
66 "When it needs to get the threshold value, Eigen calls "
67 "threshold(). "
68 "The default is NumTraits<Scalar>::epsilon()",
69 bp::return_self<>())
70 .def(
71 "setThreshold",
72 +[](Solver &self) -> Solver & {
73 return self.setThreshold(Eigen::Default);
74 },
75 bp::arg("self"),
76 "Allows to come back to the default behavior, letting Eigen use "
77 "its default formula for determining the threshold.",
78 bp::return_self<>())
79
80 .def("singularValues", &Solver::singularValues, bp::arg("self"),
81 "Returns the vector of singular values.",
82 bp::return_value_policy<bp::copy_const_reference>())
83
84 .def("solve", &solve<VectorXs>, bp::args("self", "b"),
85 "Returns the solution x of A x = b using the current "
86 "decomposition of A.")
87 .def("solve", &solve<MatrixXs>, bp::args("self", "B"),
88 "Returns the solution X of A X = B using the current "
89 "decomposition of A where B is a right hand side matrix.")
90
91 .def("threshold", &Solver::threshold, bp::arg("self"),
92 "Returns the threshold that will be used by certain methods such "
93 "as rank(). ");
94 }
95
96 private:
97 static MatrixXs matrixU(const Solver &self) { return self.matrixU(); }
98 static MatrixXs matrixV(const Solver &self) { return self.matrixV(); }
99
100 template <typename MatrixOrVector>
101 static MatrixOrVector solve(const Solver &self, const MatrixOrVector &vec) {
102 return self.solve(vec);
103 }
104};
105
106} // namespace eigenpy
107
108#endif // ifndef __eigenpy_decompositions_svdbase_hpp__