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