trapezoidal.hpp 4.9 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134
  1. /*
  2. * Copyright Nick Thompson, 2017
  3. * Use, modification and distribution are subject to the
  4. * Boost Software License, Version 1.0. (See accompanying file
  5. * LICENSE_1_0.txt or copy at http://www.boost.org/LICENSE_1_0.txt)
  6. *
  7. * Use the adaptive trapezoidal rule to estimate the integral of periodic functions over a period,
  8. * or to integrate a function whose derivative vanishes at the endpoints.
  9. *
  10. * If your function does not satisfy these conditions, and instead is simply continuous and bounded
  11. * over the whole interval, then this routine will still converge, albeit slowly. However, there
  12. * are much more efficient methods in this case, including Romberg, Simpson, and double exponential quadrature.
  13. */
  14. #ifndef BOOST_MATH_QUADRATURE_TRAPEZOIDAL_HPP
  15. #define BOOST_MATH_QUADRATURE_TRAPEZOIDAL_HPP
  16. #include <cmath>
  17. #include <limits>
  18. #include <stdexcept>
  19. #include <boost/math/constants/constants.hpp>
  20. #include <boost/math/special_functions/fpclassify.hpp>
  21. #include <boost/math/policies/error_handling.hpp>
  22. #include <boost/math/tools/cxx03_warn.hpp>
  23. namespace boost{ namespace math{ namespace quadrature {
  24. template<class F, class Real, class Policy>
  25. auto trapezoidal(F f, Real a, Real b, Real tol, std::size_t max_refinements, Real* error_estimate, Real* L1, const Policy& pol)->decltype(std::declval<F>()(std::declval<Real>()))
  26. {
  27. static const char* function = "boost::math::quadrature::trapezoidal<%1%>(F, %1%, %1%, %1%)";
  28. using std::abs;
  29. using boost::math::constants::half;
  30. // In many math texts, K represents the field of real or complex numbers.
  31. // Too bad we can't put blackboard bold into C++ source!
  32. typedef decltype(f(a)) K;
  33. if (!(boost::math::isfinite)(a))
  34. {
  35. return static_cast<K>(boost::math::policies::raise_domain_error(function, "Left endpoint of integration must be finite for adaptive trapezoidal integration but got a = %1%.\n", a, pol));
  36. }
  37. if (!(boost::math::isfinite)(b))
  38. {
  39. return static_cast<K>(boost::math::policies::raise_domain_error(function, "Right endpoint of integration must be finite for adaptive trapezoidal integration but got b = %1%.\n", b, pol));
  40. }
  41. if (a == b)
  42. {
  43. return static_cast<K>(0);
  44. }
  45. if(a > b)
  46. {
  47. return -trapezoidal(f, b, a, tol, max_refinements, error_estimate, L1, pol);
  48. }
  49. K ya = f(a);
  50. K yb = f(b);
  51. Real h = (b - a)*half<Real>();
  52. K I0 = (ya + yb)*h;
  53. Real IL0 = (abs(ya) + abs(yb))*h;
  54. K yh = f(a + h);
  55. K I1;
  56. I1 = I0*half<Real>() + yh*h;
  57. Real IL1 = IL0*half<Real>() + abs(yh)*h;
  58. // The recursion is:
  59. // I_k = 1/2 I_{k-1} + 1/2^k \sum_{j=1; j odd, j < 2^k} f(a + j(b-a)/2^k)
  60. std::size_t k = 2;
  61. // We want to go through at least 5 levels so we have sampled the function at least 20 times.
  62. // Otherwise, we could terminate prematurely and miss essential features.
  63. // This is of course possible anyway, but 20 samples seems to be a reasonable compromise.
  64. Real error = abs(I0 - I1);
  65. // I take k < 5, rather than k < 4, or some other smaller minimum number,
  66. // because I hit a truly exceptional bug where the k = 2 and k =3 refinement were bitwise equal,
  67. // but the quadrature had not yet converged.
  68. while (k < 5 || (k < max_refinements && error > tol*IL1) )
  69. {
  70. I0 = I1;
  71. IL0 = IL1;
  72. I1 = I0*half<Real>();
  73. IL1 = IL0*half<Real>();
  74. std::size_t p = static_cast<std::size_t>(1u) << k;
  75. h *= half<Real>();
  76. K sum = 0;
  77. Real absum = 0;
  78. for(std::size_t j = 1; j < p; j += 2)
  79. {
  80. K y = f(a + j*h);
  81. sum += y;
  82. absum += abs(y);
  83. }
  84. I1 += sum*h;
  85. IL1 += absum*h;
  86. ++k;
  87. error = abs(I0 - I1);
  88. }
  89. if (error_estimate)
  90. {
  91. *error_estimate = error;
  92. }
  93. if (L1)
  94. {
  95. *L1 = IL1;
  96. }
  97. return static_cast<K>(I1);
  98. }
  99. #if BOOST_WORKAROUND(BOOST_MSVC, < 1800)
  100. // Template argument deduction failure otherwise:
  101. template<class F, class Real>
  102. auto trapezoidal(F f, Real a, Real b, Real tol = 0, std::size_t max_refinements = 12, Real* error_estimate = 0, Real* L1 = 0)->decltype(std::declval<F>()(std::declval<Real>()))
  103. #elif !defined(BOOST_NO_CXX11_NULLPTR)
  104. template<class F, class Real>
  105. auto trapezoidal(F f, Real a, Real b, Real tol = boost::math::tools::root_epsilon<Real>(), std::size_t max_refinements = 12, Real* error_estimate = nullptr, Real* L1 = nullptr)->decltype(std::declval<F>()(std::declval<Real>()))
  106. #else
  107. template<class F, class Real>
  108. auto trapezoidal(F f, Real a, Real b, Real tol = boost::math::tools::root_epsilon<Real>(), std::size_t max_refinements = 12, Real* error_estimate = 0, Real* L1 = 0)->decltype(std::declval<F>()(std::declval<Real>()))
  109. #endif
  110. {
  111. #if BOOST_WORKAROUND(BOOST_MSVC, <= 1600)
  112. if (tol == 0)
  113. tol = boost::math::tools::root_epsilon<Real>();
  114. #endif
  115. return trapezoidal(f, a, b, tol, max_refinements, error_estimate, L1, boost::math::policies::policy<>());
  116. }
  117. }}}
  118. #endif