exp_sinh.hpp 4.1 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100
  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)->decltype(std::declval<F>()(std::declval<Real>())) const;
  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)->decltype(std::declval<F>()(std::declval<Real>())) const;
  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)->decltype(std::declval<F>()(std::declval<Real>())) const
  37. {
  38. typedef decltype(f(a)) K;
  39. using std::abs;
  40. using boost::math::constants::half;
  41. using boost::math::quadrature::detail::exp_sinh_detail;
  42. static const char* function = "boost::math::quadrature::exp_sinh<%1%>::integrate";
  43. // Neither limit may be a NaN:
  44. if((boost::math::isnan)(a) || (boost::math::isnan)(b))
  45. {
  46. 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()));
  47. }
  48. // Right limit is infinite:
  49. if ((boost::math::isfinite)(a) && (b >= boost::math::tools::max_value<Real>()))
  50. {
  51. // If a = 0, don't use an additional level of indirection:
  52. if (a == (Real) 0)
  53. {
  54. return m_imp->integrate(f, error, L1, function, tolerance, levels);
  55. }
  56. const auto u = [&](Real t)->K { return f(t + a); };
  57. return m_imp->integrate(u, error, L1, function, tolerance, levels);
  58. }
  59. if ((boost::math::isfinite)(b) && a <= -boost::math::tools::max_value<Real>())
  60. {
  61. const auto u = [&](Real t)->K { return f(b-t);};
  62. return m_imp->integrate(u, error, L1, function, tolerance, levels);
  63. }
  64. // Infinite limits:
  65. if ((a <= -boost::math::tools::max_value<Real>()) && (b >= boost::math::tools::max_value<Real>()))
  66. {
  67. 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()));
  68. }
  69. // If we get to here then both ends must necessarily be finite:
  70. 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()));
  71. }
  72. template<class Real, class Policy>
  73. template<class F>
  74. auto exp_sinh<Real, Policy>::integrate(const F& f, Real tolerance, Real* error, Real* L1, std::size_t* levels)->decltype(std::declval<F>()(std::declval<Real>())) const
  75. {
  76. static const char* function = "boost::math::quadrature::exp_sinh<%1%>::integrate";
  77. using std::abs;
  78. if (abs(tolerance) > 1) {
  79. std::string msg = std::string(__FILE__) + ":" + std::to_string(__LINE__) + ":" + std::string(function) + ": The tolerance provided is unusually large; did you confuse it with a domain bound?";
  80. throw std::domain_error(msg);
  81. }
  82. return m_imp->integrate(f, error, L1, function, tolerance, levels);
  83. }
  84. }}}
  85. #endif