MPRealSupport 7.2 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209
  1. // This file is part of a joint effort between Eigen, a lightweight C++ template library
  2. // for linear algebra, and MPFR C++, a C++ interface to MPFR library (http://www.holoborodko.com/pavel/)
  3. //
  4. // Copyright (C) 2010-2012 Pavel Holoborodko <pavel@holoborodko.com>
  5. // Copyright (C) 2010 Konstantin Holoborodko <konstantin@holoborodko.com>
  6. // Copyright (C) 2010 Gael Guennebaud <gael.guennebaud@inria.fr>
  7. //
  8. // This Source Code Form is subject to the terms of the Mozilla
  9. // Public License v. 2.0. If a copy of the MPL was not distributed
  10. // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
  11. #ifndef EIGEN_MPREALSUPPORT_MODULE_H
  12. #define EIGEN_MPREALSUPPORT_MODULE_H
  13. #include <Eigen/Core>
  14. #include <mpreal.h>
  15. namespace Eigen {
  16. /**
  17. * \defgroup MPRealSupport_Module MPFRC++ Support module
  18. * \code
  19. * #include <Eigen/MPRealSupport>
  20. * \endcode
  21. *
  22. * This module provides support for multi precision floating point numbers
  23. * via the <a href="http://www.holoborodko.com/pavel/mpfr">MPFR C++</a>
  24. * library which itself is built upon <a href="http://www.mpfr.org/">MPFR</a>/<a href="http://gmplib.org/">GMP</a>.
  25. *
  26. * \warning MPFR C++ is licensed under the GPL.
  27. *
  28. * You can find a copy of MPFR C++ that is known to be compatible in the unsupported/test/mpreal folder.
  29. *
  30. * Here is an example:
  31. *
  32. \code
  33. #include <iostream>
  34. #include <Eigen/MPRealSupport>
  35. #include <Eigen/LU>
  36. using namespace mpfr;
  37. using namespace Eigen;
  38. int main()
  39. {
  40. // set precision to 256 bits (double has only 53 bits)
  41. mpreal::set_default_prec(256);
  42. // Declare matrix and vector types with multi-precision scalar type
  43. typedef Matrix<mpreal,Dynamic,Dynamic> MatrixXmp;
  44. typedef Matrix<mpreal,Dynamic,1> VectorXmp;
  45. MatrixXmp A = MatrixXmp::Random(100,100);
  46. VectorXmp b = VectorXmp::Random(100);
  47. // Solve Ax=b using LU
  48. VectorXmp x = A.lu().solve(b);
  49. std::cout << "relative error: " << (A*x - b).norm() / b.norm() << std::endl;
  50. return 0;
  51. }
  52. \endcode
  53. *
  54. */
  55. template<> struct NumTraits<mpfr::mpreal>
  56. : GenericNumTraits<mpfr::mpreal>
  57. {
  58. enum {
  59. IsInteger = 0,
  60. IsSigned = 1,
  61. IsComplex = 0,
  62. RequireInitialization = 1,
  63. ReadCost = HugeCost,
  64. AddCost = HugeCost,
  65. MulCost = HugeCost
  66. };
  67. typedef mpfr::mpreal Real;
  68. typedef mpfr::mpreal NonInteger;
  69. static inline Real highest (long Precision = mpfr::mpreal::get_default_prec()) { return mpfr::maxval(Precision); }
  70. static inline Real lowest (long Precision = mpfr::mpreal::get_default_prec()) { return -mpfr::maxval(Precision); }
  71. // Constants
  72. static inline Real Pi (long Precision = mpfr::mpreal::get_default_prec()) { return mpfr::const_pi(Precision); }
  73. static inline Real Euler (long Precision = mpfr::mpreal::get_default_prec()) { return mpfr::const_euler(Precision); }
  74. static inline Real Log2 (long Precision = mpfr::mpreal::get_default_prec()) { return mpfr::const_log2(Precision); }
  75. static inline Real Catalan (long Precision = mpfr::mpreal::get_default_prec()) { return mpfr::const_catalan(Precision); }
  76. static inline Real epsilon (long Precision = mpfr::mpreal::get_default_prec()) { return mpfr::machine_epsilon(Precision); }
  77. static inline Real epsilon (const Real& x) { return mpfr::machine_epsilon(x); }
  78. #ifdef MPREAL_HAVE_DYNAMIC_STD_NUMERIC_LIMITS
  79. static inline int digits10 (long Precision = mpfr::mpreal::get_default_prec()) { return std::numeric_limits<Real>::digits10(Precision); }
  80. static inline int digits10 (const Real& x) { return std::numeric_limits<Real>::digits10(x); }
  81. #endif
  82. static inline Real dummy_precision()
  83. {
  84. mpfr_prec_t weak_prec = ((mpfr::mpreal::get_default_prec()-1) * 90) / 100;
  85. return mpfr::machine_epsilon(weak_prec);
  86. }
  87. };
  88. namespace internal {
  89. template<> inline mpfr::mpreal random<mpfr::mpreal>()
  90. {
  91. return mpfr::random();
  92. }
  93. template<> inline mpfr::mpreal random<mpfr::mpreal>(const mpfr::mpreal& a, const mpfr::mpreal& b)
  94. {
  95. return a + (b-a) * random<mpfr::mpreal>();
  96. }
  97. inline bool isMuchSmallerThan(const mpfr::mpreal& a, const mpfr::mpreal& b, const mpfr::mpreal& eps)
  98. {
  99. return mpfr::abs(a) <= mpfr::abs(b) * eps;
  100. }
  101. inline bool isApprox(const mpfr::mpreal& a, const mpfr::mpreal& b, const mpfr::mpreal& eps)
  102. {
  103. return mpfr::isEqualFuzzy(a,b,eps);
  104. }
  105. inline bool isApproxOrLessThan(const mpfr::mpreal& a, const mpfr::mpreal& b, const mpfr::mpreal& eps)
  106. {
  107. return a <= b || mpfr::isEqualFuzzy(a,b,eps);
  108. }
  109. template<> inline long double cast<mpfr::mpreal,long double>(const mpfr::mpreal& x)
  110. { return x.toLDouble(); }
  111. template<> inline double cast<mpfr::mpreal,double>(const mpfr::mpreal& x)
  112. { return x.toDouble(); }
  113. template<> inline long cast<mpfr::mpreal,long>(const mpfr::mpreal& x)
  114. { return x.toLong(); }
  115. template<> inline int cast<mpfr::mpreal,int>(const mpfr::mpreal& x)
  116. { return int(x.toLong()); }
  117. // Specialize GEBP kernel and traits for mpreal (no need for peeling, nor complicated stuff)
  118. // This also permits to directly call mpfr's routines and avoid many temporaries produced by mpreal
  119. template<>
  120. class gebp_traits<mpfr::mpreal, mpfr::mpreal, false, false>
  121. {
  122. public:
  123. typedef mpfr::mpreal ResScalar;
  124. enum {
  125. Vectorizable = false,
  126. LhsPacketSize = 1,
  127. RhsPacketSize = 1,
  128. ResPacketSize = 1,
  129. NumberOfRegisters = 1,
  130. nr = 1,
  131. mr = 1,
  132. LhsProgress = 1,
  133. RhsProgress = 1
  134. };
  135. typedef ResScalar LhsPacket;
  136. typedef ResScalar RhsPacket;
  137. typedef ResScalar ResPacket;
  138. };
  139. template<typename Index, typename DataMapper, bool ConjugateLhs, bool ConjugateRhs>
  140. struct gebp_kernel<mpfr::mpreal,mpfr::mpreal,Index,DataMapper,1,1,ConjugateLhs,ConjugateRhs>
  141. {
  142. typedef mpfr::mpreal mpreal;
  143. EIGEN_DONT_INLINE
  144. void operator()(const DataMapper& res, const mpreal* blockA, const mpreal* blockB,
  145. Index rows, Index depth, Index cols, const mpreal& alpha,
  146. Index strideA=-1, Index strideB=-1, Index offsetA=0, Index offsetB=0)
  147. {
  148. if(rows==0 || cols==0 || depth==0)
  149. return;
  150. mpreal acc1(0,mpfr_get_prec(blockA[0].mpfr_srcptr())),
  151. tmp (0,mpfr_get_prec(blockA[0].mpfr_srcptr()));
  152. if(strideA==-1) strideA = depth;
  153. if(strideB==-1) strideB = depth;
  154. for(Index i=0; i<rows; ++i)
  155. {
  156. for(Index j=0; j<cols; ++j)
  157. {
  158. const mpreal *A = blockA + i*strideA + offsetA;
  159. const mpreal *B = blockB + j*strideB + offsetB;
  160. acc1 = 0;
  161. for(Index k=0; k<depth; k++)
  162. {
  163. mpfr_mul(tmp.mpfr_ptr(), A[k].mpfr_srcptr(), B[k].mpfr_srcptr(), mpreal::get_default_rnd());
  164. mpfr_add(acc1.mpfr_ptr(), acc1.mpfr_ptr(), tmp.mpfr_ptr(), mpreal::get_default_rnd());
  165. }
  166. mpfr_mul(acc1.mpfr_ptr(), acc1.mpfr_srcptr(), alpha.mpfr_srcptr(), mpreal::get_default_rnd());
  167. mpfr_add(res(i,j).mpfr_ptr(), res(i,j).mpfr_srcptr(), acc1.mpfr_srcptr(), mpreal::get_default_rnd());
  168. }
  169. }
  170. }
  171. };
  172. } // end namespace internal
  173. }
  174. #endif // EIGEN_MPREALSUPPORT_MODULE_H