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.

143 lines
6.3 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 <limits>
  12. #include <Eigen/Eigenvalues>
  13. template<typename MatrixType> void selfadjointeigensolver(const MatrixType& m)
  14. {
  15. typedef typename MatrixType::Index Index;
  16. /* this test covers the following files:
  17. EigenSolver.h, SelfAdjointEigenSolver.h (and indirectly: Tridiagonalization.h)
  18. */
  19. Index rows = m.rows();
  20. Index cols = m.cols();
  21. typedef typename MatrixType::Scalar Scalar;
  22. typedef typename NumTraits<Scalar>::Real RealScalar;
  23. RealScalar largerEps = 10*test_precision<RealScalar>();
  24. MatrixType a = MatrixType::Random(rows,cols);
  25. MatrixType a1 = MatrixType::Random(rows,cols);
  26. MatrixType symmA = a.adjoint() * a + a1.adjoint() * a1;
  27. symmA.template triangularView<StrictlyUpper>().setZero();
  28. MatrixType b = MatrixType::Random(rows,cols);
  29. MatrixType b1 = MatrixType::Random(rows,cols);
  30. MatrixType symmB = b.adjoint() * b + b1.adjoint() * b1;
  31. symmB.template triangularView<StrictlyUpper>().setZero();
  32. SelfAdjointEigenSolver<MatrixType> eiSymm(symmA);
  33. SelfAdjointEigenSolver<MatrixType> eiDirect;
  34. eiDirect.computeDirect(symmA);
  35. // generalized eigen pb
  36. GeneralizedSelfAdjointEigenSolver<MatrixType> eiSymmGen(symmA, symmB);
  37. VERIFY_IS_EQUAL(eiSymm.info(), Success);
  38. VERIFY((symmA.template selfadjointView<Lower>() * eiSymm.eigenvectors()).isApprox(
  39. eiSymm.eigenvectors() * eiSymm.eigenvalues().asDiagonal(), largerEps));
  40. VERIFY_IS_APPROX(symmA.template selfadjointView<Lower>().eigenvalues(), eiSymm.eigenvalues());
  41. VERIFY_IS_EQUAL(eiDirect.info(), Success);
  42. VERIFY((symmA.template selfadjointView<Lower>() * eiDirect.eigenvectors()).isApprox(
  43. eiDirect.eigenvectors() * eiDirect.eigenvalues().asDiagonal(), largerEps));
  44. VERIFY_IS_APPROX(symmA.template selfadjointView<Lower>().eigenvalues(), eiDirect.eigenvalues());
  45. SelfAdjointEigenSolver<MatrixType> eiSymmNoEivecs(symmA, false);
  46. VERIFY_IS_EQUAL(eiSymmNoEivecs.info(), Success);
  47. VERIFY_IS_APPROX(eiSymm.eigenvalues(), eiSymmNoEivecs.eigenvalues());
  48. // generalized eigen problem Ax = lBx
  49. eiSymmGen.compute(symmA, symmB,Ax_lBx);
  50. VERIFY_IS_EQUAL(eiSymmGen.info(), Success);
  51. VERIFY((symmA.template selfadjointView<Lower>() * eiSymmGen.eigenvectors()).isApprox(
  52. symmB.template selfadjointView<Lower>() * (eiSymmGen.eigenvectors() * eiSymmGen.eigenvalues().asDiagonal()), largerEps));
  53. // generalized eigen problem BAx = lx
  54. eiSymmGen.compute(symmA, symmB,BAx_lx);
  55. VERIFY_IS_EQUAL(eiSymmGen.info(), Success);
  56. VERIFY((symmB.template selfadjointView<Lower>() * (symmA.template selfadjointView<Lower>() * eiSymmGen.eigenvectors())).isApprox(
  57. (eiSymmGen.eigenvectors() * eiSymmGen.eigenvalues().asDiagonal()), largerEps));
  58. // generalized eigen problem ABx = lx
  59. eiSymmGen.compute(symmA, symmB,ABx_lx);
  60. VERIFY_IS_EQUAL(eiSymmGen.info(), Success);
  61. VERIFY((symmA.template selfadjointView<Lower>() * (symmB.template selfadjointView<Lower>() * eiSymmGen.eigenvectors())).isApprox(
  62. (eiSymmGen.eigenvectors() * eiSymmGen.eigenvalues().asDiagonal()), largerEps));
  63. MatrixType sqrtSymmA = eiSymm.operatorSqrt();
  64. VERIFY_IS_APPROX(MatrixType(symmA.template selfadjointView<Lower>()), sqrtSymmA*sqrtSymmA);
  65. VERIFY_IS_APPROX(sqrtSymmA, symmA.template selfadjointView<Lower>()*eiSymm.operatorInverseSqrt());
  66. MatrixType id = MatrixType::Identity(rows, cols);
  67. VERIFY_IS_APPROX(id.template selfadjointView<Lower>().operatorNorm(), RealScalar(1));
  68. SelfAdjointEigenSolver<MatrixType> eiSymmUninitialized;
  69. VERIFY_RAISES_ASSERT(eiSymmUninitialized.info());
  70. VERIFY_RAISES_ASSERT(eiSymmUninitialized.eigenvalues());
  71. VERIFY_RAISES_ASSERT(eiSymmUninitialized.eigenvectors());
  72. VERIFY_RAISES_ASSERT(eiSymmUninitialized.operatorSqrt());
  73. VERIFY_RAISES_ASSERT(eiSymmUninitialized.operatorInverseSqrt());
  74. eiSymmUninitialized.compute(symmA, false);
  75. VERIFY_RAISES_ASSERT(eiSymmUninitialized.eigenvectors());
  76. VERIFY_RAISES_ASSERT(eiSymmUninitialized.operatorSqrt());
  77. VERIFY_RAISES_ASSERT(eiSymmUninitialized.operatorInverseSqrt());
  78. // test Tridiagonalization's methods
  79. Tridiagonalization<MatrixType> tridiag(symmA);
  80. // FIXME tridiag.matrixQ().adjoint() does not work
  81. VERIFY_IS_APPROX(MatrixType(symmA.template selfadjointView<Lower>()), tridiag.matrixQ() * tridiag.matrixT().eval() * MatrixType(tridiag.matrixQ()).adjoint());
  82. if (rows > 1)
  83. {
  84. // Test matrix with NaN
  85. symmA(0,0) = std::numeric_limits<typename MatrixType::RealScalar>::quiet_NaN();
  86. SelfAdjointEigenSolver<MatrixType> eiSymmNaN(symmA);
  87. VERIFY_IS_EQUAL(eiSymmNaN.info(), NoConvergence);
  88. }
  89. }
  90. void test_eigensolver_selfadjoint()
  91. {
  92. int s = 0;
  93. for(int i = 0; i < g_repeat; i++) {
  94. // very important to test 3x3 and 2x2 matrices since we provide special paths for them
  95. CALL_SUBTEST_1( selfadjointeigensolver(Matrix2d()) );
  96. CALL_SUBTEST_1( selfadjointeigensolver(Matrix3f()) );
  97. CALL_SUBTEST_2( selfadjointeigensolver(Matrix4d()) );
  98. s = internal::random<int>(1,EIGEN_TEST_MAX_SIZE/4);
  99. CALL_SUBTEST_3( selfadjointeigensolver(MatrixXf(s,s)) );
  100. s = internal::random<int>(1,EIGEN_TEST_MAX_SIZE/4);
  101. CALL_SUBTEST_4( selfadjointeigensolver(MatrixXd(s,s)) );
  102. s = internal::random<int>(1,EIGEN_TEST_MAX_SIZE/4);
  103. CALL_SUBTEST_5( selfadjointeigensolver(MatrixXcd(s,s)) );
  104. s = internal::random<int>(1,EIGEN_TEST_MAX_SIZE/4);
  105. CALL_SUBTEST_9( selfadjointeigensolver(Matrix<std::complex<double>,Dynamic,Dynamic,RowMajor>(s,s)) );
  106. // some trivial but implementation-wise tricky cases
  107. CALL_SUBTEST_4( selfadjointeigensolver(MatrixXd(1,1)) );
  108. CALL_SUBTEST_4( selfadjointeigensolver(MatrixXd(2,2)) );
  109. CALL_SUBTEST_6( selfadjointeigensolver(Matrix<double,1,1>()) );
  110. CALL_SUBTEST_7( selfadjointeigensolver(Matrix<double,2,2>()) );
  111. }
  112. // Test problem size constructors
  113. s = internal::random<int>(1,EIGEN_TEST_MAX_SIZE/4);
  114. CALL_SUBTEST_8(SelfAdjointEigenSolver<MatrixXf> tmp1(s));
  115. CALL_SUBTEST_8(Tridiagonalization<MatrixXf> tmp2(s));
  116. TEST_SET_BUT_UNUSED_VARIABLE(s)
  117. }