Skip to content

Commit a4b2eb4

Browse files
committed
eigenvalues: RealQZ
1 parent fe5f3ce commit a4b2eb4

File tree

7 files changed

+147
-13
lines changed

7 files changed

+147
-13
lines changed

CMakeLists.txt

Lines changed: 2 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -239,6 +239,7 @@ set(${PROJECT_NAME}_DECOMPOSITIONS_HEADERS
239239
include/eigenpy/decompositions/GeneralizedEigenSolver.hpp
240240
include/eigenpy/decompositions/GeneralizedSelfAdjointEigenSolver.hpp
241241
include/eigenpy/decompositions/HessenbergDecomposition.hpp
242+
include/eigenpy/decompositions/RealQZ.hpp
242243
include/eigenpy/decompositions/ComplexEigenSolver.hpp
243244
include/eigenpy/decompositions/ComplexSchur.hpp
244245
include/eigenpy/decompositions/PermutationMatrix.hpp
@@ -333,6 +334,7 @@ set(${PROJECT_NAME}_DECOMPOSITIONS_SOURCES
333334
src/decompositions/jacobisvd-solver.cpp
334335
src/decompositions/fullpivlu-solver.cpp
335336
src/decompositions/hessenberg-decomposition.cpp
337+
src/decompositions/real-qz.cpp
336338
src/decompositions/partialpivlu-solver.cpp
337339
src/decompositions/minres-solver.cpp
338340
src/decompositions/sparse-lu-solver.cpp

include/eigenpy/decompositions/HessenbergDecomposition.hpp

Lines changed: 17 additions & 13 deletions
Original file line numberDiff line numberDiff line change
@@ -24,26 +24,29 @@ struct HessenbergDecompositionVisitor
2424

2525
template <class PyClass>
2626
void visit(PyClass& cl) const {
27-
cl.def(bp::init<Eigen::DenseIndex>(
28-
bp::arg("size"), "Default constructor; the decomposition will be computed later. "))
29-
.def(bp::init<MatrixType>(bp::arg("matrix"),
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"),
3033
"Constructor; computes Hessenberg decomposition of given matrix. "))
3134

3235
.def("compute",
3336
&HessenbergDecompositionVisitor::compute_proxy<MatrixType>,
3437
bp::args("self", "A"),
3538
"Computes Hessenberg decomposition of given matrix. ",
3639
bp::return_self<>())
37-
.def("compute",
38-
(Solver &
39-
(Solver::*)(const Eigen::EigenBase<MatrixType>& matrix)) &
40-
Solver::compute,
41-
bp::args("self", "A"),
42-
"Computes Hessenberg decomposition of given matrix. ",
43-
bp::return_self<>())
40+
.def(
41+
"compute",
42+
(Solver & (Solver::*)(const Eigen::EigenBase<MatrixType>& matrix)) &
43+
Solver::compute,
44+
bp::args("self", "A"),
45+
"Computes Hessenberg decomposition of given matrix. ",
46+
bp::return_self<>())
4447

45-
.def("householderCoefficients", &Solver::householderCoefficients, bp::arg("self"),
46-
"Returns the Householder coefficients. ",
48+
.def("householderCoefficients", &Solver::householderCoefficients,
49+
bp::arg("self"), "Returns the Householder coefficients. ",
4750
bp::return_value_policy<bp::copy_const_reference>())
4851

4952
// matrixH
@@ -68,7 +71,8 @@ struct HessenbergDecompositionVisitor
6871

6972
private:
7073
template <typename MatrixType>
71-
static Solver& compute_proxy(Solver& self, const Eigen::EigenBase<MatrixType>& matrix) {
74+
static Solver& compute_proxy(Solver& self,
75+
const Eigen::EigenBase<MatrixType>& matrix) {
7276
return self.compute(matrix);
7377
}
7478
};
Lines changed: 93 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,93 @@
1+
/*
2+
* Copyright 2020 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+
15+
namespace eigenpy {
16+
17+
template <typename _MatrixType>
18+
struct RealQZVisitor
19+
: public boost::python::def_visitor<
20+
RealQZVisitor<_MatrixType>> {
21+
typedef _MatrixType MatrixType;
22+
typedef typename MatrixType::Scalar Scalar;
23+
typedef Eigen::RealQZ<MatrixType> Solver;
24+
25+
template <class PyClass>
26+
void visit(PyClass& cl) const {
27+
cl.def(bp::init<Eigen::DenseIndex>(
28+
bp::arg("size"), "Default constructor. "))
29+
.def(bp::init<MatrixType, MatrixType, bp::optional<bool>>(
30+
bp::args("A", "B", "computeQZ"),
31+
"Constructor; computes real QZ decomposition of given matrices. "))
32+
33+
.def("compute",
34+
&RealQZVisitor::compute_proxy<MatrixType>,
35+
bp::args("self", "A", "B"),
36+
"Computes QZ decomposition of given matrix. ",
37+
bp::return_self<>())
38+
.def("compute",
39+
(Solver &
40+
(Solver::*)(const MatrixType& A, const MatrixType& B, bool)) &
41+
Solver::compute,
42+
bp::args("self", "A", "B", "computeEigenvectors"),
43+
"Computes QZ decomposition of given matrix. ",
44+
bp::return_self<>())
45+
46+
.def("info", &Solver::info, bp::arg("self"),
47+
"NumericalIssue if the input contains INF or NaN values or "
48+
"overflow occured. Returns Success otherwise.")
49+
50+
.def("matrixQ", &Solver::matrixQ, bp::arg("self"),
51+
"Returns matrix Q in the QZ decomposition. ",
52+
bp::return_value_policy<bp::copy_const_reference>())
53+
.def("matrixS", &Solver::matrixS, bp::arg("self"),
54+
"Returns matrix S in the QZ decomposition. ",
55+
bp::return_value_policy<bp::copy_const_reference>())
56+
.def("matrixT", &Solver::matrixT, bp::arg("self"),
57+
"Returns matrix T in the QZ decomposition. ",
58+
bp::return_value_policy<bp::copy_const_reference>())
59+
.def("matrixZ", &Solver::matrixZ, bp::arg("self"),
60+
"Returns matrix Z in the QZ decomposition. ",
61+
bp::return_value_policy<bp::copy_const_reference>())
62+
63+
.def("iterations", &Solver::iterations, bp::arg("self"),
64+
"Returns number of performed QR-like iterations. ")
65+
.def("setMaxIterations", &Solver::setMaxIterations,
66+
bp::args("self", "max_iter"),
67+
"Sets the maximum number of iterations allowed.",
68+
bp::return_self<>());
69+
}
70+
71+
static void expose() {
72+
static const std::string classname =
73+
"RealQZVisitorSolver" + 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(RealQZVisitor())
80+
.def(IdVisitor<Solver>());
81+
}
82+
83+
private:
84+
template <typename MatrixType>
85+
static Solver& compute_proxy(Solver& self, const MatrixType& A,
86+
const MatrixType& B) {
87+
return self.compute(A, B);
88+
}
89+
};
90+
91+
} // namespace eigenpy
92+
93+
#endif // ifndef __eigenpy_decompositions_generalized_real_qz_hpp__

src/decompositions/decompositions.cpp

Lines changed: 2 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -13,6 +13,7 @@ void exposeGeneralizedEigenSolver();
1313
void exposeSelfAdjointEigenSolver();
1414
void exposeGeneralizedSelfAdjointEigenSolver();
1515
void exposeHessenbergDecomposition();
16+
void exposeRealQZ();
1617
void exposeComplexEigenSolver();
1718
void exposeComplexSchur();
1819
void exposeLLTSolver();
@@ -37,6 +38,7 @@ void exposeDecompositions() {
3738
exposeSelfAdjointEigenSolver();
3839
exposeGeneralizedSelfAdjointEigenSolver();
3940
exposeHessenbergDecomposition();
41+
exposeRealQZ();
4042
exposeComplexEigenSolver();
4143
exposeComplexSchur();
4244
exposeLLTSolver();

src/decompositions/real-qz.cpp

Lines changed: 13 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,13 @@
1+
2+
/*
3+
* Copyright 2024 INRIA
4+
*/
5+
6+
#include "eigenpy/decompositions/RealQZ.hpp"
7+
8+
namespace eigenpy {
9+
void exposeRealQZ() {
10+
using namespace Eigen;
11+
RealQZVisitor<MatrixXd>::expose("RealQZ");
12+
}
13+
} // namespace eigenpy

unittest/CMakeLists.txt

Lines changed: 2 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -159,6 +159,8 @@ add_python_eigenpy_lib_unit_test(
159159
"py-hessenberg-decomposition"
160160
"unittest/python/test_hessenberg_decomposition.py")
161161

162+
add_python_eigenpy_lib_unit_test("py-real-qz" "unittest/python/test_real_qz.py")
163+
162164
add_python_eigenpy_lib_unit_test("py-LLT" "unittest/python/test_LLT.py")
163165

164166
add_python_eigenpy_lib_unit_test("py-LDLT" "unittest/python/test_LDLT.py")

unittest/python/test_real_qz.py

Lines changed: 18 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,18 @@
1+
import numpy as np
2+
3+
import eigenpy
4+
5+
dim = 100
6+
rng = np.random.default_rng()
7+
A = rng.random((dim, dim))
8+
B = rng.random((dim, dim))
9+
10+
realqz = eigenpy.RealQZ(A, B)
11+
12+
Q = realqz.matrixQ()
13+
S = realqz.matrixS()
14+
Z = realqz.matrixZ()
15+
T = realqz.matrixT()
16+
17+
assert eigenpy.is_approx(A, Q @ S @ Z)
18+
assert eigenpy.is_approx(B, Q @ T @ Z)

0 commit comments

Comments
 (0)