mpreal_support.cpp 2.4 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566
  1. #include <mpreal.h> // Must be included before main.h.
  2. #include "main.h"
  3. #include <Eigen/MPRealSupport>
  4. #include <Eigen/LU>
  5. #include <Eigen/Eigenvalues>
  6. #include <sstream>
  7. using namespace mpfr;
  8. using namespace Eigen;
  9. EIGEN_DECLARE_TEST(mpreal_support)
  10. {
  11. // set precision to 256 bits (double has only 53 bits)
  12. mpreal::set_default_prec(256);
  13. typedef Matrix<mpreal,Eigen::Dynamic,Eigen::Dynamic> MatrixXmp;
  14. typedef Matrix<std::complex<mpreal>,Eigen::Dynamic,Eigen::Dynamic> MatrixXcmp;
  15. std::cerr << "epsilon = " << NumTraits<mpreal>::epsilon() << "\n";
  16. std::cerr << "dummy_precision = " << NumTraits<mpreal>::dummy_precision() << "\n";
  17. std::cerr << "highest = " << NumTraits<mpreal>::highest() << "\n";
  18. std::cerr << "lowest = " << NumTraits<mpreal>::lowest() << "\n";
  19. std::cerr << "digits10 = " << NumTraits<mpreal>::digits10() << "\n";
  20. for(int i = 0; i < g_repeat; i++) {
  21. int s = Eigen::internal::random<int>(1,100);
  22. MatrixXmp A = MatrixXmp::Random(s,s);
  23. MatrixXmp B = MatrixXmp::Random(s,s);
  24. MatrixXmp S = A.adjoint() * A;
  25. MatrixXmp X;
  26. MatrixXcmp Ac = MatrixXcmp::Random(s,s);
  27. MatrixXcmp Bc = MatrixXcmp::Random(s,s);
  28. MatrixXcmp Sc = Ac.adjoint() * Ac;
  29. MatrixXcmp Xc;
  30. // Basic stuffs
  31. VERIFY_IS_APPROX(A.real(), A);
  32. VERIFY(Eigen::internal::isApprox(A.array().abs2().sum(), A.squaredNorm()));
  33. VERIFY_IS_APPROX(A.array().exp(), exp(A.array()));
  34. VERIFY_IS_APPROX(A.array().abs2().sqrt(), A.array().abs());
  35. VERIFY_IS_APPROX(A.array().sin(), sin(A.array()));
  36. VERIFY_IS_APPROX(A.array().cos(), cos(A.array()));
  37. // Cholesky
  38. X = S.selfadjointView<Lower>().llt().solve(B);
  39. VERIFY_IS_APPROX((S.selfadjointView<Lower>()*X).eval(),B);
  40. Xc = Sc.selfadjointView<Lower>().llt().solve(Bc);
  41. VERIFY_IS_APPROX((Sc.selfadjointView<Lower>()*Xc).eval(),Bc);
  42. // partial LU
  43. X = A.lu().solve(B);
  44. VERIFY_IS_APPROX((A*X).eval(),B);
  45. // symmetric eigenvalues
  46. SelfAdjointEigenSolver<MatrixXmp> eig(S);
  47. VERIFY_IS_EQUAL(eig.info(), Success);
  48. VERIFY( (S.selfadjointView<Lower>() * eig.eigenvectors()).isApprox(eig.eigenvectors() * eig.eigenvalues().asDiagonal(), NumTraits<mpreal>::dummy_precision()*1e3) );
  49. }
  50. {
  51. MatrixXmp A(8,3); A.setRandom();
  52. // test output (interesting things happen in this code)
  53. std::stringstream stream;
  54. stream << A;
  55. }
  56. }