eigensolver_complex.cpp 6.1 KB

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