• 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) 2008 Gael Guennebaud <gael.guennebaud@inria.fr>
5 // Copyright (C) 2010,2012 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 
eigensolver(const MatrixType & m)15 template<typename MatrixType> void eigensolver(const MatrixType& m)
16 {
17   /* this test covers the following files:
18      EigenSolver.h
19   */
20   Index rows = m.rows();
21   Index cols = m.cols();
22 
23   typedef typename MatrixType::Scalar Scalar;
24   typedef typename NumTraits<Scalar>::Real RealScalar;
25   typedef Matrix<RealScalar, MatrixType::RowsAtCompileTime, 1> RealVectorType;
26   typedef typename std::complex<typename NumTraits<typename MatrixType::Scalar>::Real> Complex;
27 
28   MatrixType a = MatrixType::Random(rows,cols);
29   MatrixType a1 = MatrixType::Random(rows,cols);
30   MatrixType symmA =  a.adjoint() * a + a1.adjoint() * a1;
31 
32   EigenSolver<MatrixType> ei0(symmA);
33   VERIFY_IS_EQUAL(ei0.info(), Success);
34   VERIFY_IS_APPROX(symmA * ei0.pseudoEigenvectors(), ei0.pseudoEigenvectors() * ei0.pseudoEigenvalueMatrix());
35   VERIFY_IS_APPROX((symmA.template cast<Complex>()) * (ei0.pseudoEigenvectors().template cast<Complex>()),
36     (ei0.pseudoEigenvectors().template cast<Complex>()) * (ei0.eigenvalues().asDiagonal()));
37 
38   EigenSolver<MatrixType> ei1(a);
39   VERIFY_IS_EQUAL(ei1.info(), Success);
40   VERIFY_IS_APPROX(a * ei1.pseudoEigenvectors(), ei1.pseudoEigenvectors() * ei1.pseudoEigenvalueMatrix());
41   VERIFY_IS_APPROX(a.template cast<Complex>() * ei1.eigenvectors(),
42                    ei1.eigenvectors() * ei1.eigenvalues().asDiagonal());
43   VERIFY_IS_APPROX(ei1.eigenvectors().colwise().norm(), RealVectorType::Ones(rows).transpose());
44   VERIFY_IS_APPROX(a.eigenvalues(), ei1.eigenvalues());
45 
46   EigenSolver<MatrixType> ei2;
47   ei2.setMaxIterations(RealSchur<MatrixType>::m_maxIterationsPerRow * rows).compute(a);
48   VERIFY_IS_EQUAL(ei2.info(), Success);
49   VERIFY_IS_EQUAL(ei2.eigenvectors(), ei1.eigenvectors());
50   VERIFY_IS_EQUAL(ei2.eigenvalues(), ei1.eigenvalues());
51   if (rows > 2) {
52     ei2.setMaxIterations(1).compute(a);
53     VERIFY_IS_EQUAL(ei2.info(), NoConvergence);
54     VERIFY_IS_EQUAL(ei2.getMaxIterations(), 1);
55   }
56 
57   EigenSolver<MatrixType> eiNoEivecs(a, false);
58   VERIFY_IS_EQUAL(eiNoEivecs.info(), Success);
59   VERIFY_IS_APPROX(ei1.eigenvalues(), eiNoEivecs.eigenvalues());
60   VERIFY_IS_APPROX(ei1.pseudoEigenvalueMatrix(), eiNoEivecs.pseudoEigenvalueMatrix());
61 
62   MatrixType id = MatrixType::Identity(rows, cols);
63   VERIFY_IS_APPROX(id.operatorNorm(), RealScalar(1));
64 
65   if (rows > 2 && rows < 20)
66   {
67     // Test matrix with NaN
68     a(0,0) = std::numeric_limits<typename MatrixType::RealScalar>::quiet_NaN();
69     EigenSolver<MatrixType> eiNaN(a);
70     VERIFY_IS_EQUAL(eiNaN.info(), NoConvergence);
71   }
72 
73   // regression test for bug 1098
74   {
75     EigenSolver<MatrixType> eig(a.adjoint() * a);
76     eig.compute(a.adjoint() * a);
77   }
78 
79   // regression test for bug 478
80   {
81     a.setZero();
82     EigenSolver<MatrixType> ei3(a);
83     VERIFY_IS_EQUAL(ei3.info(), Success);
84     VERIFY_IS_MUCH_SMALLER_THAN(ei3.eigenvalues().norm(),RealScalar(1));
85     VERIFY((ei3.eigenvectors().transpose()*ei3.eigenvectors().transpose()).eval().isIdentity());
86   }
87 }
88 
eigensolver_verify_assert(const MatrixType & m)89 template<typename MatrixType> void eigensolver_verify_assert(const MatrixType& m)
90 {
91   EigenSolver<MatrixType> eig;
92   VERIFY_RAISES_ASSERT(eig.eigenvectors());
93   VERIFY_RAISES_ASSERT(eig.pseudoEigenvectors());
94   VERIFY_RAISES_ASSERT(eig.pseudoEigenvalueMatrix());
95   VERIFY_RAISES_ASSERT(eig.eigenvalues());
96 
97   MatrixType a = MatrixType::Random(m.rows(),m.cols());
98   eig.compute(a, false);
99   VERIFY_RAISES_ASSERT(eig.eigenvectors());
100   VERIFY_RAISES_ASSERT(eig.pseudoEigenvectors());
101 }
102 
test_eigensolver_generic()103 void test_eigensolver_generic()
104 {
105   int s = 0;
106   for(int i = 0; i < g_repeat; i++) {
107     CALL_SUBTEST_1( eigensolver(Matrix4f()) );
108     s = internal::random<int>(1,EIGEN_TEST_MAX_SIZE/4);
109     CALL_SUBTEST_2( eigensolver(MatrixXd(s,s)) );
110     TEST_SET_BUT_UNUSED_VARIABLE(s)
111 
112     // some trivial but implementation-wise tricky cases
113     CALL_SUBTEST_2( eigensolver(MatrixXd(1,1)) );
114     CALL_SUBTEST_2( eigensolver(MatrixXd(2,2)) );
115     CALL_SUBTEST_3( eigensolver(Matrix<double,1,1>()) );
116     CALL_SUBTEST_4( eigensolver(Matrix2d()) );
117   }
118 
119   CALL_SUBTEST_1( eigensolver_verify_assert(Matrix4f()) );
120   s = internal::random<int>(1,EIGEN_TEST_MAX_SIZE/4);
121   CALL_SUBTEST_2( eigensolver_verify_assert(MatrixXd(s,s)) );
122   CALL_SUBTEST_3( eigensolver_verify_assert(Matrix<double,1,1>()) );
123   CALL_SUBTEST_4( eigensolver_verify_assert(Matrix2d()) );
124 
125   // Test problem size constructors
126   CALL_SUBTEST_5(EigenSolver<MatrixXf> tmp(s));
127 
128   // regression test for bug 410
129   CALL_SUBTEST_2(
130   {
131      MatrixXd A(1,1);
132      A(0,0) = std::sqrt(-1.); // is Not-a-Number
133      Eigen::EigenSolver<MatrixXd> solver(A);
134      VERIFY_IS_EQUAL(solver.info(), NumericalIssue);
135   }
136   );
137 
138 #ifdef EIGEN_TEST_PART_2
139   {
140     // regression test for bug 793
141     MatrixXd a(3,3);
142     a << 0,  0,  1,
143         1,  1, 1,
144         1, 1e+200,  1;
145     Eigen::EigenSolver<MatrixXd> eig(a);
146     double scale = 1e-200; // scale to avoid overflow during the comparisons
147     VERIFY_IS_APPROX(a * eig.pseudoEigenvectors()*scale, eig.pseudoEigenvectors() * eig.pseudoEigenvalueMatrix()*scale);
148     VERIFY_IS_APPROX(a * eig.eigenvectors()*scale, eig.eigenvectors() * eig.eigenvalues().asDiagonal()*scale);
149   }
150   {
151     // check a case where all eigenvalues are null.
152     MatrixXd a(2,2);
153     a << 1,  1,
154         -1, -1;
155     Eigen::EigenSolver<MatrixXd> eig(a);
156     VERIFY_IS_APPROX(eig.pseudoEigenvectors().squaredNorm(), 2.);
157     VERIFY_IS_APPROX((a * eig.pseudoEigenvectors()).norm()+1., 1.);
158     VERIFY_IS_APPROX((eig.pseudoEigenvectors() * eig.pseudoEigenvalueMatrix()).norm()+1., 1.);
159     VERIFY_IS_APPROX((a * eig.eigenvectors()).norm()+1., 1.);
160     VERIFY_IS_APPROX((eig.eigenvectors() * eig.eigenvalues().asDiagonal()).norm()+1., 1.);
161   }
162 #endif
163 
164   TEST_SET_BUT_UNUSED_VARIABLE(s)
165 }
166