1// This file is part of Eigen, a lightweight C++ template library
2// for linear algebra.
3//
4// Copyright (C) 2008 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 <limits>
13#include <Eigen/Eigenvalues>
14
15template<typename MatrixType> void selfadjointeigensolver(const MatrixType& m)
16{
17  typedef typename MatrixType::Index Index;
18  /* this test covers the following files:
19     EigenSolver.h, SelfAdjointEigenSolver.h (and indirectly: Tridiagonalization.h)
20  */
21  Index rows = m.rows();
22  Index cols = m.cols();
23
24  typedef typename MatrixType::Scalar Scalar;
25  typedef typename NumTraits<Scalar>::Real RealScalar;
26
27  RealScalar largerEps = 10*test_precision<RealScalar>();
28
29  MatrixType a = MatrixType::Random(rows,cols);
30  MatrixType a1 = MatrixType::Random(rows,cols);
31  MatrixType symmA =  a.adjoint() * a + a1.adjoint() * a1;
32  symmA.template triangularView<StrictlyUpper>().setZero();
33
34  MatrixType b = MatrixType::Random(rows,cols);
35  MatrixType b1 = MatrixType::Random(rows,cols);
36  MatrixType symmB = b.adjoint() * b + b1.adjoint() * b1;
37  symmB.template triangularView<StrictlyUpper>().setZero();
38
39  SelfAdjointEigenSolver<MatrixType> eiSymm(symmA);
40  SelfAdjointEigenSolver<MatrixType> eiDirect;
41  eiDirect.computeDirect(symmA);
42  // generalized eigen pb
43  GeneralizedSelfAdjointEigenSolver<MatrixType> eiSymmGen(symmA, symmB);
44
45  VERIFY_IS_EQUAL(eiSymm.info(), Success);
46  VERIFY((symmA.template selfadjointView<Lower>() * eiSymm.eigenvectors()).isApprox(
47          eiSymm.eigenvectors() * eiSymm.eigenvalues().asDiagonal(), largerEps));
48  VERIFY_IS_APPROX(symmA.template selfadjointView<Lower>().eigenvalues(), eiSymm.eigenvalues());
49
50  VERIFY_IS_EQUAL(eiDirect.info(), Success);
51  VERIFY((symmA.template selfadjointView<Lower>() * eiDirect.eigenvectors()).isApprox(
52          eiDirect.eigenvectors() * eiDirect.eigenvalues().asDiagonal(), largerEps));
53  VERIFY_IS_APPROX(symmA.template selfadjointView<Lower>().eigenvalues(), eiDirect.eigenvalues());
54
55  SelfAdjointEigenSolver<MatrixType> eiSymmNoEivecs(symmA, false);
56  VERIFY_IS_EQUAL(eiSymmNoEivecs.info(), Success);
57  VERIFY_IS_APPROX(eiSymm.eigenvalues(), eiSymmNoEivecs.eigenvalues());
58
59  // generalized eigen problem Ax = lBx
60  eiSymmGen.compute(symmA, symmB,Ax_lBx);
61  VERIFY_IS_EQUAL(eiSymmGen.info(), Success);
62  VERIFY((symmA.template selfadjointView<Lower>() * eiSymmGen.eigenvectors()).isApprox(
63          symmB.template selfadjointView<Lower>() * (eiSymmGen.eigenvectors() * eiSymmGen.eigenvalues().asDiagonal()), largerEps));
64
65  // generalized eigen problem BAx = lx
66  eiSymmGen.compute(symmA, symmB,BAx_lx);
67  VERIFY_IS_EQUAL(eiSymmGen.info(), Success);
68  VERIFY((symmB.template selfadjointView<Lower>() * (symmA.template selfadjointView<Lower>() * eiSymmGen.eigenvectors())).isApprox(
69         (eiSymmGen.eigenvectors() * eiSymmGen.eigenvalues().asDiagonal()), largerEps));
70
71  // generalized eigen problem ABx = lx
72  eiSymmGen.compute(symmA, symmB,ABx_lx);
73  VERIFY_IS_EQUAL(eiSymmGen.info(), Success);
74  VERIFY((symmA.template selfadjointView<Lower>() * (symmB.template selfadjointView<Lower>() * eiSymmGen.eigenvectors())).isApprox(
75         (eiSymmGen.eigenvectors() * eiSymmGen.eigenvalues().asDiagonal()), largerEps));
76
77
78  MatrixType sqrtSymmA = eiSymm.operatorSqrt();
79  VERIFY_IS_APPROX(MatrixType(symmA.template selfadjointView<Lower>()), sqrtSymmA*sqrtSymmA);
80  VERIFY_IS_APPROX(sqrtSymmA, symmA.template selfadjointView<Lower>()*eiSymm.operatorInverseSqrt());
81
82  MatrixType id = MatrixType::Identity(rows, cols);
83  VERIFY_IS_APPROX(id.template selfadjointView<Lower>().operatorNorm(), RealScalar(1));
84
85  SelfAdjointEigenSolver<MatrixType> eiSymmUninitialized;
86  VERIFY_RAISES_ASSERT(eiSymmUninitialized.info());
87  VERIFY_RAISES_ASSERT(eiSymmUninitialized.eigenvalues());
88  VERIFY_RAISES_ASSERT(eiSymmUninitialized.eigenvectors());
89  VERIFY_RAISES_ASSERT(eiSymmUninitialized.operatorSqrt());
90  VERIFY_RAISES_ASSERT(eiSymmUninitialized.operatorInverseSqrt());
91
92  eiSymmUninitialized.compute(symmA, false);
93  VERIFY_RAISES_ASSERT(eiSymmUninitialized.eigenvectors());
94  VERIFY_RAISES_ASSERT(eiSymmUninitialized.operatorSqrt());
95  VERIFY_RAISES_ASSERT(eiSymmUninitialized.operatorInverseSqrt());
96
97  // test Tridiagonalization's methods
98  Tridiagonalization<MatrixType> tridiag(symmA);
99  // FIXME tridiag.matrixQ().adjoint() does not work
100  VERIFY_IS_APPROX(MatrixType(symmA.template selfadjointView<Lower>()), tridiag.matrixQ() * tridiag.matrixT().eval() * MatrixType(tridiag.matrixQ()).adjoint());
101
102  if (rows > 1)
103  {
104    // Test matrix with NaN
105    symmA(0,0) = std::numeric_limits<typename MatrixType::RealScalar>::quiet_NaN();
106    SelfAdjointEigenSolver<MatrixType> eiSymmNaN(symmA);
107    VERIFY_IS_EQUAL(eiSymmNaN.info(), NoConvergence);
108  }
109}
110
111void test_eigensolver_selfadjoint()
112{
113  int s = 0;
114  for(int i = 0; i < g_repeat; i++) {
115    // very important to test 3x3 and 2x2 matrices since we provide special paths for them
116    CALL_SUBTEST_1( selfadjointeigensolver(Matrix2d()) );
117    CALL_SUBTEST_1( selfadjointeigensolver(Matrix3f()) );
118    CALL_SUBTEST_2( selfadjointeigensolver(Matrix4d()) );
119    s = internal::random<int>(1,EIGEN_TEST_MAX_SIZE/4);
120    CALL_SUBTEST_3( selfadjointeigensolver(MatrixXf(s,s)) );
121    s = internal::random<int>(1,EIGEN_TEST_MAX_SIZE/4);
122    CALL_SUBTEST_4( selfadjointeigensolver(MatrixXd(s,s)) );
123    s = internal::random<int>(1,EIGEN_TEST_MAX_SIZE/4);
124    CALL_SUBTEST_5( selfadjointeigensolver(MatrixXcd(s,s)) );
125
126    s = internal::random<int>(1,EIGEN_TEST_MAX_SIZE/4);
127    CALL_SUBTEST_9( selfadjointeigensolver(Matrix<std::complex<double>,Dynamic,Dynamic,RowMajor>(s,s)) );
128
129    // some trivial but implementation-wise tricky cases
130    CALL_SUBTEST_4( selfadjointeigensolver(MatrixXd(1,1)) );
131    CALL_SUBTEST_4( selfadjointeigensolver(MatrixXd(2,2)) );
132    CALL_SUBTEST_6( selfadjointeigensolver(Matrix<double,1,1>()) );
133    CALL_SUBTEST_7( selfadjointeigensolver(Matrix<double,2,2>()) );
134  }
135
136  // Test problem size constructors
137  s = internal::random<int>(1,EIGEN_TEST_MAX_SIZE/4);
138  CALL_SUBTEST_8(SelfAdjointEigenSolver<MatrixXf> tmp1(s));
139  CALL_SUBTEST_8(Tridiagonalization<MatrixXf> tmp2(s));
140
141  TEST_SET_BUT_UNUSED_VARIABLE(s)
142}
143
144