[136] | 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 |
|
---|
| 16 | /* Check that two column vectors are approximately equal upto permutations,
|
---|
| 17 | by checking that the k-th power sums are equal for k = 1, ..., vec1.rows() */
|
---|
| 18 | template<typename VectorType>
|
---|
| 19 | void verify_is_approx_upto_permutation(const VectorType& vec1, const VectorType& vec2)
|
---|
| 20 | {
|
---|
| 21 | typedef typename NumTraits<typename VectorType::Scalar>::Real RealScalar;
|
---|
| 22 |
|
---|
| 23 | VERIFY(vec1.cols() == 1);
|
---|
| 24 | VERIFY(vec2.cols() == 1);
|
---|
| 25 | VERIFY(vec1.rows() == vec2.rows());
|
---|
| 26 | for (int k = 1; k <= vec1.rows(); ++k)
|
---|
| 27 | {
|
---|
| 28 | VERIFY_IS_APPROX(vec1.array().pow(RealScalar(k)).sum(), vec2.array().pow(RealScalar(k)).sum());
|
---|
| 29 | }
|
---|
| 30 | }
|
---|
| 31 |
|
---|
| 32 |
|
---|
| 33 | template<typename MatrixType> void eigensolver(const MatrixType& m)
|
---|
| 34 | {
|
---|
| 35 | typedef typename MatrixType::Index Index;
|
---|
| 36 | /* this test covers the following files:
|
---|
| 37 | ComplexEigenSolver.h, and indirectly ComplexSchur.h
|
---|
| 38 | */
|
---|
| 39 | Index rows = m.rows();
|
---|
| 40 | Index cols = m.cols();
|
---|
| 41 |
|
---|
| 42 | typedef typename MatrixType::Scalar Scalar;
|
---|
| 43 | typedef typename NumTraits<Scalar>::Real RealScalar;
|
---|
| 44 |
|
---|
| 45 | MatrixType a = MatrixType::Random(rows,cols);
|
---|
| 46 | MatrixType symmA = a.adjoint() * a;
|
---|
| 47 |
|
---|
| 48 | ComplexEigenSolver<MatrixType> ei0(symmA);
|
---|
| 49 | VERIFY_IS_EQUAL(ei0.info(), Success);
|
---|
| 50 | VERIFY_IS_APPROX(symmA * ei0.eigenvectors(), ei0.eigenvectors() * ei0.eigenvalues().asDiagonal());
|
---|
| 51 |
|
---|
| 52 | ComplexEigenSolver<MatrixType> ei1(a);
|
---|
| 53 | VERIFY_IS_EQUAL(ei1.info(), Success);
|
---|
| 54 | VERIFY_IS_APPROX(a * ei1.eigenvectors(), ei1.eigenvectors() * ei1.eigenvalues().asDiagonal());
|
---|
| 55 | // Note: If MatrixType is real then a.eigenvalues() uses EigenSolver and thus
|
---|
| 56 | // another algorithm so results may differ slightly
|
---|
| 57 | verify_is_approx_upto_permutation(a.eigenvalues(), ei1.eigenvalues());
|
---|
| 58 |
|
---|
| 59 | ComplexEigenSolver<MatrixType> ei2;
|
---|
| 60 | ei2.setMaxIterations(ComplexSchur<MatrixType>::m_maxIterationsPerRow * rows).compute(a);
|
---|
| 61 | VERIFY_IS_EQUAL(ei2.info(), Success);
|
---|
| 62 | VERIFY_IS_EQUAL(ei2.eigenvectors(), ei1.eigenvectors());
|
---|
| 63 | VERIFY_IS_EQUAL(ei2.eigenvalues(), ei1.eigenvalues());
|
---|
| 64 | if (rows > 2) {
|
---|
| 65 | ei2.setMaxIterations(1).compute(a);
|
---|
| 66 | VERIFY_IS_EQUAL(ei2.info(), NoConvergence);
|
---|
| 67 | VERIFY_IS_EQUAL(ei2.getMaxIterations(), 1);
|
---|
| 68 | }
|
---|
| 69 |
|
---|
| 70 | ComplexEigenSolver<MatrixType> eiNoEivecs(a, false);
|
---|
| 71 | VERIFY_IS_EQUAL(eiNoEivecs.info(), Success);
|
---|
| 72 | VERIFY_IS_APPROX(ei1.eigenvalues(), eiNoEivecs.eigenvalues());
|
---|
| 73 |
|
---|
| 74 | // Regression test for issue #66
|
---|
| 75 | MatrixType z = MatrixType::Zero(rows,cols);
|
---|
| 76 | ComplexEigenSolver<MatrixType> eiz(z);
|
---|
| 77 | VERIFY((eiz.eigenvalues().cwiseEqual(0)).all());
|
---|
| 78 |
|
---|
| 79 | MatrixType id = MatrixType::Identity(rows, cols);
|
---|
| 80 | VERIFY_IS_APPROX(id.operatorNorm(), RealScalar(1));
|
---|
| 81 |
|
---|
| 82 | if (rows > 1)
|
---|
| 83 | {
|
---|
| 84 | // Test matrix with NaN
|
---|
| 85 | a(0,0) = std::numeric_limits<typename MatrixType::RealScalar>::quiet_NaN();
|
---|
| 86 | ComplexEigenSolver<MatrixType> eiNaN(a);
|
---|
| 87 | VERIFY_IS_EQUAL(eiNaN.info(), NoConvergence);
|
---|
| 88 | }
|
---|
| 89 | }
|
---|
| 90 |
|
---|
| 91 | template<typename MatrixType> void eigensolver_verify_assert(const MatrixType& m)
|
---|
| 92 | {
|
---|
| 93 | ComplexEigenSolver<MatrixType> eig;
|
---|
| 94 | VERIFY_RAISES_ASSERT(eig.eigenvectors());
|
---|
| 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 | }
|
---|
| 101 |
|
---|
| 102 | void test_eigensolver_complex()
|
---|
| 103 | {
|
---|
| 104 | int s = 0;
|
---|
| 105 | for(int i = 0; i < g_repeat; i++) {
|
---|
| 106 | CALL_SUBTEST_1( eigensolver(Matrix4cf()) );
|
---|
| 107 | s = internal::random<int>(1,EIGEN_TEST_MAX_SIZE/4);
|
---|
| 108 | CALL_SUBTEST_2( eigensolver(MatrixXcd(s,s)) );
|
---|
| 109 | CALL_SUBTEST_3( eigensolver(Matrix<std::complex<float>, 1, 1>()) );
|
---|
| 110 | CALL_SUBTEST_4( eigensolver(Matrix3f()) );
|
---|
| 111 | }
|
---|
| 112 | CALL_SUBTEST_1( eigensolver_verify_assert(Matrix4cf()) );
|
---|
| 113 | s = internal::random<int>(1,EIGEN_TEST_MAX_SIZE/4);
|
---|
| 114 | CALL_SUBTEST_2( eigensolver_verify_assert(MatrixXcd(s,s)) );
|
---|
| 115 | CALL_SUBTEST_3( eigensolver_verify_assert(Matrix<std::complex<float>, 1, 1>()) );
|
---|
| 116 | CALL_SUBTEST_4( eigensolver_verify_assert(Matrix3f()) );
|
---|
| 117 |
|
---|
| 118 | // Test problem size constructors
|
---|
| 119 | CALL_SUBTEST_5(ComplexEigenSolver<MatrixXf> tmp(s));
|
---|
| 120 |
|
---|
| 121 | TEST_SET_BUT_UNUSED_VARIABLE(s)
|
---|
| 122 | }
|
---|