eigenpy 3.12.0
Bindings between Numpy and Eigen using Boost.Python
Loading...
Searching...
No Matches
HessenbergDecomposition.hpp
1/*
2 * Copyright 2025 INRIA
3 */
4
5#ifndef __eigenpy_decompositions_hessenberg_decomposition_hpp__
6#define __eigenpy_decompositions_hessenberg_decomposition_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>
19 : public boost::python::def_visitor<
20 HessenbergDecompositionVisitor<_MatrixType>> {
21 typedef _MatrixType MatrixType;
22 typedef typename MatrixType::Scalar Scalar;
23 typedef Eigen::HessenbergDecomposition<MatrixType> Solver;
24
25 template <class PyClass>
26 void visit(PyClass& cl) const {
27 cl
28 .def(bp::init<Eigen::DenseIndex>(
29 bp::arg("size"),
30 "Default constructor; the decomposition will be computed later. "))
31 .def(bp::init<MatrixType>(
32 bp::arg("matrix"),
33 "Constructor; computes Hessenberg decomposition of given matrix. "))
34
35 .def(
36 "compute",
37 (Solver & (Solver::*)(const Eigen::EigenBase<MatrixType>& matrix)) &
38 Solver::compute,
39 bp::args("self", "A"),
40 "Computes Hessenberg decomposition of given matrix. ",
41 bp::return_self<>())
42
43 .def("householderCoefficients", &Solver::householderCoefficients,
44 bp::arg("self"), "Returns the Householder coefficients. ",
45 bp::return_value_policy<bp::copy_const_reference>())
46
47 .def(
48 "matrixQ",
49 +[](const Solver& c) -> MatrixType { return c.matrixQ(); },
50 "Reconstructs the orthogonal matrix Q in the decomposition.")
51 .def(
52 "matrixH",
53 +[](const Solver& c) -> MatrixType { return c.matrixH(); },
54 "Constructs the Hessenberg matrix H in the decomposition.")
55
56 .def("packedMatrix", &Solver::packedMatrix, bp::arg("self"),
57 "Returns the internal representation of the decomposition. ",
58 bp::return_value_policy<bp::copy_const_reference>());
59 }
60
61 static void expose() {
62 static const std::string classname =
63 "HessenbergDecomposition" + scalar_name<Scalar>::shortname();
64 expose(classname);
65 }
66
67 static void expose(const std::string& name) {
68 bp::class_<Solver>(name.c_str(), bp::no_init)
70 .def(IdVisitor<Solver>());
71 }
72};
73
74} // namespace eigenpy
75
76#endif // ifndef __eigenpy_decompositions_hessenberg_decomposition_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