1 // This file is part of Eigen, a lightweight C++ template library 2 // for linear algebra. 3 // 4 // Copyright (C) 2012 Alexey Korepanov <kaikaikai (at) yandex.ru> 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 15 template<typename MatrixType> void real_qz(const MatrixType& m) 16 { 17 /* this test covers the following files: 18 RealQZ.h 19 */ 20 using std::abs; 21 typedef typename MatrixType::Index Index; 22 typedef typename MatrixType::Scalar Scalar; 23 24 Index dim = m.cols(); 25 26 MatrixType A = MatrixType::Random(dim,dim), 27 B = MatrixType::Random(dim,dim); 28 29 30 // Regression test for bug 985: Randomly set rows or columns to zero 31 Index k=internal::random<Index>(0, dim-1); 32 switch(internal::random<int>(0,10)) { 33 case 0: 34 A.row(k).setZero(); break; 35 case 1: 36 A.col(k).setZero(); break; 37 case 2: 38 B.row(k).setZero(); break; 39 case 3: 40 B.col(k).setZero(); break; 41 default: 42 break; 43 } 44 45 RealQZ<MatrixType> qz(dim); 46 // TODO enable full-prealocation of required memory, this probably requires an in-place mode for HessenbergDecomposition 47 //Eigen::internal::set_is_malloc_allowed(false); 48 qz.compute(A,B); 49 //Eigen::internal::set_is_malloc_allowed(true); 50 51 VERIFY_IS_EQUAL(qz.info(), Success); 52 // check for zeros 53 bool all_zeros = true; 54 for (Index i=0; i<A.cols(); i++) 55 for (Index j=0; j<i; j++) { 56 if (abs(qz.matrixT()(i,j))!=Scalar(0.0)) 57 { 58 std::cerr << "Error: T(" << i << "," << j << ") = " << qz.matrixT()(i,j) << std::endl; 59 all_zeros = false; 60 } 61 if (j<i-1 && abs(qz.matrixS()(i,j))!=Scalar(0.0)) 62 { 63 std::cerr << "Error: S(" << i << "," << j << ") = " << qz.matrixS()(i,j) << std::endl; 64 all_zeros = false; 65 } 66 if (j==i-1 && j>0 && abs(qz.matrixS()(i,j))!=Scalar(0.0) && abs(qz.matrixS()(i-1,j-1))!=Scalar(0.0)) 67 { 68 std::cerr << "Error: S(" << i << "," << j << ") = " << qz.matrixS()(i,j) << " && S(" << i-1 << "," << j-1 << ") = " << qz.matrixS()(i-1,j-1) << std::endl; 69 all_zeros = false; 70 } 71 } 72 VERIFY_IS_EQUAL(all_zeros, true); 73 VERIFY_IS_APPROX(qz.matrixQ()*qz.matrixS()*qz.matrixZ(), A); 74 VERIFY_IS_APPROX(qz.matrixQ()*qz.matrixT()*qz.matrixZ(), B); 75 VERIFY_IS_APPROX(qz.matrixQ()*qz.matrixQ().adjoint(), MatrixType::Identity(dim,dim)); 76 VERIFY_IS_APPROX(qz.matrixZ()*qz.matrixZ().adjoint(), MatrixType::Identity(dim,dim)); 77 } 78 79 void test_real_qz() 80 { 81 int s = 0; 82 for(int i = 0; i < g_repeat; i++) { 83 CALL_SUBTEST_1( real_qz(Matrix4f()) ); 84 s = internal::random<int>(1,EIGEN_TEST_MAX_SIZE/4); 85 CALL_SUBTEST_2( real_qz(MatrixXd(s,s)) ); 86 87 // some trivial but implementation-wise tricky cases 88 CALL_SUBTEST_2( real_qz(MatrixXd(1,1)) ); 89 CALL_SUBTEST_2( real_qz(MatrixXd(2,2)) ); 90 CALL_SUBTEST_3( real_qz(Matrix<double,1,1>()) ); 91 CALL_SUBTEST_4( real_qz(Matrix2d()) ); 92 } 93 94 TEST_SET_BUT_UNUSED_VARIABLE(s) 95 } 96