exp_sinh_quadrature_test.cpp 24 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596
  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. #define BOOST_TEST_MODULE exp_sinh_quadrature_test
  7. #include <complex>
  8. #include <boost/multiprecision/cpp_complex.hpp>
  9. #include <boost/math/concepts/real_concept.hpp>
  10. #include <boost/test/included/unit_test.hpp>
  11. #include <boost/test/tools/floating_point_comparison.hpp>
  12. #include <boost/math/quadrature/exp_sinh.hpp>
  13. #include <boost/math/special_functions/sinc.hpp>
  14. #include <boost/math/special_functions/bessel.hpp>
  15. #include <boost/multiprecision/cpp_bin_float.hpp>
  16. #include <boost/multiprecision/cpp_dec_float.hpp>
  17. #include <boost/math/special_functions/next.hpp>
  18. #include <boost/math/special_functions/gamma.hpp>
  19. #include <boost/math/special_functions/sinc.hpp>
  20. #include <boost/type_traits/is_class.hpp>
  21. #ifdef BOOST_HAS_FLOAT128
  22. #include <boost/multiprecision/complex128.hpp>
  23. #endif
  24. using std::exp;
  25. using std::cos;
  26. using std::tan;
  27. using std::log;
  28. using std::sqrt;
  29. using std::abs;
  30. using std::sinh;
  31. using std::cosh;
  32. using std::pow;
  33. using std::atan;
  34. using boost::multiprecision::cpp_bin_float_50;
  35. using boost::multiprecision::cpp_bin_float_100;
  36. using boost::multiprecision::cpp_bin_float_quad;
  37. using boost::math::constants::pi;
  38. using boost::math::constants::half_pi;
  39. using boost::math::constants::two_div_pi;
  40. using boost::math::constants::half;
  41. using boost::math::constants::third;
  42. using boost::math::constants::half;
  43. using boost::math::constants::third;
  44. using boost::math::constants::catalan;
  45. using boost::math::constants::ln_two;
  46. using boost::math::constants::root_two;
  47. using boost::math::constants::root_two_pi;
  48. using boost::math::constants::root_pi;
  49. using boost::math::quadrature::exp_sinh;
  50. #if !defined(TEST1) && !defined(TEST2) && !defined(TEST3) && !defined(TEST4) && !defined(TEST5) && !defined(TEST6) && !defined(TEST7) && !defined(TEST8)
  51. # define TEST1
  52. # define TEST2
  53. # define TEST3
  54. # define TEST4
  55. # define TEST5
  56. # define TEST6
  57. # define TEST7
  58. # define TEST8
  59. #endif
  60. #ifdef BOOST_MSVC
  61. #pragma warning (disable:4127)
  62. #endif
  63. //
  64. // Coefficient generation code:
  65. //
  66. template <class T>
  67. void print_levels(const T& v, const char* suffix)
  68. {
  69. std::cout << "{\n";
  70. for (unsigned i = 0; i < v.size(); ++i)
  71. {
  72. std::cout << " { ";
  73. for (unsigned j = 0; j < v[i].size(); ++j)
  74. {
  75. std::cout << v[i][j] << suffix << ", ";
  76. }
  77. std::cout << "},\n";
  78. }
  79. std::cout << " };\n";
  80. }
  81. template <class T>
  82. void print_levels(const std::pair<T, T>& p, const char* suffix = "")
  83. {
  84. std::cout << " static const std::vector<std::vector<Real> > abscissa = ";
  85. print_levels(p.first, suffix);
  86. std::cout << " static const std::vector<std::vector<Real> > weights = ";
  87. print_levels(p.second, suffix);
  88. }
  89. template <class Real, class TargetType>
  90. std::pair<std::vector<std::vector<Real>>, std::vector<std::vector<Real>> > generate_constants(unsigned max_rows)
  91. {
  92. using boost::math::constants::half_pi;
  93. using boost::math::constants::two_div_pi;
  94. using boost::math::constants::pi;
  95. auto g = [](Real t)->Real { return exp(half_pi<Real>()*sinh(t)); };
  96. auto w = [](Real t)->Real { return cosh(t)*half_pi<Real>()*exp(half_pi<Real>()*sinh(t)); };
  97. std::vector<std::vector<Real>> abscissa, weights;
  98. std::vector<Real> temp;
  99. Real tmp = (Real(boost::math::tools::log_min_value<TargetType>()) + log(Real(boost::math::tools::epsilon<TargetType>())))*0.5f;
  100. Real t_min = asinh(two_div_pi<Real>()*tmp);
  101. // truncate t_min to an exact binary value:
  102. t_min = floor(t_min * 128) / 128;
  103. std::cout << "m_t_min = " << t_min << ";\n";
  104. // t_max is chosen to make g'(t_max) ~ sqrt(max) (g' grows faster than g).
  105. // This will allow some flexibility on the users part; they can at least square a number function without overflow.
  106. // But there is no unique choice; the further out we can evaluate the function, the better we can do on slowly decaying integrands.
  107. const Real t_max = log(2 * two_div_pi<Real>()*log(2 * two_div_pi<Real>()*sqrt(Real(boost::math::tools::max_value<TargetType>()))));
  108. Real h = 1;
  109. for (Real t = t_min; t < t_max; t += h)
  110. {
  111. temp.push_back(g(t));
  112. }
  113. abscissa.push_back(temp);
  114. temp.clear();
  115. for (Real t = t_min; t < t_max; t += h)
  116. {
  117. temp.push_back(w(t * h));
  118. }
  119. weights.push_back(temp);
  120. temp.clear();
  121. for (unsigned row = 1; row < max_rows; ++row)
  122. {
  123. h /= 2;
  124. for (Real t = t_min + h; t < t_max; t += 2 * h)
  125. temp.push_back(g(t));
  126. abscissa.push_back(temp);
  127. temp.clear();
  128. }
  129. h = 1;
  130. for (unsigned row = 1; row < max_rows; ++row)
  131. {
  132. h /= 2;
  133. for (Real t = t_min + h; t < t_max; t += 2 * h)
  134. temp.push_back(w(t));
  135. weights.push_back(temp);
  136. temp.clear();
  137. }
  138. return std::make_pair(abscissa, weights);
  139. }
  140. template <class Real>
  141. const exp_sinh<Real>& get_integrator()
  142. {
  143. static const exp_sinh<Real> integrator(14);
  144. return integrator;
  145. }
  146. template <class Real>
  147. Real get_convergence_tolerance()
  148. {
  149. return boost::math::tools::root_epsilon<Real>();
  150. }
  151. template<class Real>
  152. void test_right_limit_infinite()
  153. {
  154. std::cout << "Testing right limit infinite for tanh_sinh in 'A Comparison of Three High Precision Quadrature Schemes' on type " << boost::typeindex::type_id<Real>().pretty_name() << "\n";
  155. Real tol = 10 * boost::math::tools::epsilon<Real>();
  156. Real Q;
  157. Real Q_expected;
  158. Real error;
  159. Real L1;
  160. auto integrator = get_integrator<Real>();
  161. // Example 12
  162. const auto f2 = [](const Real& t)->Real { return exp(-t)/sqrt(t); };
  163. Q = integrator.integrate(f2, get_convergence_tolerance<Real>(), &error, &L1);
  164. Q_expected = root_pi<Real>();
  165. Real tol_mult = 1;
  166. // Multiprecision type have higher error rates, probably evaluation of f() is less accurate:
  167. if (std::numeric_limits<Real>::digits10 > std::numeric_limits<long double>::digits10)
  168. tol_mult = 12;
  169. else if (std::numeric_limits<Real>::digits10 > std::numeric_limits<double>::digits10)
  170. tol_mult = 5;
  171. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, tol * tol_mult);
  172. // The integrand is strictly positive, so it coincides with the value of the integral:
  173. BOOST_CHECK_CLOSE_FRACTION(L1, Q_expected, tol * tol_mult);
  174. auto f3 = [](Real t)->Real { Real z = exp(-t); if (z == 0) { return z; } return z*cos(t); };
  175. Q = integrator.integrate(f3, get_convergence_tolerance<Real>(), &error, &L1);
  176. Q_expected = half<Real>();
  177. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, tol);
  178. Q = integrator.integrate(f3, 10, std::numeric_limits<Real>::has_infinity ? std::numeric_limits<Real>::infinity() : boost::math::tools::max_value<Real>(), get_convergence_tolerance<Real>(), &error, &L1);
  179. Q_expected = boost::lexical_cast<Real>("-6.6976341310426674140007086979326069121526743314567805278252392932e-6");
  180. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, 10 * tol);
  181. // Integrating through zero risks precision loss:
  182. Q = integrator.integrate(f3, -10, std::numeric_limits<Real>::has_infinity ? std::numeric_limits<Real>::infinity() : boost::math::tools::max_value<Real>(), get_convergence_tolerance<Real>(), &error, &L1);
  183. Q_expected = boost::lexical_cast<Real>("-15232.3213626280525704332288302799653087046646639974940243044623285817777006");
  184. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, std::numeric_limits<Real>::digits10 > 30 ? 1000 * tol : tol);
  185. auto f4 = [](Real t)->Real { return 1/(1+t*t); };
  186. Q = integrator.integrate(f4, 1, std::numeric_limits<Real>::has_infinity ? std::numeric_limits<Real>::infinity() : boost::math::tools::max_value<Real>(), get_convergence_tolerance<Real>(), &error, &L1);
  187. Q_expected = pi<Real>()/4;
  188. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, tol);
  189. BOOST_CHECK_CLOSE_FRACTION(L1, Q_expected, tol);
  190. Q = integrator.integrate(f4, 20, std::numeric_limits<Real>::has_infinity ? std::numeric_limits<Real>::infinity() : boost::math::tools::max_value<Real>(), get_convergence_tolerance<Real>(), &error, &L1);
  191. Q_expected = boost::lexical_cast<Real>("0.0499583957219427614100062870348448814912770804235071744108534548299835954767");
  192. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, tol);
  193. BOOST_CHECK_CLOSE_FRACTION(L1, Q_expected, tol);
  194. Q = integrator.integrate(f4, 500, std::numeric_limits<Real>::has_infinity ? std::numeric_limits<Real>::infinity() : boost::math::tools::max_value<Real>(), get_convergence_tolerance<Real>(), &error, &L1);
  195. Q_expected = boost::lexical_cast<Real>("0.0019999973333397333150476759363217553199063513829126652556286269630");
  196. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, tol);
  197. BOOST_CHECK_CLOSE_FRACTION(L1, Q_expected, tol);
  198. }
  199. template<class Real>
  200. void test_left_limit_infinite()
  201. {
  202. std::cout << "Testing left limit infinite for 1/(1+t^2) on type " << boost::typeindex::type_id<Real>().pretty_name() << "\n";
  203. Real tol = 10 * boost::math::tools::epsilon<Real>();
  204. Real Q;
  205. Real Q_expected;
  206. Real error;
  207. Real L1;
  208. auto integrator = get_integrator<Real>();
  209. // Example 11:
  210. auto f1 = [](const Real& t)->Real { return 1/(1+t*t);};
  211. Q = integrator.integrate(f1, std::numeric_limits<Real>::has_infinity ? -std::numeric_limits<Real>::infinity() : -boost::math::tools::max_value<Real>(), 0, get_convergence_tolerance<Real>(), &error, &L1);
  212. Q_expected = half_pi<Real>();
  213. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, tol);
  214. BOOST_CHECK_CLOSE_FRACTION(L1, Q_expected, tol);
  215. Q = integrator.integrate(f1, std::numeric_limits<Real>::has_infinity ? -std::numeric_limits<Real>::infinity() : -boost::math::tools::max_value<Real>(), -20, get_convergence_tolerance<Real>(), &error, &L1);
  216. Q_expected = boost::lexical_cast<Real>("0.0499583957219427614100062870348448814912770804235071744108534548299835954767");
  217. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, tol);
  218. BOOST_CHECK_CLOSE_FRACTION(L1, Q_expected, tol);
  219. Q = integrator.integrate(f1, std::numeric_limits<Real>::has_infinity ? -std::numeric_limits<Real>::infinity() : -boost::math::tools::max_value<Real>(), -500, get_convergence_tolerance<Real>(), &error, &L1);
  220. Q_expected = boost::lexical_cast<Real>("0.0019999973333397333150476759363217553199063513829126652556286269630");
  221. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, tol);
  222. BOOST_CHECK_CLOSE_FRACTION(L1, Q_expected, tol);
  223. }
  224. // Some examples of tough integrals from NR, section 4.5.4:
  225. template<class Real>
  226. void test_nr_examples()
  227. {
  228. using std::sin;
  229. using std::cos;
  230. using std::pow;
  231. using std::exp;
  232. using std::sqrt;
  233. std::cout << "Testing type " << boost::typeindex::type_id<Real>().pretty_name() << "\n";
  234. Real tol = 10 * boost::math::tools::epsilon<Real>();
  235. std::cout << std::setprecision(std::numeric_limits<Real>::digits10);
  236. Real Q;
  237. Real Q_expected;
  238. Real L1;
  239. Real error;
  240. auto integrator = get_integrator<Real>();
  241. auto f0 = [] (Real)->Real { return (Real) 0; };
  242. Q = integrator.integrate(f0, get_convergence_tolerance<Real>(), &error, &L1);
  243. Q_expected = 0;
  244. BOOST_CHECK_CLOSE_FRACTION(Q, 0.0f, tol);
  245. BOOST_CHECK_CLOSE_FRACTION(L1, 0.0f, tol);
  246. auto f = [](const Real& x)->Real { return 1/(1+x*x); };
  247. Q = integrator.integrate(f, get_convergence_tolerance<Real>(), &error, &L1);
  248. Q_expected = half_pi<Real>();
  249. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, tol);
  250. BOOST_CHECK_CLOSE_FRACTION(L1, Q_expected, tol);
  251. auto f1 = [](Real x)->Real {
  252. Real z1 = exp(-x);
  253. if (z1 == 0)
  254. {
  255. return (Real) 0;
  256. }
  257. Real z2 = pow(x, -3*half<Real>())*z1;
  258. if (z2 == 0)
  259. {
  260. return (Real) 0;
  261. }
  262. return sin(x*half<Real>())*z2;
  263. };
  264. Q = integrator.integrate(f1, get_convergence_tolerance<Real>(), &error, &L1);
  265. Q_expected = sqrt(pi<Real>()*(sqrt((Real) 5) - 2));
  266. // The integrand is oscillatory; the accuracy is low.
  267. Real tol_mul = 1;
  268. if (std::numeric_limits<Real>::digits10 > 40)
  269. tol_mul = 500000;
  270. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, tol_mul * tol);
  271. auto f2 = [](Real x)->Real { return x > boost::math::tools::log_max_value<Real>() ? Real(0) : Real(pow(x, -(Real) 2/(Real) 7)*exp(-x*x)); };
  272. Q = integrator.integrate(f2, get_convergence_tolerance<Real>(), &error, &L1);
  273. Q_expected = half<Real>()*boost::math::tgamma((Real) 5/ (Real) 14);
  274. tol_mul = 1;
  275. if (std::numeric_limits<Real>::is_specialized == false)
  276. tol_mul = 6;
  277. else if (std::numeric_limits<Real>::digits10 > 40)
  278. tol_mul = 100;
  279. else
  280. tol_mul = 3;
  281. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, tol_mul * tol);
  282. BOOST_CHECK_CLOSE_FRACTION(L1, Q_expected, tol_mul * tol);
  283. auto f3 = [](Real x)->Real { return (Real) 1/ (sqrt(x)*(1+x)); };
  284. Q = integrator.integrate(f3, get_convergence_tolerance<Real>(), &error, &L1);
  285. Q_expected = pi<Real>();
  286. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, 10*boost::math::tools::epsilon<Real>());
  287. BOOST_CHECK_CLOSE_FRACTION(L1, Q_expected, 10*boost::math::tools::epsilon<Real>());
  288. auto f4 = [](const Real& t)->Real { return t > boost::math::tools::log_max_value<Real>() ? Real(0) : Real(exp(-t*t*half<Real>())); };
  289. Q = integrator.integrate(f4, get_convergence_tolerance<Real>(), &error, &L1);
  290. Q_expected = root_two_pi<Real>()/2;
  291. tol_mul = 1;
  292. if (std::numeric_limits<Real>::digits10 > 40)
  293. tol_mul = 5000;
  294. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, tol_mul * tol);
  295. BOOST_CHECK_CLOSE_FRACTION(L1, Q_expected, tol_mul * tol);
  296. auto f5 = [](const Real& t)->Real { return 1/cosh(t);};
  297. Q = integrator.integrate(f5, get_convergence_tolerance<Real>(), &error, &L1);
  298. Q_expected = half_pi<Real>();
  299. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, tol * 12); // Fails at float precision without higher error rate
  300. BOOST_CHECK_CLOSE_FRACTION(L1, Q_expected, tol * 12);
  301. }
  302. // Definite integrals found in the CRC Handbook of Mathematical Formulas
  303. template<class Real>
  304. void test_crc()
  305. {
  306. using std::sin;
  307. using std::pow;
  308. using std::exp;
  309. using std::sqrt;
  310. using std::log;
  311. using std::cos;
  312. std::cout << "Testing integral from CRC handbook on type " << boost::typeindex::type_id<Real>().pretty_name() << "\n";
  313. Real tol = 10 * boost::math::tools::epsilon<Real>();
  314. std::cout << std::setprecision(std::numeric_limits<Real>::digits10);
  315. Real Q;
  316. Real Q_expected;
  317. Real L1;
  318. Real error;
  319. auto integrator = get_integrator<Real>();
  320. auto f0 = [](const Real& x)->Real { return x > boost::math::tools::log_max_value<Real>() ? Real(0) : Real(log(x)*exp(-x)); };
  321. Q = integrator.integrate(f0, get_convergence_tolerance<Real>(), &error, &L1);
  322. Q_expected = -boost::math::constants::euler<Real>();
  323. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, tol);
  324. // Test the integral representation of the gamma function:
  325. auto f1 = [](Real t)->Real { Real x = exp(-t);
  326. if(x == 0)
  327. {
  328. return (Real) 0;
  329. }
  330. return pow(t, (Real) 12 - 1)*x;
  331. };
  332. Q = integrator.integrate(f1, get_convergence_tolerance<Real>(), &error, &L1);
  333. Q_expected = boost::math::tgamma(12.0f);
  334. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, tol);
  335. // Integral representation of the modified bessel function:
  336. // K_5(12)
  337. auto f2 = [](Real t)->Real {
  338. Real x = 12*cosh(t);
  339. if (x > boost::math::tools::log_max_value<Real>())
  340. {
  341. return (Real) 0;
  342. }
  343. return exp(-x)*cosh(5*t);
  344. };
  345. Q = integrator.integrate(f2, get_convergence_tolerance<Real>(), &error, &L1);
  346. Q_expected = boost::math::cyl_bessel_k<int, Real>(5, 12);
  347. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, tol);
  348. // Laplace transform of cos(at)
  349. Real a = 20;
  350. Real s = 1;
  351. auto f3 = [&](Real t)->Real {
  352. Real x = s * t;
  353. if (x > boost::math::tools::log_max_value<Real>())
  354. {
  355. return (Real) 0;
  356. }
  357. return cos(a * t) * exp(-x);
  358. };
  359. // Since the integrand is oscillatory, we increase the tolerance:
  360. Real tol_mult = 10;
  361. // Multiprecision type have higher error rates, probably evaluation of f() is less accurate:
  362. if (!boost::is_class<Real>::value)
  363. {
  364. // For high oscillation frequency, the quadrature sum is ill-conditioned.
  365. Q = integrator.integrate(f3, get_convergence_tolerance<Real>(), &error, &L1);
  366. Q_expected = s/(a*a+s*s);
  367. if (std::numeric_limits<Real>::digits10 > std::numeric_limits<double>::digits10)
  368. tol_mult = 5000; // we should really investigate this more??
  369. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, tol_mult*tol);
  370. }
  371. //
  372. // This one doesn't pass for real_concept..
  373. //
  374. if (std::numeric_limits<Real>::is_specialized)
  375. {
  376. // Laplace transform of J_0(t):
  377. auto f4 = [&](Real t)->Real {
  378. Real x = s * t;
  379. if (x > boost::math::tools::log_max_value<Real>())
  380. {
  381. return (Real)0;
  382. }
  383. return boost::math::cyl_bessel_j(0, t) * exp(-x);
  384. };
  385. Q = integrator.integrate(f4, get_convergence_tolerance<Real>(), &error, &L1);
  386. Q_expected = 1 / sqrt(1 + s*s);
  387. tol_mult = 3;
  388. // Multiprecision type have higher error rates, probably evaluation of f() is less accurate:
  389. if (std::numeric_limits<Real>::digits10 > std::numeric_limits<long double>::digits10)
  390. tol_mult = 750;
  391. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, tol_mult * tol);
  392. }
  393. auto f6 = [](const Real& t)->Real { return t > boost::math::tools::log_max_value<Real>() ? Real(0) : Real(exp(-t*t)*log(t));};
  394. Q = integrator.integrate(f6, get_convergence_tolerance<Real>(), &error, &L1);
  395. Q_expected = -boost::math::constants::root_pi<Real>()*(boost::math::constants::euler<Real>() + 2*ln_two<Real>())/4;
  396. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, tol);
  397. // CRC Section 5.5, integral 591
  398. // The parameter p allows us to control the strength of the singularity.
  399. // Rapid convergence is not guaranteed for this function, as the branch cut makes it non-analytic on a disk.
  400. // This converges only when our test type has an extended exponent range as all the area of the integral
  401. // occurs so close to 0 (or 1) that we need abscissa values exceptionally small to find it.
  402. // "There's a lot of room at the bottom".
  403. // This version is transformed via argument substitution (exp(-x) for x) so that the integral is spread
  404. // over (0, INF).
  405. tol *= boost::math::tools::digits<Real>() > 100 ? 100000 : 75;
  406. for (Real pn = 99; pn > 0; pn -= 10) {
  407. Real p = pn / 100;
  408. auto f = [&](Real x)->Real
  409. {
  410. return x > 1000 * boost::math::tools::log_max_value<Real>() ? Real(0) : Real(exp(-x * (1 - p) + p * log(-boost::math::expm1(-x))));
  411. };
  412. Q = integrator.integrate(f, get_convergence_tolerance<Real>(), &error, &L1);
  413. Q_expected = 1 / boost::math::sinc_pi(p*pi<Real>());
  414. /*
  415. std::cout << std::setprecision(std::numeric_limits<Real>::max_digits10) << p << std::endl;
  416. std::cout << std::setprecision(std::numeric_limits<Real>::max_digits10) << Q << std::endl;
  417. std::cout << std::setprecision(std::numeric_limits<Real>::max_digits10) << Q_expected << std::endl;
  418. std::cout << fabs((Q - Q_expected) / Q_expected) << std::endl;
  419. */
  420. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, tol);
  421. }
  422. // and for p < 1:
  423. for (Real p = -0.99; p < 0; p += 0.1) {
  424. auto f = [&](Real x)->Real
  425. {
  426. return x > 1000 * boost::math::tools::log_max_value<Real>() ? Real(0) : Real(exp(-p * log(-boost::math::expm1(-x)) - (1 + p) * x));
  427. };
  428. Q = integrator.integrate(f, get_convergence_tolerance<Real>(), &error, &L1);
  429. Q_expected = 1 / boost::math::sinc_pi(p*pi<Real>());
  430. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, tol);
  431. }
  432. }
  433. template<class Complex>
  434. void test_complex_modified_bessel()
  435. {
  436. std::cout << "Testing complex modified Bessel function on type " << boost::typeindex::type_id<Complex>().pretty_name() << "\n";
  437. typedef typename Complex::value_type Real;
  438. Real tol = 100 * boost::math::tools::epsilon<Real>();
  439. Real error;
  440. Real L1;
  441. auto integrator = get_integrator<Real>();
  442. // Integral Representation of Modified Complex Bessel function:
  443. // https://en.wikipedia.org/wiki/Bessel_function#Modified_Bessel_functions
  444. Complex z{2, 3};
  445. const auto f = [&z](const Real& t)->Complex
  446. {
  447. using std::cosh;
  448. using std::exp;
  449. Real cosht = cosh(t);
  450. if (cosht > boost::math::tools::log_max_value<Real>())
  451. {
  452. return Complex{0, 0};
  453. }
  454. Complex arg = -z*cosht;
  455. Complex res = exp(arg);
  456. return res;
  457. };
  458. Complex K0 = integrator.integrate(f, get_convergence_tolerance<Real>(), &error, &L1);
  459. // Mathematica code: N[BesselK[0, 2 + 3 I], 140]
  460. Real K0_x_expected = boost::lexical_cast<Real>("-0.08296852656762551490517953520589186885781541203818846830385526187936132191822538822296497597191327722262903004145527496422090506197776994");
  461. Real K0_y_expected = boost::lexical_cast<Real>("0.027949603635183423629723306332336002340909030265538548521150904238352846705644065168365102147901993976999717171115546662967229050834575193041");
  462. BOOST_CHECK_CLOSE_FRACTION(K0.real(), K0_x_expected, tol);
  463. BOOST_CHECK_CLOSE_FRACTION(K0.imag(), K0_y_expected, tol);
  464. }
  465. BOOST_AUTO_TEST_CASE(exp_sinh_quadrature_test)
  466. {
  467. //
  468. // Uncomment to generate the coefficients:
  469. //
  470. /*
  471. std::cout << std::scientific << std::setprecision(8);
  472. print_levels(generate_constants<cpp_bin_float_100, float>(8), "f");
  473. std::cout << std::setprecision(18);
  474. print_levels(generate_constants<cpp_bin_float_100, double>(8), "");
  475. std::cout << std::setprecision(35);
  476. print_levels(generate_constants<cpp_bin_float_100, cpp_bin_float_quad>(8), "L");
  477. */
  478. #ifdef TEST1
  479. test_left_limit_infinite<float>();
  480. test_right_limit_infinite<float>();
  481. test_nr_examples<float>();
  482. test_crc<float>();
  483. #endif
  484. #ifdef TEST2
  485. test_left_limit_infinite<double>();
  486. test_right_limit_infinite<double>();
  487. test_nr_examples<double>();
  488. test_crc<double>();
  489. #endif
  490. #ifndef BOOST_MATH_NO_LONG_DOUBLE_MATH_FUNCTIONS
  491. #ifdef TEST3
  492. test_left_limit_infinite<long double>();
  493. test_right_limit_infinite<long double>();
  494. test_nr_examples<long double>();
  495. test_crc<long double>();
  496. #endif
  497. #endif
  498. #ifdef TEST4
  499. test_left_limit_infinite<cpp_bin_float_quad>();
  500. test_right_limit_infinite<cpp_bin_float_quad>();
  501. test_nr_examples<cpp_bin_float_quad>();
  502. test_crc<cpp_bin_float_quad>();
  503. #endif
  504. #ifndef BOOST_MATH_NO_LONG_DOUBLE_MATH_FUNCTIONS
  505. #ifdef TEST5
  506. test_left_limit_infinite<boost::math::concepts::real_concept>();
  507. test_right_limit_infinite<boost::math::concepts::real_concept>();
  508. test_nr_examples<boost::math::concepts::real_concept>();
  509. test_crc<boost::math::concepts::real_concept>();
  510. #endif
  511. #endif
  512. #ifdef TEST6
  513. test_left_limit_infinite<boost::multiprecision::cpp_bin_float_50>();
  514. test_right_limit_infinite<boost::multiprecision::cpp_bin_float_50>();
  515. test_nr_examples<boost::multiprecision::cpp_bin_float_50>();
  516. test_crc<boost::multiprecision::cpp_bin_float_50>();
  517. #endif
  518. #ifdef TEST7
  519. test_left_limit_infinite<boost::multiprecision::cpp_dec_float_50>();
  520. test_right_limit_infinite<boost::multiprecision::cpp_dec_float_50>();
  521. test_nr_examples<boost::multiprecision::cpp_dec_float_50>();
  522. //
  523. // This one causes stack overflows on the CI machine, but not locally,
  524. // assume it's due to resticted resources on the server, and <shrug> for now...
  525. //
  526. #if ! BOOST_WORKAROUND(BOOST_MSVC, == 1900)
  527. test_crc<boost::multiprecision::cpp_dec_float_50>();
  528. #endif
  529. #endif
  530. #ifdef TEST8
  531. test_complex_modified_bessel<std::complex<float>>();
  532. test_complex_modified_bessel<std::complex<double>>();
  533. test_complex_modified_bessel<std::complex<long double>>();
  534. #ifdef BOOST_HAS_FLOAT128
  535. test_complex_modified_bessel<boost::multiprecision::complex128>();
  536. #endif
  537. test_complex_modified_bessel<boost::multiprecision::cpp_complex_quad>();
  538. #endif
  539. }