1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2008-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 <limits>
13 #include <Eigen/Eigenvalues>
14 #include <Eigen/LU>
15
find_pivot(typename MatrixType::Scalar tol,MatrixType & diffs,Index col=0)16 template<typename MatrixType> bool find_pivot(typename MatrixType::Scalar tol, MatrixType &diffs, Index col=0)
17 {
18 bool match = diffs.diagonal().sum() <= tol;
19 if(match || col==diffs.cols())
20 {
21 return match;
22 }
23 else
24 {
25 Index n = diffs.cols();
26 std::vector<std::pair<Index,Index> > transpositions;
27 for(Index i=col; i<n; ++i)
28 {
29 Index best_index(0);
30 if(diffs.col(col).segment(col,n-i).minCoeff(&best_index) > tol)
31 break;
32
33 best_index += col;
34
35 diffs.row(col).swap(diffs.row(best_index));
36 if(find_pivot(tol,diffs,col+1)) return true;
37 diffs.row(col).swap(diffs.row(best_index));
38
39 // move current pivot to the end
40 diffs.row(n-(i-col)-1).swap(diffs.row(best_index));
41 transpositions.push_back(std::pair<Index,Index>(n-(i-col)-1,best_index));
42 }
43 // restore
44 for(Index k=transpositions.size()-1; k>=0; --k)
45 diffs.row(transpositions[k].first).swap(diffs.row(transpositions[k].second));
46 }
47 return false;
48 }
49
50 /* Check that two column vectors are approximately equal upto permutations.
51 * Initially, this method checked that the k-th power sums are equal for all k = 1, ..., vec1.rows(),
52 * however this strategy is numerically inacurate because of numerical cancellation issues.
53 */
54 template<typename VectorType>
verify_is_approx_upto_permutation(const VectorType & vec1,const VectorType & vec2)55 void verify_is_approx_upto_permutation(const VectorType& vec1, const VectorType& vec2)
56 {
57 typedef typename VectorType::Scalar Scalar;
58 typedef typename NumTraits<Scalar>::Real RealScalar;
59
60 VERIFY(vec1.cols() == 1);
61 VERIFY(vec2.cols() == 1);
62 VERIFY(vec1.rows() == vec2.rows());
63
64 Index n = vec1.rows();
65 RealScalar tol = test_precision<RealScalar>()*test_precision<RealScalar>()*numext::maxi(vec1.squaredNorm(),vec2.squaredNorm());
66 Matrix<RealScalar,Dynamic,Dynamic> diffs = (vec1.rowwise().replicate(n) - vec2.rowwise().replicate(n).transpose()).cwiseAbs2();
67
68 VERIFY( find_pivot(tol, diffs) );
69 }
70
71
eigensolver(const MatrixType & m)72 template<typename MatrixType> void eigensolver(const MatrixType& m)
73 {
74 typedef typename MatrixType::Index Index;
75 /* this test covers the following files:
76 ComplexEigenSolver.h, and indirectly ComplexSchur.h
77 */
78 Index rows = m.rows();
79 Index cols = m.cols();
80
81 typedef typename MatrixType::Scalar Scalar;
82 typedef typename NumTraits<Scalar>::Real RealScalar;
83
84 MatrixType a = MatrixType::Random(rows,cols);
85 MatrixType symmA = a.adjoint() * a;
86
87 ComplexEigenSolver<MatrixType> ei0(symmA);
88 VERIFY_IS_EQUAL(ei0.info(), Success);
89 VERIFY_IS_APPROX(symmA * ei0.eigenvectors(), ei0.eigenvectors() * ei0.eigenvalues().asDiagonal());
90
91 ComplexEigenSolver<MatrixType> ei1(a);
92 VERIFY_IS_EQUAL(ei1.info(), Success);
93 VERIFY_IS_APPROX(a * ei1.eigenvectors(), ei1.eigenvectors() * ei1.eigenvalues().asDiagonal());
94 // Note: If MatrixType is real then a.eigenvalues() uses EigenSolver and thus
95 // another algorithm so results may differ slightly
96 verify_is_approx_upto_permutation(a.eigenvalues(), ei1.eigenvalues());
97
98 ComplexEigenSolver<MatrixType> ei2;
99 ei2.setMaxIterations(ComplexSchur<MatrixType>::m_maxIterationsPerRow * rows).compute(a);
100 VERIFY_IS_EQUAL(ei2.info(), Success);
101 VERIFY_IS_EQUAL(ei2.eigenvectors(), ei1.eigenvectors());
102 VERIFY_IS_EQUAL(ei2.eigenvalues(), ei1.eigenvalues());
103 if (rows > 2) {
104 ei2.setMaxIterations(1).compute(a);
105 VERIFY_IS_EQUAL(ei2.info(), NoConvergence);
106 VERIFY_IS_EQUAL(ei2.getMaxIterations(), 1);
107 }
108
109 ComplexEigenSolver<MatrixType> eiNoEivecs(a, false);
110 VERIFY_IS_EQUAL(eiNoEivecs.info(), Success);
111 VERIFY_IS_APPROX(ei1.eigenvalues(), eiNoEivecs.eigenvalues());
112
113 // Regression test for issue #66
114 MatrixType z = MatrixType::Zero(rows,cols);
115 ComplexEigenSolver<MatrixType> eiz(z);
116 VERIFY((eiz.eigenvalues().cwiseEqual(0)).all());
117
118 MatrixType id = MatrixType::Identity(rows, cols);
119 VERIFY_IS_APPROX(id.operatorNorm(), RealScalar(1));
120
121 if (rows > 1 && rows < 20)
122 {
123 // Test matrix with NaN
124 a(0,0) = std::numeric_limits<typename MatrixType::RealScalar>::quiet_NaN();
125 ComplexEigenSolver<MatrixType> eiNaN(a);
126 VERIFY_IS_EQUAL(eiNaN.info(), NoConvergence);
127 }
128
129 // regression test for bug 1098
130 {
131 ComplexEigenSolver<MatrixType> eig(a.adjoint() * a);
132 eig.compute(a.adjoint() * a);
133 }
134
135 // regression test for bug 478
136 {
137 a.setZero();
138 ComplexEigenSolver<MatrixType> ei3(a);
139 VERIFY_IS_EQUAL(ei3.info(), Success);
140 VERIFY_IS_MUCH_SMALLER_THAN(ei3.eigenvalues().norm(),RealScalar(1));
141 VERIFY((ei3.eigenvectors().transpose()*ei3.eigenvectors().transpose()).eval().isIdentity());
142 }
143 }
144
eigensolver_verify_assert(const MatrixType & m)145 template<typename MatrixType> void eigensolver_verify_assert(const MatrixType& m)
146 {
147 ComplexEigenSolver<MatrixType> eig;
148 VERIFY_RAISES_ASSERT(eig.eigenvectors());
149 VERIFY_RAISES_ASSERT(eig.eigenvalues());
150
151 MatrixType a = MatrixType::Random(m.rows(),m.cols());
152 eig.compute(a, false);
153 VERIFY_RAISES_ASSERT(eig.eigenvectors());
154 }
155
test_eigensolver_complex()156 void test_eigensolver_complex()
157 {
158 int s = 0;
159 for(int i = 0; i < g_repeat; i++) {
160 CALL_SUBTEST_1( eigensolver(Matrix4cf()) );
161 s = internal::random<int>(1,EIGEN_TEST_MAX_SIZE/4);
162 CALL_SUBTEST_2( eigensolver(MatrixXcd(s,s)) );
163 CALL_SUBTEST_3( eigensolver(Matrix<std::complex<float>, 1, 1>()) );
164 CALL_SUBTEST_4( eigensolver(Matrix3f()) );
165 TEST_SET_BUT_UNUSED_VARIABLE(s)
166 }
167 CALL_SUBTEST_1( eigensolver_verify_assert(Matrix4cf()) );
168 s = internal::random<int>(1,EIGEN_TEST_MAX_SIZE/4);
169 CALL_SUBTEST_2( eigensolver_verify_assert(MatrixXcd(s,s)) );
170 CALL_SUBTEST_3( eigensolver_verify_assert(Matrix<std::complex<float>, 1, 1>()) );
171 CALL_SUBTEST_4( eigensolver_verify_assert(Matrix3f()) );
172
173 // Test problem size constructors
174 CALL_SUBTEST_5(ComplexEigenSolver<MatrixXf> tmp(s));
175
176 TEST_SET_BUT_UNUSED_VARIABLE(s)
177 }
178