eigenpy 3.12.0
Bindings between Numpy and Eigen using Boost.Python
Loading...
Searching...
No Matches
RealQZ.hpp
1/*
2 * Copyright 2025 INRIA
3 */
4
5#ifndef __eigenpy_decompositions_generalized_real_qz_hpp__
6#define __eigenpy_decompositions_generalized_real_qz_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<RealQZVisitor<_MatrixType>> {
20 typedef _MatrixType MatrixType;
21 typedef typename MatrixType::Scalar Scalar;
22 typedef Eigen::RealQZ<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, MatrixType, bp::optional<bool>>(
29 bp::args("A", "B", "computeQZ"),
30 "Constructor; computes real QZ decomposition of given matrices. "))
31
32 .def("compute", &RealQZVisitor::compute_proxy<MatrixType>,
33 bp::args("self", "A", "B"),
34 "Computes QZ decomposition of given matrix. ",
35 bp::return_self<>())
36 .def("compute",
37 (Solver &
38 (Solver::*)(const MatrixType& A, const MatrixType& B, bool)) &
39 Solver::compute,
40 bp::args("self", "A", "B", "computeEigenvectors"),
41 "Computes QZ decomposition of given matrix. ", bp::return_self<>())
42
43 .def("info", &Solver::info, bp::arg("self"),
44 "NumericalIssue if the input contains INF or NaN values or "
45 "overflow occured. Returns Success otherwise.")
46
47 .def("matrixQ", &Solver::matrixQ, bp::arg("self"),
48 "Returns matrix Q in the QZ decomposition. ",
49 bp::return_value_policy<bp::copy_const_reference>())
50 .def("matrixS", &Solver::matrixS, bp::arg("self"),
51 "Returns matrix S in the QZ decomposition. ",
52 bp::return_value_policy<bp::copy_const_reference>())
53 .def("matrixT", &Solver::matrixT, bp::arg("self"),
54 "Returns matrix T in the QZ decomposition. ",
55 bp::return_value_policy<bp::copy_const_reference>())
56 .def("matrixZ", &Solver::matrixZ, bp::arg("self"),
57 "Returns matrix Z in the QZ decomposition. ",
58 bp::return_value_policy<bp::copy_const_reference>())
59
60 .def("iterations", &Solver::iterations, bp::arg("self"),
61 "Returns number of performed QR-like iterations. ")
62 .def("setMaxIterations", &Solver::setMaxIterations,
63 bp::args("self", "max_iter"),
64 "Sets the maximum number of iterations allowed.",
65 bp::return_self<>());
66 }
67
68 static void expose() {
69 static const std::string classname =
70 "RealQZVisitor" + scalar_name<Scalar>::shortname();
71 expose(classname);
72 }
73
74 static void expose(const std::string& name) {
75 bp::class_<Solver>(name.c_str(), bp::no_init)
76 .def(RealQZVisitor())
77 .def(IdVisitor<Solver>());
78 }
79
80 private:
81 template <typename MatrixType>
82 static Solver& compute_proxy(Solver& self, const MatrixType& A,
83 const MatrixType& B) {
84 return self.compute(A, B);
85 }
86};
87
88} // namespace eigenpy
89
90#endif // ifndef __eigenpy_decompositions_generalized_real_qz_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