blob: 7269452598ffe5c7faf3852022e73b9d5983d275 [file] [log] [blame]
Brian Silverman72890c22015-09-19 14:37:37 -04001// 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
Austin Schuh189376f2018-12-20 22:11:15 +110016template<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 */
Brian Silverman72890c22015-09-19 14:37:37 -040054template<typename VectorType>
55void verify_is_approx_upto_permutation(const VectorType& vec1, const VectorType& vec2)
56{
Austin Schuh189376f2018-12-20 22:11:15 +110057 typedef typename VectorType::Scalar Scalar;
58 typedef typename NumTraits<Scalar>::Real RealScalar;
Brian Silverman72890c22015-09-19 14:37:37 -040059
60 VERIFY(vec1.cols() == 1);
61 VERIFY(vec2.cols() == 1);
62 VERIFY(vec1.rows() == vec2.rows());
Austin Schuh189376f2018-12-20 22:11:15 +110063
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) );
Brian Silverman72890c22015-09-19 14:37:37 -040069}
70
71
72template<typename MatrixType> void eigensolver(const MatrixType& m)
73{
Brian Silverman72890c22015-09-19 14:37:37 -040074 /* this test covers the following files:
75 ComplexEigenSolver.h, and indirectly ComplexSchur.h
76 */
77 Index rows = m.rows();
78 Index cols = m.cols();
79
80 typedef typename MatrixType::Scalar Scalar;
81 typedef typename NumTraits<Scalar>::Real RealScalar;
82
83 MatrixType a = MatrixType::Random(rows,cols);
84 MatrixType symmA = a.adjoint() * a;
85
86 ComplexEigenSolver<MatrixType> ei0(symmA);
87 VERIFY_IS_EQUAL(ei0.info(), Success);
88 VERIFY_IS_APPROX(symmA * ei0.eigenvectors(), ei0.eigenvectors() * ei0.eigenvalues().asDiagonal());
89
90 ComplexEigenSolver<MatrixType> ei1(a);
91 VERIFY_IS_EQUAL(ei1.info(), Success);
92 VERIFY_IS_APPROX(a * ei1.eigenvectors(), ei1.eigenvectors() * ei1.eigenvalues().asDiagonal());
93 // Note: If MatrixType is real then a.eigenvalues() uses EigenSolver and thus
94 // another algorithm so results may differ slightly
95 verify_is_approx_upto_permutation(a.eigenvalues(), ei1.eigenvalues());
96
97 ComplexEigenSolver<MatrixType> ei2;
98 ei2.setMaxIterations(ComplexSchur<MatrixType>::m_maxIterationsPerRow * rows).compute(a);
99 VERIFY_IS_EQUAL(ei2.info(), Success);
100 VERIFY_IS_EQUAL(ei2.eigenvectors(), ei1.eigenvectors());
101 VERIFY_IS_EQUAL(ei2.eigenvalues(), ei1.eigenvalues());
102 if (rows > 2) {
103 ei2.setMaxIterations(1).compute(a);
104 VERIFY_IS_EQUAL(ei2.info(), NoConvergence);
105 VERIFY_IS_EQUAL(ei2.getMaxIterations(), 1);
106 }
107
108 ComplexEigenSolver<MatrixType> eiNoEivecs(a, false);
109 VERIFY_IS_EQUAL(eiNoEivecs.info(), Success);
110 VERIFY_IS_APPROX(ei1.eigenvalues(), eiNoEivecs.eigenvalues());
111
112 // Regression test for issue #66
113 MatrixType z = MatrixType::Zero(rows,cols);
114 ComplexEigenSolver<MatrixType> eiz(z);
115 VERIFY((eiz.eigenvalues().cwiseEqual(0)).all());
116
117 MatrixType id = MatrixType::Identity(rows, cols);
118 VERIFY_IS_APPROX(id.operatorNorm(), RealScalar(1));
119
Austin Schuh189376f2018-12-20 22:11:15 +1100120 if (rows > 1 && rows < 20)
Brian Silverman72890c22015-09-19 14:37:37 -0400121 {
122 // Test matrix with NaN
123 a(0,0) = std::numeric_limits<typename MatrixType::RealScalar>::quiet_NaN();
124 ComplexEigenSolver<MatrixType> eiNaN(a);
125 VERIFY_IS_EQUAL(eiNaN.info(), NoConvergence);
126 }
Austin Schuh189376f2018-12-20 22:11:15 +1100127
128 // regression test for bug 1098
129 {
130 ComplexEigenSolver<MatrixType> eig(a.adjoint() * a);
131 eig.compute(a.adjoint() * a);
132 }
133
134 // regression test for bug 478
135 {
136 a.setZero();
137 ComplexEigenSolver<MatrixType> ei3(a);
138 VERIFY_IS_EQUAL(ei3.info(), Success);
139 VERIFY_IS_MUCH_SMALLER_THAN(ei3.eigenvalues().norm(),RealScalar(1));
140 VERIFY((ei3.eigenvectors().transpose()*ei3.eigenvectors().transpose()).eval().isIdentity());
141 }
Brian Silverman72890c22015-09-19 14:37:37 -0400142}
143
144template<typename MatrixType> void eigensolver_verify_assert(const MatrixType& m)
145{
146 ComplexEigenSolver<MatrixType> eig;
147 VERIFY_RAISES_ASSERT(eig.eigenvectors());
148 VERIFY_RAISES_ASSERT(eig.eigenvalues());
149
150 MatrixType a = MatrixType::Random(m.rows(),m.cols());
151 eig.compute(a, false);
152 VERIFY_RAISES_ASSERT(eig.eigenvectors());
153}
154
155void test_eigensolver_complex()
156{
157 int s = 0;
158 for(int i = 0; i < g_repeat; i++) {
159 CALL_SUBTEST_1( eigensolver(Matrix4cf()) );
160 s = internal::random<int>(1,EIGEN_TEST_MAX_SIZE/4);
161 CALL_SUBTEST_2( eigensolver(MatrixXcd(s,s)) );
162 CALL_SUBTEST_3( eigensolver(Matrix<std::complex<float>, 1, 1>()) );
163 CALL_SUBTEST_4( eigensolver(Matrix3f()) );
Austin Schuh189376f2018-12-20 22:11:15 +1100164 TEST_SET_BUT_UNUSED_VARIABLE(s)
Brian Silverman72890c22015-09-19 14:37:37 -0400165 }
166 CALL_SUBTEST_1( eigensolver_verify_assert(Matrix4cf()) );
167 s = internal::random<int>(1,EIGEN_TEST_MAX_SIZE/4);
168 CALL_SUBTEST_2( eigensolver_verify_assert(MatrixXcd(s,s)) );
169 CALL_SUBTEST_3( eigensolver_verify_assert(Matrix<std::complex<float>, 1, 1>()) );
170 CALL_SUBTEST_4( eigensolver_verify_assert(Matrix3f()) );
171
172 // Test problem size constructors
173 CALL_SUBTEST_5(ComplexEigenSolver<MatrixXf> tmp(s));
174
175 TEST_SET_BUT_UNUSED_VARIABLE(s)
176}