eigenpy 3.12.0
Bindings between Numpy and Eigen using Boost.Python
Loading...
Searching...
No Matches
Tridiagonalization.hpp
1/*
2 * Copyright 2025 INRIA
3 */
4
5#ifndef __eigenpy_decompositions_tridiagonalization_hpp__
6#define __eigenpy_decompositions_tridiagonalization_hpp__
7
8#include <Eigen/Core>
9#include <Eigen/Eigenvalues>
10
11#include "eigenpy/eigen-to-python.hpp"
12#include "eigenpy/eigenpy.hpp"
13#include "eigenpy/utils/scalar-name.hpp"
14
15namespace eigenpy {
16
17template <typename _MatrixType>
18struct TridiagonalizationVisitor : public boost::python::def_visitor<
19 TridiagonalizationVisitor<_MatrixType>> {
20 typedef _MatrixType MatrixType;
21 typedef typename MatrixType::Scalar Scalar;
22 typedef Eigen::Tridiagonalization<MatrixType> Solver;
23 typedef Eigen::VectorXd VectorType;
24
25 template <class PyClass>
26 void visit(PyClass &cl) const {
27 cl.def(
28 bp::init<Eigen::DenseIndex>(bp::arg("size"), "Default constructor. "))
29 .def(bp::init<MatrixType>(bp::arg("matrix"),
30 "Constructor; computes tridiagonal "
31 "decomposition of given matrix. "))
32
33 .def(
34 "compute",
35 (Solver & (Solver::*)(const Eigen::EigenBase<MatrixType> &matrix)) &
36 Solver::compute,
37 bp::args("self", "matrix"),
38 "Computes tridiagonal decomposition of given matrix. ",
39 bp::return_self<>())
40
41 .def("householderCoefficients", &Solver::householderCoefficients,
42 bp::arg("self"), "Returns the Householder coefficients. ")
43 .def("packedMatrix", &Solver::packedMatrix, bp::arg("self"),
44 "Returns the internal representation of the decomposition. ",
45 bp::return_value_policy<bp::copy_const_reference>())
46
47 .def(
48 "matrixQ",
49 +[](const Solver &c) -> MatrixType { return c.matrixQ(); },
50 "Returns the unitary matrix Q in the decomposition.")
51 .def(
52 "matrixT",
53 +[](const Solver &c) -> MatrixType { return c.matrixT(); },
54 "Returns an expression of the tridiagonal matrix T in the "
55 "decomposition.")
56
57 .def(
58 "diagonal",
59 +[](const Solver &c) -> VectorType { return c.diagonal(); },
60 bp::arg("self"),
61 "Returns the diagonal of the tridiagonal matrix T in the "
62 "decomposition. ")
63
64 .def(
65 "subDiagonal",
66 +[](const Solver &c) -> VectorType { return c.subDiagonal(); },
67 bp::arg("self"),
68 "Returns the subdiagonal of the tridiagonal matrix T in the "
69 "decomposition.");
70 }
71
72 static void expose() {
73 static const std::string classname =
74 "TridiagonalizationVisitor" + scalar_name<Scalar>::shortname();
75 expose(classname);
76 }
77
78 static void expose(const std::string &name) {
79 bp::class_<Solver>(name.c_str(), bp::no_init)
81 .def(IdVisitor<Solver>());
82 }
83};
84
85} // namespace eigenpy
86
87#endif // ifndef __eigenpy_decompositions_tridiagonalization_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