• Home
  • Line#
  • Scopes#
  • Navigate#
  • Raw
  • Download
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2009 Gael Guennebaud <gael.guennebaud@inria.fr>
5 // Copyright (C) 2010 Jitse Niesen <jitse@maths.leeds.ac.uk>
6 //
7 // This Source Code Form is subject to the terms of the Mozilla
8 // Public License v. 2.0. If a copy of the MPL was not distributed
9 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
10 
11 #include "main.h"
12 #include <Eigen/Eigenvalues>
13 
hessenberg(int size=Size)14 template<typename Scalar,int Size> void hessenberg(int size = Size)
15 {
16   typedef Matrix<Scalar,Size,Size> MatrixType;
17 
18   // Test basic functionality: A = U H U* and H is Hessenberg
19   for(int counter = 0; counter < g_repeat; ++counter) {
20     MatrixType m = MatrixType::Random(size,size);
21     HessenbergDecomposition<MatrixType> hess(m);
22     MatrixType Q = hess.matrixQ();
23     MatrixType H = hess.matrixH();
24     VERIFY_IS_APPROX(m, Q * H * Q.adjoint());
25     for(int row = 2; row < size; ++row) {
26       for(int col = 0; col < row-1; ++col) {
27 	VERIFY(H(row,col) == (typename MatrixType::Scalar)0);
28       }
29     }
30   }
31 
32   // Test whether compute() and constructor returns same result
33   MatrixType A = MatrixType::Random(size, size);
34   HessenbergDecomposition<MatrixType> cs1;
35   cs1.compute(A);
36   HessenbergDecomposition<MatrixType> cs2(A);
37   VERIFY_IS_EQUAL(cs1.matrixH().eval(), cs2.matrixH().eval());
38   MatrixType cs1Q = cs1.matrixQ();
39   MatrixType cs2Q = cs2.matrixQ();
40   VERIFY_IS_EQUAL(cs1Q, cs2Q);
41 
42   // Test assertions for when used uninitialized
43   HessenbergDecomposition<MatrixType> hessUninitialized;
44   VERIFY_RAISES_ASSERT( hessUninitialized.matrixH() );
45   VERIFY_RAISES_ASSERT( hessUninitialized.matrixQ() );
46   VERIFY_RAISES_ASSERT( hessUninitialized.householderCoefficients() );
47   VERIFY_RAISES_ASSERT( hessUninitialized.packedMatrix() );
48 
49   // TODO: Add tests for packedMatrix() and householderCoefficients()
50 }
51 
test_hessenberg()52 void test_hessenberg()
53 {
54   CALL_SUBTEST_1(( hessenberg<std::complex<double>,1>() ));
55   CALL_SUBTEST_2(( hessenberg<std::complex<double>,2>() ));
56   CALL_SUBTEST_3(( hessenberg<std::complex<float>,4>() ));
57   CALL_SUBTEST_4(( hessenberg<float,Dynamic>(internal::random<int>(1,EIGEN_TEST_MAX_SIZE)) ));
58   CALL_SUBTEST_5(( hessenberg<std::complex<double>,Dynamic>(internal::random<int>(1,EIGEN_TEST_MAX_SIZE)) ));
59 
60   // Test problem size constructors
61   CALL_SUBTEST_6(HessenbergDecomposition<MatrixXf>(10));
62 }
63