schur_real.cpp (3964B)
1 // This file is part of Eigen, a lightweight C++ template library 2 // for linear algebra. 3 // 4 // Copyright (C) 2010,2012 Jitse Niesen <jitse@maths.leeds.ac.uk> 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 #include "main.h" 11 #include <limits> 12 #include <Eigen/Eigenvalues> 13 14 template<typename MatrixType> void verifyIsQuasiTriangular(const MatrixType& T) 15 { 16 const Index size = T.cols(); 17 typedef typename MatrixType::Scalar Scalar; 18 19 // Check T is lower Hessenberg 20 for(int row = 2; row < size; ++row) { 21 for(int col = 0; col < row - 1; ++col) { 22 VERIFY(T(row,col) == Scalar(0)); 23 } 24 } 25 26 // Check that any non-zero on the subdiagonal is followed by a zero and is 27 // part of a 2x2 diagonal block with imaginary eigenvalues. 28 for(int row = 1; row < size; ++row) { 29 if (T(row,row-1) != Scalar(0)) { 30 VERIFY(row == size-1 || T(row+1,row) == 0); 31 Scalar tr = T(row-1,row-1) + T(row,row); 32 Scalar det = T(row-1,row-1) * T(row,row) - T(row-1,row) * T(row,row-1); 33 VERIFY(4 * det > tr * tr); 34 } 35 } 36 } 37 38 template<typename MatrixType> void schur(int size = MatrixType::ColsAtCompileTime) 39 { 40 // Test basic functionality: T is quasi-triangular and A = U T U* 41 for(int counter = 0; counter < g_repeat; ++counter) { 42 MatrixType A = MatrixType::Random(size, size); 43 RealSchur<MatrixType> schurOfA(A); 44 VERIFY_IS_EQUAL(schurOfA.info(), Success); 45 MatrixType U = schurOfA.matrixU(); 46 MatrixType T = schurOfA.matrixT(); 47 verifyIsQuasiTriangular(T); 48 VERIFY_IS_APPROX(A, U * T * U.transpose()); 49 } 50 51 // Test asserts when not initialized 52 RealSchur<MatrixType> rsUninitialized; 53 VERIFY_RAISES_ASSERT(rsUninitialized.matrixT()); 54 VERIFY_RAISES_ASSERT(rsUninitialized.matrixU()); 55 VERIFY_RAISES_ASSERT(rsUninitialized.info()); 56 57 // Test whether compute() and constructor returns same result 58 MatrixType A = MatrixType::Random(size, size); 59 RealSchur<MatrixType> rs1; 60 rs1.compute(A); 61 RealSchur<MatrixType> rs2(A); 62 VERIFY_IS_EQUAL(rs1.info(), Success); 63 VERIFY_IS_EQUAL(rs2.info(), Success); 64 VERIFY_IS_EQUAL(rs1.matrixT(), rs2.matrixT()); 65 VERIFY_IS_EQUAL(rs1.matrixU(), rs2.matrixU()); 66 67 // Test maximum number of iterations 68 RealSchur<MatrixType> rs3; 69 rs3.setMaxIterations(RealSchur<MatrixType>::m_maxIterationsPerRow * size).compute(A); 70 VERIFY_IS_EQUAL(rs3.info(), Success); 71 VERIFY_IS_EQUAL(rs3.matrixT(), rs1.matrixT()); 72 VERIFY_IS_EQUAL(rs3.matrixU(), rs1.matrixU()); 73 if (size > 2) { 74 rs3.setMaxIterations(1).compute(A); 75 VERIFY_IS_EQUAL(rs3.info(), NoConvergence); 76 VERIFY_IS_EQUAL(rs3.getMaxIterations(), 1); 77 } 78 79 MatrixType Atriangular = A; 80 Atriangular.template triangularView<StrictlyLower>().setZero(); 81 rs3.setMaxIterations(1).compute(Atriangular); // triangular matrices do not need any iterations 82 VERIFY_IS_EQUAL(rs3.info(), Success); 83 VERIFY_IS_APPROX(rs3.matrixT(), Atriangular); // approx because of scaling... 84 VERIFY_IS_EQUAL(rs3.matrixU(), MatrixType::Identity(size, size)); 85 86 // Test computation of only T, not U 87 RealSchur<MatrixType> rsOnlyT(A, false); 88 VERIFY_IS_EQUAL(rsOnlyT.info(), Success); 89 VERIFY_IS_EQUAL(rs1.matrixT(), rsOnlyT.matrixT()); 90 VERIFY_RAISES_ASSERT(rsOnlyT.matrixU()); 91 92 if (size > 2 && size < 20) 93 { 94 // Test matrix with NaN 95 A(0,0) = std::numeric_limits<typename MatrixType::Scalar>::quiet_NaN(); 96 RealSchur<MatrixType> rsNaN(A); 97 VERIFY_IS_EQUAL(rsNaN.info(), NoConvergence); 98 } 99 } 100 101 EIGEN_DECLARE_TEST(schur_real) 102 { 103 CALL_SUBTEST_1(( schur<Matrix4f>() )); 104 CALL_SUBTEST_2(( schur<MatrixXd>(internal::random<int>(1,EIGEN_TEST_MAX_SIZE/4)) )); 105 CALL_SUBTEST_3(( schur<Matrix<float, 1, 1> >() )); 106 CALL_SUBTEST_4(( schur<Matrix<double, 3, 3, Eigen::RowMajor> >() )); 107 108 // Test problem size constructors 109 CALL_SUBTEST_5(RealSchur<MatrixXf>(10)); 110 }