1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2012-2016 Gael Guennebaud <[email protected]>
5 //
6 // This Source Code Form is subject to the terms of the Mozilla
7 // Public License v. 2.0. If a copy of the MPL was not distributed
8 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
9
10 #define EIGEN_RUNTIME_NO_MALLOC
11 #include "main.h"
12 #include <limits>
13 #include <Eigen/Eigenvalues>
14 #include <Eigen/LU>
15
generalized_eigensolver_real(const MatrixType & m)16 template<typename MatrixType> void generalized_eigensolver_real(const MatrixType& m)
17 {
18 /* this test covers the following files:
19 GeneralizedEigenSolver.h
20 */
21 Index rows = m.rows();
22 Index cols = m.cols();
23
24 typedef typename MatrixType::Scalar Scalar;
25 typedef std::complex<Scalar> ComplexScalar;
26 typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType;
27
28 MatrixType a = MatrixType::Random(rows,cols);
29 MatrixType b = MatrixType::Random(rows,cols);
30 MatrixType a1 = MatrixType::Random(rows,cols);
31 MatrixType b1 = MatrixType::Random(rows,cols);
32 MatrixType spdA = a.adjoint() * a + a1.adjoint() * a1;
33 MatrixType spdB = b.adjoint() * b + b1.adjoint() * b1;
34
35 // lets compare to GeneralizedSelfAdjointEigenSolver
36 {
37 GeneralizedSelfAdjointEigenSolver<MatrixType> symmEig(spdA, spdB);
38 GeneralizedEigenSolver<MatrixType> eig(spdA, spdB);
39
40 VERIFY_IS_EQUAL(eig.eigenvalues().imag().cwiseAbs().maxCoeff(), 0);
41
42 VectorType realEigenvalues = eig.eigenvalues().real();
43 std::sort(realEigenvalues.data(), realEigenvalues.data()+realEigenvalues.size());
44 VERIFY_IS_APPROX(realEigenvalues, symmEig.eigenvalues());
45
46 // check eigenvectors
47 typename GeneralizedEigenSolver<MatrixType>::EigenvectorsType D = eig.eigenvalues().asDiagonal();
48 typename GeneralizedEigenSolver<MatrixType>::EigenvectorsType V = eig.eigenvectors();
49 VERIFY_IS_APPROX(spdA*V, spdB*V*D);
50 }
51
52 // non symmetric case:
53 {
54 GeneralizedEigenSolver<MatrixType> eig(rows);
55 // TODO enable full-prealocation of required memory, this probably requires an in-place mode for HessenbergDecomposition
56 //Eigen::internal::set_is_malloc_allowed(false);
57 eig.compute(a,b);
58 //Eigen::internal::set_is_malloc_allowed(true);
59 for(Index k=0; k<cols; ++k)
60 {
61 Matrix<ComplexScalar,Dynamic,Dynamic> tmp = (eig.betas()(k)*a).template cast<ComplexScalar>() - eig.alphas()(k)*b;
62 if(tmp.size()>1 && tmp.norm()>(std::numeric_limits<Scalar>::min)())
63 tmp /= tmp.norm();
64 VERIFY_IS_MUCH_SMALLER_THAN( std::abs(tmp.determinant()), Scalar(1) );
65 }
66 // check eigenvectors
67 typename GeneralizedEigenSolver<MatrixType>::EigenvectorsType D = eig.eigenvalues().asDiagonal();
68 typename GeneralizedEigenSolver<MatrixType>::EigenvectorsType V = eig.eigenvectors();
69 VERIFY_IS_APPROX(a*V, b*V*D);
70 }
71
72 // regression test for bug 1098
73 {
74 GeneralizedSelfAdjointEigenSolver<MatrixType> eig1(a.adjoint() * a,b.adjoint() * b);
75 eig1.compute(a.adjoint() * a,b.adjoint() * b);
76 GeneralizedEigenSolver<MatrixType> eig2(a.adjoint() * a,b.adjoint() * b);
77 eig2.compute(a.adjoint() * a,b.adjoint() * b);
78 }
79
80 // check without eigenvectors
81 {
82 GeneralizedEigenSolver<MatrixType> eig1(spdA, spdB, true);
83 GeneralizedEigenSolver<MatrixType> eig2(spdA, spdB, false);
84 VERIFY_IS_APPROX(eig1.eigenvalues(), eig2.eigenvalues());
85 }
86 }
87
EIGEN_DECLARE_TEST(eigensolver_generalized_real)88 EIGEN_DECLARE_TEST(eigensolver_generalized_real)
89 {
90 for(int i = 0; i < g_repeat; i++) {
91 int s = 0;
92 CALL_SUBTEST_1( generalized_eigensolver_real(Matrix4f()) );
93 s = internal::random<int>(1,EIGEN_TEST_MAX_SIZE/4);
94 CALL_SUBTEST_2( generalized_eigensolver_real(MatrixXd(s,s)) );
95
96 // some trivial but implementation-wise special cases
97 CALL_SUBTEST_2( generalized_eigensolver_real(MatrixXd(1,1)) );
98 CALL_SUBTEST_2( generalized_eigensolver_real(MatrixXd(2,2)) );
99 CALL_SUBTEST_3( generalized_eigensolver_real(Matrix<double,1,1>()) );
100 CALL_SUBTEST_4( generalized_eigensolver_real(Matrix2d()) );
101 TEST_SET_BUT_UNUSED_VARIABLE(s)
102 }
103 }
104