eigenpy 3.12.0
Bindings between Numpy and Eigen using Boost.Python
Loading...
Searching...
No Matches
BDCSVD.hpp
1/*
2 * Copyright 2025 INRIA
3 */
4
5#ifndef __eigenpy_decompositions_bdcsvd_hpp__
6#define __eigenpy_decompositions_bdcsvd_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 _MatrixType>
20 : public boost::python::def_visitor<BDCSVDVisitor<_MatrixType>> {
21 typedef _MatrixType MatrixType;
22 typedef Eigen::BDCSVD<MatrixType> Solver;
23 typedef typename MatrixType::Scalar Scalar;
24
25 template <class PyClass>
26 void visit(PyClass &cl) const {
27 cl.def(bp::init<>(bp::arg("self"), "Default constructor"))
28 .def(bp::init<Eigen::DenseIndex, Eigen::DenseIndex,
29 bp::optional<unsigned int>>(
30 bp::args("self", "rows", "cols", "computationOptions "),
31 "Default Constructor with memory preallocation. "))
32 .def(bp::init<MatrixType, bp::optional<unsigned int>>(
33 bp::args("self", "matrix", "computationOptions "),
34 "Constructor performing the decomposition of given matrix."))
35
36 .def("cols", &Solver::cols, bp::arg("self"),
37 "Returns the number of columns. ")
38 .def("compute",
39 (Solver & (Solver::*)(const MatrixType &matrix)) & Solver::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 (Solver & (Solver::*)(const MatrixType &matrix,
48 unsigned int computationOptions)) &
49 Solver::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", &Solver::rows, bp::arg("self"),
55 "Returns the number of rows. ")
56 .def("setSwitchSize", &Solver::setSwitchSize, bp::args("self", "s"))
57
59 }
60
61 static void expose() {
62 static const std::string classname =
63 "BDCSVD_" + scalar_name<Scalar>::shortname();
64 expose(classname);
65 }
66
67 static void expose(const std::string &name) {
68 bp::class_<Solver, boost::noncopyable>(
69 name.c_str(),
70 "Class Bidiagonal Divide and Conquer SVD.\n\n"
71 "This class first reduces the input matrix to bi-diagonal form using "
72 "class "
73 "UpperBidiagonalization, and then performs a divide-and-conquer "
74 "diagonalization. "
75 "Small blocks are diagonalized using class JacobiSVD. You can control "
76 "the "
77 "switching size with the setSwitchSize() method, default is 16. For "
78 "small matrice "
79 "(<16), it is thus preferable to directly use JacobiSVD. For larger "
80 "ones, BDCSVD "
81 "is highly recommended and can several order of magnitude faster.\n\n"
82 "Warming: this algorithm is unlikely to provide accurate result when "
83 "compiled with "
84 "unsafe math optimizations. For instance, this concerns Intel's "
85 "compiler (ICC), which "
86 "performs such optimization by default unless you compile with the "
87 "-fp-model precise "
88 "option. Likewise, the -ffast-math option of GCC or clang will "
89 "significantly degrade the "
90 "accuracy.",
91 bp::no_init)
92 .def(BDCSVDVisitor())
93 .def(IdVisitor<Solver>());
94 }
95};
96
97} // namespace eigenpy
98
99#endif // ifndef __eigenpy_decompositions_bdcsvd_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