eigensolver_generic.cpp 9.2 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247
  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,2012 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. template<typename EigType,typename MatType>
  14. void check_eigensolver_for_given_mat(const EigType &eig, const MatType& a)
  15. {
  16. typedef typename NumTraits<typename MatType::Scalar>::Real RealScalar;
  17. typedef Matrix<RealScalar, MatType::RowsAtCompileTime, 1> RealVectorType;
  18. typedef typename std::complex<RealScalar> Complex;
  19. Index n = a.rows();
  20. VERIFY_IS_EQUAL(eig.info(), Success);
  21. VERIFY_IS_APPROX(a * eig.pseudoEigenvectors(), eig.pseudoEigenvectors() * eig.pseudoEigenvalueMatrix());
  22. VERIFY_IS_APPROX(a.template cast<Complex>() * eig.eigenvectors(),
  23. eig.eigenvectors() * eig.eigenvalues().asDiagonal());
  24. VERIFY_IS_APPROX(eig.eigenvectors().colwise().norm(), RealVectorType::Ones(n).transpose());
  25. VERIFY_IS_APPROX(a.eigenvalues(), eig.eigenvalues());
  26. }
  27. template<typename MatrixType> void eigensolver(const MatrixType& m)
  28. {
  29. /* this test covers the following files:
  30. EigenSolver.h
  31. */
  32. Index rows = m.rows();
  33. Index cols = m.cols();
  34. typedef typename MatrixType::Scalar Scalar;
  35. typedef typename NumTraits<Scalar>::Real RealScalar;
  36. typedef typename std::complex<RealScalar> Complex;
  37. MatrixType a = MatrixType::Random(rows,cols);
  38. MatrixType a1 = MatrixType::Random(rows,cols);
  39. MatrixType symmA = a.adjoint() * a + a1.adjoint() * a1;
  40. EigenSolver<MatrixType> ei0(symmA);
  41. VERIFY_IS_EQUAL(ei0.info(), Success);
  42. VERIFY_IS_APPROX(symmA * ei0.pseudoEigenvectors(), ei0.pseudoEigenvectors() * ei0.pseudoEigenvalueMatrix());
  43. VERIFY_IS_APPROX((symmA.template cast<Complex>()) * (ei0.pseudoEigenvectors().template cast<Complex>()),
  44. (ei0.pseudoEigenvectors().template cast<Complex>()) * (ei0.eigenvalues().asDiagonal()));
  45. EigenSolver<MatrixType> ei1(a);
  46. CALL_SUBTEST( check_eigensolver_for_given_mat(ei1,a) );
  47. EigenSolver<MatrixType> ei2;
  48. ei2.setMaxIterations(RealSchur<MatrixType>::m_maxIterationsPerRow * rows).compute(a);
  49. VERIFY_IS_EQUAL(ei2.info(), Success);
  50. VERIFY_IS_EQUAL(ei2.eigenvectors(), ei1.eigenvectors());
  51. VERIFY_IS_EQUAL(ei2.eigenvalues(), ei1.eigenvalues());
  52. if (rows > 2) {
  53. ei2.setMaxIterations(1).compute(a);
  54. VERIFY_IS_EQUAL(ei2.info(), NoConvergence);
  55. VERIFY_IS_EQUAL(ei2.getMaxIterations(), 1);
  56. }
  57. EigenSolver<MatrixType> eiNoEivecs(a, false);
  58. VERIFY_IS_EQUAL(eiNoEivecs.info(), Success);
  59. VERIFY_IS_APPROX(ei1.eigenvalues(), eiNoEivecs.eigenvalues());
  60. VERIFY_IS_APPROX(ei1.pseudoEigenvalueMatrix(), eiNoEivecs.pseudoEigenvalueMatrix());
  61. MatrixType id = MatrixType::Identity(rows, cols);
  62. VERIFY_IS_APPROX(id.operatorNorm(), RealScalar(1));
  63. if (rows > 2 && rows < 20)
  64. {
  65. // Test matrix with NaN
  66. a(0,0) = std::numeric_limits<typename MatrixType::RealScalar>::quiet_NaN();
  67. EigenSolver<MatrixType> eiNaN(a);
  68. VERIFY_IS_NOT_EQUAL(eiNaN.info(), Success);
  69. }
  70. // regression test for bug 1098
  71. {
  72. EigenSolver<MatrixType> eig(a.adjoint() * a);
  73. eig.compute(a.adjoint() * a);
  74. }
  75. // regression test for bug 478
  76. {
  77. a.setZero();
  78. EigenSolver<MatrixType> ei3(a);
  79. VERIFY_IS_EQUAL(ei3.info(), Success);
  80. VERIFY_IS_MUCH_SMALLER_THAN(ei3.eigenvalues().norm(),RealScalar(1));
  81. VERIFY((ei3.eigenvectors().transpose()*ei3.eigenvectors().transpose()).eval().isIdentity());
  82. }
  83. }
  84. template<typename MatrixType> void eigensolver_verify_assert(const MatrixType& m)
  85. {
  86. EigenSolver<MatrixType> eig;
  87. VERIFY_RAISES_ASSERT(eig.eigenvectors());
  88. VERIFY_RAISES_ASSERT(eig.pseudoEigenvectors());
  89. VERIFY_RAISES_ASSERT(eig.pseudoEigenvalueMatrix());
  90. VERIFY_RAISES_ASSERT(eig.eigenvalues());
  91. MatrixType a = MatrixType::Random(m.rows(),m.cols());
  92. eig.compute(a, false);
  93. VERIFY_RAISES_ASSERT(eig.eigenvectors());
  94. VERIFY_RAISES_ASSERT(eig.pseudoEigenvectors());
  95. }
  96. template<typename CoeffType>
  97. Matrix<typename CoeffType::Scalar,Dynamic,Dynamic>
  98. make_companion(const CoeffType& coeffs)
  99. {
  100. Index n = coeffs.size()-1;
  101. Matrix<typename CoeffType::Scalar,Dynamic,Dynamic> res(n,n);
  102. res.setZero();
  103. res.row(0) = -coeffs.tail(n) / coeffs(0);
  104. res.diagonal(-1).setOnes();
  105. return res;
  106. }
  107. template<int>
  108. void eigensolver_generic_extra()
  109. {
  110. {
  111. // regression test for bug 793
  112. MatrixXd a(3,3);
  113. a << 0, 0, 1,
  114. 1, 1, 1,
  115. 1, 1e+200, 1;
  116. Eigen::EigenSolver<MatrixXd> eig(a);
  117. double scale = 1e-200; // scale to avoid overflow during the comparisons
  118. VERIFY_IS_APPROX(a * eig.pseudoEigenvectors()*scale, eig.pseudoEigenvectors() * eig.pseudoEigenvalueMatrix()*scale);
  119. VERIFY_IS_APPROX(a * eig.eigenvectors()*scale, eig.eigenvectors() * eig.eigenvalues().asDiagonal()*scale);
  120. }
  121. {
  122. // check a case where all eigenvalues are null.
  123. MatrixXd a(2,2);
  124. a << 1, 1,
  125. -1, -1;
  126. Eigen::EigenSolver<MatrixXd> eig(a);
  127. VERIFY_IS_APPROX(eig.pseudoEigenvectors().squaredNorm(), 2.);
  128. VERIFY_IS_APPROX((a * eig.pseudoEigenvectors()).norm()+1., 1.);
  129. VERIFY_IS_APPROX((eig.pseudoEigenvectors() * eig.pseudoEigenvalueMatrix()).norm()+1., 1.);
  130. VERIFY_IS_APPROX((a * eig.eigenvectors()).norm()+1., 1.);
  131. VERIFY_IS_APPROX((eig.eigenvectors() * eig.eigenvalues().asDiagonal()).norm()+1., 1.);
  132. }
  133. // regression test for bug 933
  134. {
  135. {
  136. VectorXd coeffs(5); coeffs << 1, -3, -175, -225, 2250;
  137. MatrixXd C = make_companion(coeffs);
  138. EigenSolver<MatrixXd> eig(C);
  139. CALL_SUBTEST( check_eigensolver_for_given_mat(eig,C) );
  140. }
  141. {
  142. // this test is tricky because it requires high accuracy in smallest eigenvalues
  143. VectorXd coeffs(5); coeffs << 6.154671e-15, -1.003870e-10, -9.819570e-01, 3.995715e+03, 2.211511e+08;
  144. MatrixXd C = make_companion(coeffs);
  145. EigenSolver<MatrixXd> eig(C);
  146. CALL_SUBTEST( check_eigensolver_for_given_mat(eig,C) );
  147. Index n = C.rows();
  148. for(Index i=0;i<n;++i)
  149. {
  150. typedef std::complex<double> Complex;
  151. MatrixXcd ac = C.cast<Complex>();
  152. ac.diagonal().array() -= eig.eigenvalues()(i);
  153. VectorXd sv = ac.jacobiSvd().singularValues();
  154. // comparing to sv(0) is not enough here to catch the "bug",
  155. // the hard-coded 1.0 is important!
  156. VERIFY_IS_MUCH_SMALLER_THAN(sv(n-1), 1.0);
  157. }
  158. }
  159. }
  160. // regression test for bug 1557
  161. {
  162. // this test is interesting because it contains zeros on the diagonal.
  163. MatrixXd A_bug1557(3,3);
  164. A_bug1557 << 0, 0, 0, 1, 0, 0.5887907064808635127, 0, 1, 0;
  165. EigenSolver<MatrixXd> eig(A_bug1557);
  166. CALL_SUBTEST( check_eigensolver_for_given_mat(eig,A_bug1557) );
  167. }
  168. // regression test for bug 1174
  169. {
  170. Index n = 12;
  171. MatrixXf A_bug1174(n,n);
  172. A_bug1174 << 262144, 0, 0, 262144, 786432, 0, 0, 0, 0, 0, 0, 786432,
  173. 262144, 0, 0, 262144, 786432, 0, 0, 0, 0, 0, 0, 786432,
  174. 262144, 0, 0, 262144, 786432, 0, 0, 0, 0, 0, 0, 786432,
  175. 262144, 0, 0, 262144, 786432, 0, 0, 0, 0, 0, 0, 786432,
  176. 0, 262144, 262144, 0, 0, 262144, 262144, 262144, 262144, 262144, 262144, 0,
  177. 0, 262144, 262144, 0, 0, 262144, 262144, 262144, 262144, 262144, 262144, 0,
  178. 0, 262144, 262144, 0, 0, 262144, 262144, 262144, 262144, 262144, 262144, 0,
  179. 0, 262144, 262144, 0, 0, 262144, 262144, 262144, 262144, 262144, 262144, 0,
  180. 0, 262144, 262144, 0, 0, 262144, 262144, 262144, 262144, 262144, 262144, 0,
  181. 0, 262144, 262144, 0, 0, 262144, 262144, 262144, 262144, 262144, 262144, 0,
  182. 0, 262144, 262144, 0, 0, 262144, 262144, 262144, 262144, 262144, 262144, 0,
  183. 0, 262144, 262144, 0, 0, 262144, 262144, 262144, 262144, 262144, 262144, 0;
  184. EigenSolver<MatrixXf> eig(A_bug1174);
  185. CALL_SUBTEST( check_eigensolver_for_given_mat(eig,A_bug1174) );
  186. }
  187. }
  188. EIGEN_DECLARE_TEST(eigensolver_generic)
  189. {
  190. int s = 0;
  191. for(int i = 0; i < g_repeat; i++) {
  192. CALL_SUBTEST_1( eigensolver(Matrix4f()) );
  193. s = internal::random<int>(1,EIGEN_TEST_MAX_SIZE/4);
  194. CALL_SUBTEST_2( eigensolver(MatrixXd(s,s)) );
  195. TEST_SET_BUT_UNUSED_VARIABLE(s)
  196. // some trivial but implementation-wise tricky cases
  197. CALL_SUBTEST_2( eigensolver(MatrixXd(1,1)) );
  198. CALL_SUBTEST_2( eigensolver(MatrixXd(2,2)) );
  199. CALL_SUBTEST_3( eigensolver(Matrix<double,1,1>()) );
  200. CALL_SUBTEST_4( eigensolver(Matrix2d()) );
  201. }
  202. CALL_SUBTEST_1( eigensolver_verify_assert(Matrix4f()) );
  203. s = internal::random<int>(1,EIGEN_TEST_MAX_SIZE/4);
  204. CALL_SUBTEST_2( eigensolver_verify_assert(MatrixXd(s,s)) );
  205. CALL_SUBTEST_3( eigensolver_verify_assert(Matrix<double,1,1>()) );
  206. CALL_SUBTEST_4( eigensolver_verify_assert(Matrix2d()) );
  207. // Test problem size constructors
  208. CALL_SUBTEST_5(EigenSolver<MatrixXf> tmp(s));
  209. // regression test for bug 410
  210. CALL_SUBTEST_2(
  211. {
  212. MatrixXd A(1,1);
  213. A(0,0) = std::sqrt(-1.); // is Not-a-Number
  214. Eigen::EigenSolver<MatrixXd> solver(A);
  215. VERIFY_IS_EQUAL(solver.info(), NumericalIssue);
  216. }
  217. );
  218. CALL_SUBTEST_2( eigensolver_generic_extra<0>() );
  219. TEST_SET_BUT_UNUSED_VARIABLE(s)
  220. }