eigenpy 3.12.0
Bindings between Numpy and Eigen using Boost.Python
Loading...
Searching...
No Matches
JacobiSVD.hpp
1/*
2 * Copyright 2025 INRIA
3 */
4
5#ifndef __eigenpy_decompositions_jacobisvd_hpp__
6#define __eigenpy_decompositions_jacobisvd_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#include "eigenpy/decompositions/SVDBase.hpp"
15
16namespace eigenpy {
17
18template <typename JacobiSVD>
20 : public boost::python::def_visitor<JacobiSVDVisitor<JacobiSVD>> {
21 typedef typename JacobiSVD::MatrixType MatrixType;
22 typedef typename MatrixType::Scalar Scalar;
23
24 template <class PyClass>
25 void visit(PyClass &cl) const {
26 cl.def(bp::init<>(bp::arg("self"), "Default constructor"))
27 .def(bp::init<Eigen::DenseIndex, Eigen::DenseIndex,
28 bp::optional<unsigned int>>(
29 bp::args("self", "rows", "cols", "computationOptions "),
30 "Default Constructor with memory preallocation."))
31 .def(bp::init<MatrixType, bp::optional<unsigned int>>(
32 bp::args("self", "matrix", "computationOptions "),
33 "Constructor performing the decomposition of given matrix."))
34
35 .def("cols", &JacobiSVD::cols, bp::arg("self"),
36 "Returns the number of columns. ")
37 .def("compute",
38 (JacobiSVD & (JacobiSVD::*)(const MatrixType &matrix)) &
39 JacobiSVD::compute,
40 bp::args("self", "matrix"),
41 "Method performing the decomposition of given matrix. Computes "
42 "Thin/Full "
43 "unitaries U/V if specified using the Options template parameter "
44 "or the class constructor. ",
45 bp::return_self<>())
46 .def("compute",
47 (JacobiSVD & (JacobiSVD::*)(const MatrixType &matrix,
48 unsigned int computationOptions)) &
49 JacobiSVD::compute,
50 bp::args("self", "matrix", "computationOptions"),
51 "Method performing the decomposition of given matrix, as "
52 "specified by the computationOptions parameter. ",
53 bp::return_self<>())
54 .def("rows", &JacobiSVD::rows, bp::arg("self"),
55 "Returns the number of rows.")
56
58 }
59
60 static void expose() {
61 static const std::string classname =
62 "JacobiSVD_" + scalar_name<Scalar>::shortname();
63 expose(classname);
64 }
65
66 static void expose(const std::string &name) {
67 bp::class_<JacobiSVD, boost::noncopyable>(
68 name.c_str(),
69 "Two-sided Jacobi SVD decomposition of a rectangular matrix. \n\n"
70 "SVD decomposition consists in decomposing any n-by-p matrix A as a "
71 "product "
72 "A=USV∗ where U is a n-by-n unitary, V is a p-by-p unitary, and S is a "
73 "n-by-p r"
74 "eal positive matrix which is zero outside of its main diagonal; the "
75 "diagonal "
76 "entries of S are known as the singular values of A and the columns of "
77 "U and V "
78 "are known as the left and right singular vectors of A respectively. "
79 "\n\n"
80 "Singular values are always sorted in decreasing order. \n\n "
81 "This JacobiSVD decomposition computes only the singular values by "
82 "default. "
83 "If you want U or V, you need to ask for them explicitly. \n\n"
84 "You can ask for only thin U or V to be computed, meaning the "
85 "following. "
86 "In case of a rectangular n-by-p matrix, letting m be the smaller "
87 "value among "
88 "n and p, there are only m singular vectors; the remaining columns of "
89 "U and V "
90 "do not correspond to actual singular vectors. Asking for thin U or V "
91 "means asking "
92 "for only their m first columns to be formed. So U is then a n-by-m "
93 "matrix, and V "
94 "is then a p-by-m matrix. Notice that thin U and V are all you need "
95 "for (least squares) "
96 "solving.",
97 bp::no_init)
98 .def(JacobiSVDVisitor())
100 }
101};
102
103} // namespace eigenpy
104
105#endif // ifndef __eigenpy_decompositions_jacobisvd_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