You can not select more than 25 topics Topics must start with a letter or number, can include dashes ('-') and can be up to 35 characters long.

223 lines
9.5 KiB

  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 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 "svd_fill.h"
  12. #include <limits>
  13. #include <Eigen/Eigenvalues>
  14. template<typename MatrixType> void selfadjointeigensolver_essential_check(const MatrixType& m)
  15. {
  16. typedef typename MatrixType::Scalar Scalar;
  17. typedef typename NumTraits<Scalar>::Real RealScalar;
  18. RealScalar eival_eps = (std::min)(test_precision<RealScalar>(), NumTraits<Scalar>::dummy_precision()*20000);
  19. SelfAdjointEigenSolver<MatrixType> eiSymm(m);
  20. VERIFY_IS_EQUAL(eiSymm.info(), Success);
  21. VERIFY_IS_APPROX(m.template selfadjointView<Lower>() * eiSymm.eigenvectors(),
  22. eiSymm.eigenvectors() * eiSymm.eigenvalues().asDiagonal());
  23. VERIFY_IS_APPROX(m.template selfadjointView<Lower>().eigenvalues(), eiSymm.eigenvalues());
  24. VERIFY_IS_UNITARY(eiSymm.eigenvectors());
  25. if(m.cols()<=4)
  26. {
  27. SelfAdjointEigenSolver<MatrixType> eiDirect;
  28. eiDirect.computeDirect(m);
  29. VERIFY_IS_EQUAL(eiDirect.info(), Success);
  30. VERIFY_IS_APPROX(eiSymm.eigenvalues(), eiDirect.eigenvalues());
  31. if(! eiSymm.eigenvalues().isApprox(eiDirect.eigenvalues(), eival_eps) )
  32. {
  33. std::cerr << "reference eigenvalues: " << eiSymm.eigenvalues().transpose() << "\n"
  34. << "obtained eigenvalues: " << eiDirect.eigenvalues().transpose() << "\n"
  35. << "diff: " << (eiSymm.eigenvalues()-eiDirect.eigenvalues()).transpose() << "\n"
  36. << "error (eps): " << (eiSymm.eigenvalues()-eiDirect.eigenvalues()).norm() / eiSymm.eigenvalues().norm() << " (" << eival_eps << ")\n";
  37. }
  38. VERIFY(eiSymm.eigenvalues().isApprox(eiDirect.eigenvalues(), eival_eps));
  39. VERIFY_IS_APPROX(m.template selfadjointView<Lower>() * eiDirect.eigenvectors(),
  40. eiDirect.eigenvectors() * eiDirect.eigenvalues().asDiagonal());
  41. VERIFY_IS_APPROX(m.template selfadjointView<Lower>().eigenvalues(), eiDirect.eigenvalues());
  42. VERIFY_IS_UNITARY(eiDirect.eigenvectors());
  43. }
  44. }
  45. template<typename MatrixType> void selfadjointeigensolver(const MatrixType& m)
  46. {
  47. typedef typename MatrixType::Index Index;
  48. /* this test covers the following files:
  49. EigenSolver.h, SelfAdjointEigenSolver.h (and indirectly: Tridiagonalization.h)
  50. */
  51. Index rows = m.rows();
  52. Index cols = m.cols();
  53. typedef typename MatrixType::Scalar Scalar;
  54. typedef typename NumTraits<Scalar>::Real RealScalar;
  55. RealScalar largerEps = 10*test_precision<RealScalar>();
  56. MatrixType a = MatrixType::Random(rows,cols);
  57. MatrixType a1 = MatrixType::Random(rows,cols);
  58. MatrixType symmA = a.adjoint() * a + a1.adjoint() * a1;
  59. MatrixType symmC = symmA;
  60. svd_fill_random(symmA,Symmetric);
  61. symmA.template triangularView<StrictlyUpper>().setZero();
  62. symmC.template triangularView<StrictlyUpper>().setZero();
  63. MatrixType b = MatrixType::Random(rows,cols);
  64. MatrixType b1 = MatrixType::Random(rows,cols);
  65. MatrixType symmB = b.adjoint() * b + b1.adjoint() * b1;
  66. symmB.template triangularView<StrictlyUpper>().setZero();
  67. CALL_SUBTEST( selfadjointeigensolver_essential_check(symmA) );
  68. SelfAdjointEigenSolver<MatrixType> eiSymm(symmA);
  69. // generalized eigen pb
  70. GeneralizedSelfAdjointEigenSolver<MatrixType> eiSymmGen(symmC, symmB);
  71. SelfAdjointEigenSolver<MatrixType> eiSymmNoEivecs(symmA, false);
  72. VERIFY_IS_EQUAL(eiSymmNoEivecs.info(), Success);
  73. VERIFY_IS_APPROX(eiSymm.eigenvalues(), eiSymmNoEivecs.eigenvalues());
  74. // generalized eigen problem Ax = lBx
  75. eiSymmGen.compute(symmC, symmB,Ax_lBx);
  76. VERIFY_IS_EQUAL(eiSymmGen.info(), Success);
  77. VERIFY((symmC.template selfadjointView<Lower>() * eiSymmGen.eigenvectors()).isApprox(
  78. symmB.template selfadjointView<Lower>() * (eiSymmGen.eigenvectors() * eiSymmGen.eigenvalues().asDiagonal()), largerEps));
  79. // generalized eigen problem BAx = lx
  80. eiSymmGen.compute(symmC, symmB,BAx_lx);
  81. VERIFY_IS_EQUAL(eiSymmGen.info(), Success);
  82. VERIFY((symmB.template selfadjointView<Lower>() * (symmC.template selfadjointView<Lower>() * eiSymmGen.eigenvectors())).isApprox(
  83. (eiSymmGen.eigenvectors() * eiSymmGen.eigenvalues().asDiagonal()), largerEps));
  84. // generalized eigen problem ABx = lx
  85. eiSymmGen.compute(symmC, symmB,ABx_lx);
  86. VERIFY_IS_EQUAL(eiSymmGen.info(), Success);
  87. VERIFY((symmC.template selfadjointView<Lower>() * (symmB.template selfadjointView<Lower>() * eiSymmGen.eigenvectors())).isApprox(
  88. (eiSymmGen.eigenvectors() * eiSymmGen.eigenvalues().asDiagonal()), largerEps));
  89. eiSymm.compute(symmC);
  90. MatrixType sqrtSymmA = eiSymm.operatorSqrt();
  91. VERIFY_IS_APPROX(MatrixType(symmC.template selfadjointView<Lower>()), sqrtSymmA*sqrtSymmA);
  92. VERIFY_IS_APPROX(sqrtSymmA, symmC.template selfadjointView<Lower>()*eiSymm.operatorInverseSqrt());
  93. MatrixType id = MatrixType::Identity(rows, cols);
  94. VERIFY_IS_APPROX(id.template selfadjointView<Lower>().operatorNorm(), RealScalar(1));
  95. SelfAdjointEigenSolver<MatrixType> eiSymmUninitialized;
  96. VERIFY_RAISES_ASSERT(eiSymmUninitialized.info());
  97. VERIFY_RAISES_ASSERT(eiSymmUninitialized.eigenvalues());
  98. VERIFY_RAISES_ASSERT(eiSymmUninitialized.eigenvectors());
  99. VERIFY_RAISES_ASSERT(eiSymmUninitialized.operatorSqrt());
  100. VERIFY_RAISES_ASSERT(eiSymmUninitialized.operatorInverseSqrt());
  101. eiSymmUninitialized.compute(symmA, false);
  102. VERIFY_RAISES_ASSERT(eiSymmUninitialized.eigenvectors());
  103. VERIFY_RAISES_ASSERT(eiSymmUninitialized.operatorSqrt());
  104. VERIFY_RAISES_ASSERT(eiSymmUninitialized.operatorInverseSqrt());
  105. // test Tridiagonalization's methods
  106. Tridiagonalization<MatrixType> tridiag(symmC);
  107. VERIFY_IS_APPROX(tridiag.diagonal(), tridiag.matrixT().diagonal());
  108. VERIFY_IS_APPROX(tridiag.subDiagonal(), tridiag.matrixT().template diagonal<-1>());
  109. Matrix<RealScalar,Dynamic,Dynamic> T = tridiag.matrixT();
  110. if(rows>1 && cols>1) {
  111. // FIXME check that upper and lower part are 0:
  112. //VERIFY(T.topRightCorner(rows-2, cols-2).template triangularView<Upper>().isZero());
  113. }
  114. VERIFY_IS_APPROX(tridiag.diagonal(), T.diagonal());
  115. VERIFY_IS_APPROX(tridiag.subDiagonal(), T.template diagonal<1>());
  116. VERIFY_IS_APPROX(MatrixType(symmC.template selfadjointView<Lower>()), tridiag.matrixQ() * tridiag.matrixT().eval() * MatrixType(tridiag.matrixQ()).adjoint());
  117. VERIFY_IS_APPROX(MatrixType(symmC.template selfadjointView<Lower>()), tridiag.matrixQ() * tridiag.matrixT() * tridiag.matrixQ().adjoint());
  118. // Test computation of eigenvalues from tridiagonal matrix
  119. if(rows > 1)
  120. {
  121. SelfAdjointEigenSolver<MatrixType> eiSymmTridiag;
  122. eiSymmTridiag.computeFromTridiagonal(tridiag.matrixT().diagonal(), tridiag.matrixT().diagonal(-1), ComputeEigenvectors);
  123. VERIFY_IS_APPROX(eiSymm.eigenvalues(), eiSymmTridiag.eigenvalues());
  124. VERIFY_IS_APPROX(tridiag.matrixT(), eiSymmTridiag.eigenvectors().real() * eiSymmTridiag.eigenvalues().asDiagonal() * eiSymmTridiag.eigenvectors().real().transpose());
  125. }
  126. if (rows > 1 && rows < 20)
  127. {
  128. // Test matrix with NaN
  129. symmC(0,0) = std::numeric_limits<typename MatrixType::RealScalar>::quiet_NaN();
  130. SelfAdjointEigenSolver<MatrixType> eiSymmNaN(symmC);
  131. VERIFY_IS_EQUAL(eiSymmNaN.info(), NoConvergence);
  132. }
  133. // regression test for bug 1098
  134. {
  135. SelfAdjointEigenSolver<MatrixType> eig(a.adjoint() * a);
  136. eig.compute(a.adjoint() * a);
  137. }
  138. }
  139. void bug_854()
  140. {
  141. Matrix3d m;
  142. m << 850.961, 51.966, 0,
  143. 51.966, 254.841, 0,
  144. 0, 0, 0;
  145. selfadjointeigensolver_essential_check(m);
  146. }
  147. void bug_1014()
  148. {
  149. Matrix3d m;
  150. m << 0.11111111111111114658, 0, 0,
  151. 0, 0.11111111111111109107, 0,
  152. 0, 0, 0.11111111111111107719;
  153. selfadjointeigensolver_essential_check(m);
  154. }
  155. void test_eigensolver_selfadjoint()
  156. {
  157. int s = 0;
  158. for(int i = 0; i < g_repeat; i++) {
  159. // trivial test for 1x1 matrices:
  160. CALL_SUBTEST_1( selfadjointeigensolver(Matrix<float, 1, 1>()));
  161. CALL_SUBTEST_1( selfadjointeigensolver(Matrix<double, 1, 1>()));
  162. // very important to test 3x3 and 2x2 matrices since we provide special paths for them
  163. CALL_SUBTEST_12( selfadjointeigensolver(Matrix2f()) );
  164. CALL_SUBTEST_12( selfadjointeigensolver(Matrix2d()) );
  165. CALL_SUBTEST_13( selfadjointeigensolver(Matrix3f()) );
  166. CALL_SUBTEST_13( selfadjointeigensolver(Matrix3d()) );
  167. CALL_SUBTEST_2( selfadjointeigensolver(Matrix4d()) );
  168. s = internal::random<int>(1,EIGEN_TEST_MAX_SIZE/4);
  169. CALL_SUBTEST_3( selfadjointeigensolver(MatrixXf(s,s)) );
  170. CALL_SUBTEST_4( selfadjointeigensolver(MatrixXd(s,s)) );
  171. CALL_SUBTEST_5( selfadjointeigensolver(MatrixXcd(s,s)) );
  172. CALL_SUBTEST_9( selfadjointeigensolver(Matrix<std::complex<double>,Dynamic,Dynamic,RowMajor>(s,s)) );
  173. TEST_SET_BUT_UNUSED_VARIABLE(s)
  174. // some trivial but implementation-wise tricky cases
  175. CALL_SUBTEST_4( selfadjointeigensolver(MatrixXd(1,1)) );
  176. CALL_SUBTEST_4( selfadjointeigensolver(MatrixXd(2,2)) );
  177. CALL_SUBTEST_6( selfadjointeigensolver(Matrix<double,1,1>()) );
  178. CALL_SUBTEST_7( selfadjointeigensolver(Matrix<double,2,2>()) );
  179. }
  180. CALL_SUBTEST_13( bug_854() );
  181. CALL_SUBTEST_13( bug_1014() );
  182. // Test problem size constructors
  183. s = internal::random<int>(1,EIGEN_TEST_MAX_SIZE/4);
  184. CALL_SUBTEST_8(SelfAdjointEigenSolver<MatrixXf> tmp1(s));
  185. CALL_SUBTEST_8(Tridiagonalization<MatrixXf> tmp2(s));
  186. TEST_SET_BUT_UNUSED_VARIABLE(s)
  187. }