eigensolver_generalized_real.cpp 4.0 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103
  1. // This file is part of Eigen, a lightweight C++ template library
  2. // for linear algebra.
  3. //
  4. // Copyright (C) 2012-2016 Gael Guennebaud <gael.guennebaud@inria.fr>
  5. //
  6. // This Source Code Form is subject to the terms of the Mozilla
  7. // Public License v. 2.0. If a copy of the MPL was not distributed
  8. // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
  9. #define EIGEN_RUNTIME_NO_MALLOC
  10. #include "main.h"
  11. #include <limits>
  12. #include <Eigen/Eigenvalues>
  13. #include <Eigen/LU>
  14. template<typename MatrixType> void generalized_eigensolver_real(const MatrixType& m)
  15. {
  16. /* this test covers the following files:
  17. GeneralizedEigenSolver.h
  18. */
  19. Index rows = m.rows();
  20. Index cols = m.cols();
  21. typedef typename MatrixType::Scalar Scalar;
  22. typedef std::complex<Scalar> ComplexScalar;
  23. typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType;
  24. MatrixType a = MatrixType::Random(rows,cols);
  25. MatrixType b = MatrixType::Random(rows,cols);
  26. MatrixType a1 = MatrixType::Random(rows,cols);
  27. MatrixType b1 = MatrixType::Random(rows,cols);
  28. MatrixType spdA = a.adjoint() * a + a1.adjoint() * a1;
  29. MatrixType spdB = b.adjoint() * b + b1.adjoint() * b1;
  30. // lets compare to GeneralizedSelfAdjointEigenSolver
  31. {
  32. GeneralizedSelfAdjointEigenSolver<MatrixType> symmEig(spdA, spdB);
  33. GeneralizedEigenSolver<MatrixType> eig(spdA, spdB);
  34. VERIFY_IS_EQUAL(eig.eigenvalues().imag().cwiseAbs().maxCoeff(), 0);
  35. VectorType realEigenvalues = eig.eigenvalues().real();
  36. std::sort(realEigenvalues.data(), realEigenvalues.data()+realEigenvalues.size());
  37. VERIFY_IS_APPROX(realEigenvalues, symmEig.eigenvalues());
  38. // check eigenvectors
  39. typename GeneralizedEigenSolver<MatrixType>::EigenvectorsType D = eig.eigenvalues().asDiagonal();
  40. typename GeneralizedEigenSolver<MatrixType>::EigenvectorsType V = eig.eigenvectors();
  41. VERIFY_IS_APPROX(spdA*V, spdB*V*D);
  42. }
  43. // non symmetric case:
  44. {
  45. GeneralizedEigenSolver<MatrixType> eig(rows);
  46. // TODO enable full-prealocation of required memory, this probably requires an in-place mode for HessenbergDecomposition
  47. //Eigen::internal::set_is_malloc_allowed(false);
  48. eig.compute(a,b);
  49. //Eigen::internal::set_is_malloc_allowed(true);
  50. for(Index k=0; k<cols; ++k)
  51. {
  52. Matrix<ComplexScalar,Dynamic,Dynamic> tmp = (eig.betas()(k)*a).template cast<ComplexScalar>() - eig.alphas()(k)*b;
  53. if(tmp.size()>1 && tmp.norm()>(std::numeric_limits<Scalar>::min)())
  54. tmp /= tmp.norm();
  55. VERIFY_IS_MUCH_SMALLER_THAN( std::abs(tmp.determinant()), Scalar(1) );
  56. }
  57. // check eigenvectors
  58. typename GeneralizedEigenSolver<MatrixType>::EigenvectorsType D = eig.eigenvalues().asDiagonal();
  59. typename GeneralizedEigenSolver<MatrixType>::EigenvectorsType V = eig.eigenvectors();
  60. VERIFY_IS_APPROX(a*V, b*V*D);
  61. }
  62. // regression test for bug 1098
  63. {
  64. GeneralizedSelfAdjointEigenSolver<MatrixType> eig1(a.adjoint() * a,b.adjoint() * b);
  65. eig1.compute(a.adjoint() * a,b.adjoint() * b);
  66. GeneralizedEigenSolver<MatrixType> eig2(a.adjoint() * a,b.adjoint() * b);
  67. eig2.compute(a.adjoint() * a,b.adjoint() * b);
  68. }
  69. // check without eigenvectors
  70. {
  71. GeneralizedEigenSolver<MatrixType> eig1(spdA, spdB, true);
  72. GeneralizedEigenSolver<MatrixType> eig2(spdA, spdB, false);
  73. VERIFY_IS_APPROX(eig1.eigenvalues(), eig2.eigenvalues());
  74. }
  75. }
  76. EIGEN_DECLARE_TEST(eigensolver_generalized_real)
  77. {
  78. for(int i = 0; i < g_repeat; i++) {
  79. int s = 0;
  80. CALL_SUBTEST_1( generalized_eigensolver_real(Matrix4f()) );
  81. s = internal::random<int>(1,EIGEN_TEST_MAX_SIZE/4);
  82. CALL_SUBTEST_2( generalized_eigensolver_real(MatrixXd(s,s)) );
  83. // some trivial but implementation-wise special cases
  84. CALL_SUBTEST_2( generalized_eigensolver_real(MatrixXd(1,1)) );
  85. CALL_SUBTEST_2( generalized_eigensolver_real(MatrixXd(2,2)) );
  86. CALL_SUBTEST_3( generalized_eigensolver_real(Matrix<double,1,1>()) );
  87. CALL_SUBTEST_4( generalized_eigensolver_real(Matrix2d()) );
  88. TEST_SET_BUT_UNUSED_VARIABLE(s)
  89. }
  90. }