exp_sinh.hpp 4.2 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101
  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 <cmath>
  16. #include <limits>
  17. #include <memory>
  18. #include <string>
  19. #include <boost/math/quadrature/detail/exp_sinh_detail.hpp>
  20. namespace boost{ namespace math{ namespace quadrature {
  21. template<class Real, class Policy = policies::policy<> >
  22. class exp_sinh
  23. {
  24. public:
  25. exp_sinh(size_t max_refinements = 9)
  26. : m_imp(std::make_shared<detail::exp_sinh_detail<Real, Policy>>(max_refinements)) {}
  27. template<class F>
  28. 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>()));
  29. template<class F>
  30. 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>()));
  31. private:
  32. std::shared_ptr<detail::exp_sinh_detail<Real, Policy>> m_imp;
  33. };
  34. template<class Real, class Policy>
  35. template<class F>
  36. 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>()))
  37. {
  38. typedef decltype(f(a)) K;
  39. static_assert(!std::is_integral<K>::value,
  40. "The return type cannot be integral, it must be either a real or complex floating point type.");
  41. using std::abs;
  42. using boost::math::constants::half;
  43. using boost::math::quadrature::detail::exp_sinh_detail;
  44. static const char* function = "boost::math::quadrature::exp_sinh<%1%>::integrate";
  45. // Neither limit may be a NaN:
  46. if((boost::math::isnan)(a) || (boost::math::isnan)(b))
  47. {
  48. 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()));
  49. }
  50. // Right limit is infinite:
  51. if ((boost::math::isfinite)(a) && (b >= boost::math::tools::max_value<Real>()))
  52. {
  53. // If a = 0, don't use an additional level of indirection:
  54. if (a == static_cast<Real>(0))
  55. {
  56. return m_imp->integrate(f, error, L1, function, tolerance, levels);
  57. }
  58. const auto u = [&](Real t)->K { return f(t + a); };
  59. return m_imp->integrate(u, error, L1, function, tolerance, levels);
  60. }
  61. if ((boost::math::isfinite)(b) && a <= -boost::math::tools::max_value<Real>())
  62. {
  63. const auto u = [&](Real t)->K { return f(b-t);};
  64. return m_imp->integrate(u, error, L1, function, tolerance, levels);
  65. }
  66. // Infinite limits:
  67. if ((a <= -boost::math::tools::max_value<Real>()) && (b >= boost::math::tools::max_value<Real>()))
  68. {
  69. 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()));
  70. }
  71. // If we get to here then both ends must necessarily be finite:
  72. 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()));
  73. }
  74. template<class Real, class Policy>
  75. template<class F>
  76. 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>()))
  77. {
  78. static const char* function = "boost::math::quadrature::exp_sinh<%1%>::integrate";
  79. using std::abs;
  80. if (abs(tolerance) > 1) {
  81. return policies::raise_domain_error(function, "The tolerance provided (%1%) is unusually large; did you confuse it with a domain bound?", tolerance, Policy());
  82. }
  83. return m_imp->integrate(f, error, L1, function, tolerance, levels);
  84. }
  85. }}}
  86. #endif