eigenpy 3.12.0
Bindings between Numpy and Eigen using Boost.Python
Loading...
Searching...
No Matches
RealSchur.hpp
1/*
2 * Copyright 2025 INRIA
3 */
4
5#ifndef __eigenpy_decompositions_generalized_real_schur_hpp__
6#define __eigenpy_decompositions_generalized_real_schur_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<RealSchurVisitor<_MatrixType>> {
20 typedef _MatrixType MatrixType;
21 typedef typename MatrixType::Scalar Scalar;
22 typedef Eigen::RealSchur<MatrixType> Solver;
23
24 template <class PyClass>
25 void visit(PyClass& cl) const {
26 cl.def(
27 bp::init<Eigen::DenseIndex>(bp::arg("size"), "Default constructor. "))
28 .def(bp::init<MatrixType, bp::optional<bool>>(
29 bp::args("matrix", "computeU"),
30 "Constructor; computes real Schur decomposition of given matrix. "))
31
32 .def("compute", &RealSchurVisitor::compute_proxy<MatrixType>,
33 bp::args("self", "matrix"),
34 "Computes Schur decomposition of given matrix. ",
35 bp::return_self<>())
36 .def("compute",
37 (Solver &
38 (Solver::*)(const Eigen::EigenBase<MatrixType>& matrix, bool)) &
39 Solver::compute,
40 bp::args("self", "matrix", "computeEigenvectors"),
41 "Computes Schur decomposition of given matrix. ",
42 bp::return_self<>())
43
44 .def("computeFromHessenberg",
45 (Solver & (Solver::*)(const MatrixType& matrixH,
46 const MatrixType& matrixQ, bool)) &
47 Solver::computeFromHessenberg,
48 bp::args("self", "matrixH", "matrixQ", "computeU"),
49 "Compute Schur decomposition from a given Hessenberg matrix. ",
50 bp::return_self<>())
51
52 .def("info", &Solver::info, bp::arg("self"),
53 "NumericalIssue if the input contains INF or NaN values or "
54 "overflow occured. Returns Success otherwise.")
55
56 .def("matrixT", &Solver::matrixT, bp::arg("self"),
57 "Returns the quasi-triangular matrix in the Schur decomposition.",
58 bp::return_value_policy<bp::copy_const_reference>())
59 .def("matrixU", &Solver::matrixU, bp::arg("self"),
60 "Returns the orthogonal matrix in the Schur decomposition. ",
61 bp::return_value_policy<bp::copy_const_reference>())
62
63 .def("setMaxIterations", &Solver::setMaxIterations,
64 bp::args("self", "max_iter"),
65 "Sets the maximum number of iterations allowed.",
66 bp::return_self<>())
67 .def("getMaxIterations", &Solver::getMaxIterations, bp::arg("self"),
68 "Returns the maximum number of iterations.");
69 }
70
71 static void expose() {
72 static const std::string classname =
73 "RealSchurVisitor" + scalar_name<Scalar>::shortname();
74 expose(classname);
75 }
76
77 static void expose(const std::string& name) {
78 bp::class_<Solver>(name.c_str(), bp::no_init)
79 .def(RealSchurVisitor())
80 .def(IdVisitor<Solver>());
81 }
82
83 private:
84 template <typename MatrixType>
85 static Solver& compute_proxy(Solver& self, const MatrixType& A) {
86 return self.compute(A);
87 }
88};
89
90} // namespace eigenpy
91
92#endif // ifndef __eigenpy_decompositions_generalized_real_schur_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