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.

203 lines
7.3 KiB

  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. * You can find a copy of MPFR C++ that is known to be compatible in the unsupported/test/mpreal folder.
  27. *
  28. * Here is an example:
  29. *
  30. \code
  31. #include <iostream>
  32. #include <Eigen/MPRealSupport>
  33. #include <Eigen/LU>
  34. using namespace mpfr;
  35. using namespace Eigen;
  36. int main()
  37. {
  38. // set precision to 256 bits (double has only 53 bits)
  39. mpreal::set_default_prec(256);
  40. // Declare matrix and vector types with multi-precision scalar type
  41. typedef Matrix<mpreal,Dynamic,Dynamic> MatrixXmp;
  42. typedef Matrix<mpreal,Dynamic,1> VectorXmp;
  43. MatrixXmp A = MatrixXmp::Random(100,100);
  44. VectorXmp b = VectorXmp::Random(100);
  45. // Solve Ax=b using LU
  46. VectorXmp x = A.lu().solve(b);
  47. std::cout << "relative error: " << (A*x - b).norm() / b.norm() << std::endl;
  48. return 0;
  49. }
  50. \endcode
  51. *
  52. */
  53. template<> struct NumTraits<mpfr::mpreal>
  54. : GenericNumTraits<mpfr::mpreal>
  55. {
  56. enum {
  57. IsInteger = 0,
  58. IsSigned = 1,
  59. IsComplex = 0,
  60. RequireInitialization = 1,
  61. ReadCost = 10,
  62. AddCost = 10,
  63. MulCost = 40
  64. };
  65. typedef mpfr::mpreal Real;
  66. typedef mpfr::mpreal NonInteger;
  67. inline static Real highest (long Precision = mpfr::mpreal::get_default_prec()) { return mpfr::maxval(Precision); }
  68. inline static Real lowest (long Precision = mpfr::mpreal::get_default_prec()) { return -mpfr::maxval(Precision); }
  69. // Constants
  70. inline static Real Pi (long Precision = mpfr::mpreal::get_default_prec()) { return mpfr::const_pi(Precision); }
  71. inline static Real Euler (long Precision = mpfr::mpreal::get_default_prec()) { return mpfr::const_euler(Precision); }
  72. inline static Real Log2 (long Precision = mpfr::mpreal::get_default_prec()) { return mpfr::const_log2(Precision); }
  73. inline static Real Catalan (long Precision = mpfr::mpreal::get_default_prec()) { return mpfr::const_catalan(Precision); }
  74. inline static Real epsilon (long Precision = mpfr::mpreal::get_default_prec()) { return mpfr::machine_epsilon(Precision); }
  75. inline static Real epsilon (const Real& x) { return mpfr::machine_epsilon(x); }
  76. inline static Real dummy_precision()
  77. {
  78. unsigned int weak_prec = ((mpfr::mpreal::get_default_prec()-1) * 90) / 100;
  79. return mpfr::machine_epsilon(weak_prec);
  80. }
  81. };
  82. namespace internal {
  83. template<> inline mpfr::mpreal random<mpfr::mpreal>()
  84. {
  85. return mpfr::random();
  86. }
  87. template<> inline mpfr::mpreal random<mpfr::mpreal>(const mpfr::mpreal& a, const mpfr::mpreal& b)
  88. {
  89. return a + (b-a) * random<mpfr::mpreal>();
  90. }
  91. inline bool isMuchSmallerThan(const mpfr::mpreal& a, const mpfr::mpreal& b, const mpfr::mpreal& eps)
  92. {
  93. return mpfr::abs(a) <= mpfr::abs(b) * eps;
  94. }
  95. inline bool isApprox(const mpfr::mpreal& a, const mpfr::mpreal& b, const mpfr::mpreal& eps)
  96. {
  97. return mpfr::isEqualFuzzy(a,b,eps);
  98. }
  99. inline bool isApproxOrLessThan(const mpfr::mpreal& a, const mpfr::mpreal& b, const mpfr::mpreal& eps)
  100. {
  101. return a <= b || mpfr::isEqualFuzzy(a,b,eps);
  102. }
  103. template<> inline long double cast<mpfr::mpreal,long double>(const mpfr::mpreal& x)
  104. { return x.toLDouble(); }
  105. template<> inline double cast<mpfr::mpreal,double>(const mpfr::mpreal& x)
  106. { return x.toDouble(); }
  107. template<> inline long cast<mpfr::mpreal,long>(const mpfr::mpreal& x)
  108. { return x.toLong(); }
  109. template<> inline int cast<mpfr::mpreal,int>(const mpfr::mpreal& x)
  110. { return int(x.toLong()); }
  111. // Specialize GEBP kernel and traits for mpreal (no need for peeling, nor complicated stuff)
  112. // This also permits to directly call mpfr's routines and avoid many temporaries produced by mpreal
  113. template<>
  114. class gebp_traits<mpfr::mpreal, mpfr::mpreal, false, false>
  115. {
  116. public:
  117. typedef mpfr::mpreal ResScalar;
  118. enum {
  119. nr = 2, // must be 2 for proper packing...
  120. mr = 1,
  121. WorkSpaceFactor = nr,
  122. LhsProgress = 1,
  123. RhsProgress = 1
  124. };
  125. };
  126. template<typename Index, int mr, int nr, bool ConjugateLhs, bool ConjugateRhs>
  127. struct gebp_kernel<mpfr::mpreal,mpfr::mpreal,Index,mr,nr,ConjugateLhs,ConjugateRhs>
  128. {
  129. typedef mpfr::mpreal mpreal;
  130. EIGEN_DONT_INLINE
  131. void operator()(mpreal* res, Index resStride, const mpreal* blockA, const mpreal* blockB, Index rows, Index depth, Index cols, mpreal alpha,
  132. Index strideA=-1, Index strideB=-1, Index offsetA=0, Index offsetB=0, mpreal* /*unpackedB*/ = 0)
  133. {
  134. mpreal acc1, acc2, tmp;
  135. if(strideA==-1) strideA = depth;
  136. if(strideB==-1) strideB = depth;
  137. for(Index j=0; j<cols; j+=nr)
  138. {
  139. Index actual_nr = (std::min<Index>)(nr,cols-j);
  140. mpreal *C1 = res + j*resStride;
  141. mpreal *C2 = res + (j+1)*resStride;
  142. for(Index i=0; i<rows; i++)
  143. {
  144. mpreal *B = const_cast<mpreal*>(blockB) + j*strideB + offsetB*actual_nr;
  145. mpreal *A = const_cast<mpreal*>(blockA) + i*strideA + offsetA;
  146. acc1 = 0;
  147. acc2 = 0;
  148. for(Index k=0; k<depth; k++)
  149. {
  150. mpfr_mul(tmp.mpfr_ptr(), A[k].mpfr_ptr(), B[0].mpfr_ptr(), mpreal::get_default_rnd());
  151. mpfr_add(acc1.mpfr_ptr(), acc1.mpfr_ptr(), tmp.mpfr_ptr(), mpreal::get_default_rnd());
  152. if(actual_nr==2) {
  153. mpfr_mul(tmp.mpfr_ptr(), A[k].mpfr_ptr(), B[1].mpfr_ptr(), mpreal::get_default_rnd());
  154. mpfr_add(acc2.mpfr_ptr(), acc2.mpfr_ptr(), tmp.mpfr_ptr(), mpreal::get_default_rnd());
  155. }
  156. B+=actual_nr;
  157. }
  158. mpfr_mul(acc1.mpfr_ptr(), acc1.mpfr_ptr(), alpha.mpfr_ptr(), mpreal::get_default_rnd());
  159. mpfr_add(C1[i].mpfr_ptr(), C1[i].mpfr_ptr(), acc1.mpfr_ptr(), mpreal::get_default_rnd());
  160. if(actual_nr==2) {
  161. mpfr_mul(acc2.mpfr_ptr(), acc2.mpfr_ptr(), alpha.mpfr_ptr(), mpreal::get_default_rnd());
  162. mpfr_add(C2[i].mpfr_ptr(), C2[i].mpfr_ptr(), acc2.mpfr_ptr(), mpreal::get_default_rnd());
  163. }
  164. }
  165. }
  166. }
  167. };
  168. } // end namespace internal
  169. }
  170. #endif // EIGEN_MPREALSUPPORT_MODULE_H