exp_sinh.hpp 7.4 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181
  1. // Copyright Nick Thompson, 2017
  2. // Use, modification and distribution are subject to the
  3. // Boost Software License, Version 1.0.
  4. // (See accompanying file LICENSE_1_0.txt
  5. // or copy at http://www.boost.org/LICENSE_1_0.txt)
  6. /*
  7. * This class performs exp-sinh quadrature on half infinite intervals.
  8. *
  9. * References:
  10. *
  11. * 1) Tanaka, Ken'ichiro, et al. "Function classes for double exponential integration formulas." Numerische Mathematik 111.4 (2009): 631-655.
  12. */
  13. #ifndef BOOST_MATH_QUADRATURE_EXP_SINH_HPP
  14. #define BOOST_MATH_QUADRATURE_EXP_SINH_HPP
  15. #include <boost/math/tools/config.hpp>
  16. #include <boost/math/quadrature/detail/exp_sinh_detail.hpp>
  17. #ifndef BOOST_MATH_HAS_NVRTC
  18. #include <cmath>
  19. #include <limits>
  20. #include <memory>
  21. #include <string>
  22. namespace boost{ namespace math{ namespace quadrature {
  23. template<class Real, class Policy = policies::policy<> >
  24. class exp_sinh
  25. {
  26. public:
  27. exp_sinh(size_t max_refinements = 9)
  28. : m_imp(std::make_shared<detail::exp_sinh_detail<Real, Policy>>(max_refinements)) {}
  29. template<class F>
  30. auto integrate(const F& f, Real a, Real b, Real tol = boost::math::tools::root_epsilon<Real>(), Real* error = nullptr, Real* L1 = nullptr, std::size_t* levels = nullptr) const ->decltype(std::declval<F>()(std::declval<Real>()));
  31. template<class F>
  32. auto integrate(const F& f, Real tol = boost::math::tools::root_epsilon<Real>(), Real* error = nullptr, Real* L1 = nullptr, std::size_t* levels = nullptr) const ->decltype(std::declval<F>()(std::declval<Real>()));
  33. private:
  34. std::shared_ptr<detail::exp_sinh_detail<Real, Policy>> m_imp;
  35. };
  36. template<class Real, class Policy>
  37. template<class F>
  38. auto exp_sinh<Real, Policy>::integrate(const F& f, Real a, Real b, Real tolerance, Real* error, Real* L1, std::size_t* levels) const ->decltype(std::declval<F>()(std::declval<Real>()))
  39. {
  40. typedef decltype(f(a)) K;
  41. static_assert(!std::is_integral<K>::value,
  42. "The return type cannot be integral, it must be either a real or complex floating point type.");
  43. using std::abs;
  44. using boost::math::constants::half;
  45. using boost::math::quadrature::detail::exp_sinh_detail;
  46. static const char* function = "boost::math::quadrature::exp_sinh<%1%>::integrate";
  47. // Neither limit may be a NaN:
  48. if((boost::math::isnan)(a) || (boost::math::isnan)(b))
  49. {
  50. return static_cast<K>(policies::raise_domain_error(function, "NaN supplied as one limit of integration - sorry I don't know what to do", a, Policy()));
  51. }
  52. // Right limit is infinite:
  53. if ((boost::math::isfinite)(a) && (b >= boost::math::tools::max_value<Real>()))
  54. {
  55. // If a = 0, don't use an additional level of indirection:
  56. if (a == static_cast<Real>(0))
  57. {
  58. return m_imp->integrate(f, error, L1, function, tolerance, levels);
  59. }
  60. const auto u = [&](Real t)->K { return f(t + a); };
  61. return m_imp->integrate(u, error, L1, function, tolerance, levels);
  62. }
  63. if ((boost::math::isfinite)(b) && a <= -boost::math::tools::max_value<Real>())
  64. {
  65. const auto u = [&](Real t)->K { return f(b-t);};
  66. return m_imp->integrate(u, error, L1, function, tolerance, levels);
  67. }
  68. // Infinite limits:
  69. if ((a <= -boost::math::tools::max_value<Real>()) && (b >= boost::math::tools::max_value<Real>()))
  70. {
  71. return static_cast<K>(policies::raise_domain_error(function, "Use sinh_sinh quadrature for integration over the whole real line; exp_sinh is for half infinite integrals.", a, Policy()));
  72. }
  73. // If we get to here then both ends must necessarily be finite:
  74. return static_cast<K>(policies::raise_domain_error(function, "Use tanh_sinh quadrature for integration over finite domains; exp_sinh is for half infinite integrals.", a, Policy()));
  75. }
  76. template<class Real, class Policy>
  77. template<class F>
  78. auto exp_sinh<Real, Policy>::integrate(const F& f, Real tolerance, Real* error, Real* L1, std::size_t* levels) const ->decltype(std::declval<F>()(std::declval<Real>()))
  79. {
  80. static const char* function = "boost::math::quadrature::exp_sinh<%1%>::integrate";
  81. using std::abs;
  82. if (abs(tolerance) > 1) {
  83. return policies::raise_domain_error(function, "The tolerance provided (%1%) is unusually large; did you confuse it with a domain bound?", tolerance, Policy());
  84. }
  85. return m_imp->integrate(f, error, L1, function, tolerance, levels);
  86. }
  87. }}}
  88. #endif // BOOST_MATH_HAS_NVRTC
  89. #ifdef BOOST_MATH_ENABLE_CUDA
  90. #include <boost/math/tools/type_traits.hpp>
  91. #include <boost/math/tools/cstdint.hpp>
  92. #include <boost/math/tools/precision.hpp>
  93. #include <boost/math/policies/error_handling.hpp>
  94. #include <boost/math/constants/constants.hpp>
  95. namespace boost {
  96. namespace math {
  97. namespace quadrature {
  98. template <class F, class Real, class Policy = policies::policy<> >
  99. __device__ auto exp_sinh_integrate(const F& f, Real a, Real b, Real tolerance, Real* error, Real* L1, boost::math::size_t* levels)
  100. {
  101. BOOST_MATH_STD_USING
  102. using K = decltype(f(a));
  103. static_assert(!boost::math::is_integral<K>::value,
  104. "The return type cannot be integral, it must be either a real or complex floating point type.");
  105. constexpr auto function = "boost::math::quadrature::exp_sinh<%1%>::integrate";
  106. // Neither limit may be a NaN:
  107. if((boost::math::isnan)(a) || (boost::math::isnan)(b))
  108. {
  109. return static_cast<K>(policies::raise_domain_error(function, "NaN supplied as one limit of integration - sorry I don't know what to do", a, Policy()));
  110. }
  111. // Right limit is infinite:
  112. if ((boost::math::isfinite)(a) && (b >= boost::math::tools::max_value<Real>()))
  113. {
  114. // If a = 0, don't use an additional level of indirection:
  115. if (a == static_cast<Real>(0))
  116. {
  117. return detail::exp_sinh_integrate_impl(f, tolerance, error, L1, levels);
  118. }
  119. const auto u = [&](Real t)->K { return f(t + a); };
  120. return detail::exp_sinh_integrate_impl(u, tolerance, error, L1, levels);
  121. }
  122. if ((boost::math::isfinite)(b) && a <= -boost::math::tools::max_value<Real>())
  123. {
  124. const auto u = [&](Real t)->K { return f(b-t);};
  125. return detail::exp_sinh_integrate_impl(u, tolerance, error, L1, levels);
  126. }
  127. // Infinite limits:
  128. if ((a <= -boost::math::tools::max_value<Real>()) && (b >= boost::math::tools::max_value<Real>()))
  129. {
  130. return static_cast<K>(policies::raise_domain_error(function, "Use sinh_sinh quadrature for integration over the whole real line; exp_sinh is for half infinite integrals.", a, Policy()));
  131. }
  132. // If we get to here then both ends must necessarily be finite:
  133. return static_cast<K>(policies::raise_domain_error(function, "Use tanh_sinh quadrature for integration over finite domains; exp_sinh is for half infinite integrals.", a, Policy()));
  134. }
  135. template <class F, class Real, class Policy = policies::policy<> >
  136. __device__ auto exp_sinh_integrate(const F& f, Real tolerance, Real* error, Real* L1, boost::math::size_t* levels)
  137. {
  138. BOOST_MATH_STD_USING
  139. constexpr auto function = "boost::math::quadrature::exp_sinh<%1%>::integrate";
  140. if (abs(tolerance) > 1) {
  141. return policies::raise_domain_error(function, "The tolerance provided (%1%) is unusually large; did you confuse it with a domain bound?", tolerance, Policy());
  142. }
  143. return detail::exp_sinh_integrate_impl(f, tolerance, error, L1, levels);
  144. }
  145. } // namespace quadrature
  146. } // namespace math
  147. } // namespace boost
  148. #endif // BOOST_MATH_ENABLE_CUDA
  149. #endif // BOOST_MATH_QUADRATURE_EXP_SINH_HPP