gauss_quadrature_test.cpp 17 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515
  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 tanh_sinh_quadrature_test
  7. #include <complex>
  8. //#include <boost/multiprecision/mpc.hpp>
  9. #include <boost/config.hpp>
  10. #include <boost/detail/workaround.hpp>
  11. #if !defined(BOOST_NO_CXX11_DECLTYPE) && !defined(BOOST_NO_CXX11_TRAILING_RESULT_TYPES) && !defined(BOOST_NO_SFINAE_EXPR)
  12. #include <boost/math/concepts/real_concept.hpp>
  13. #include <boost/test/included/unit_test.hpp>
  14. #include <boost/test/tools/floating_point_comparison.hpp>
  15. #include <boost/math/quadrature/gauss.hpp>
  16. #include <boost/math/special_functions/sinc.hpp>
  17. #include <boost/multiprecision/cpp_bin_float.hpp>
  18. #include <boost/multiprecision/cpp_complex.hpp>
  19. #ifdef BOOST_HAS_FLOAT128
  20. #include <boost/multiprecision/complex128.hpp>
  21. #endif
  22. #ifdef _MSC_VER
  23. #pragma warning(disable:4127) // Conditional expression is constant
  24. #endif
  25. #if !defined(TEST1) && !defined(TEST2) && !defined(TEST3)
  26. # define TEST1
  27. # define TEST2
  28. # define TEST3
  29. #endif
  30. using std::expm1;
  31. using std::atan;
  32. using std::tan;
  33. using std::log;
  34. using std::log1p;
  35. using std::asinh;
  36. using std::atanh;
  37. using std::sqrt;
  38. using std::isnormal;
  39. using std::abs;
  40. using std::sinh;
  41. using std::tanh;
  42. using std::cosh;
  43. using std::pow;
  44. using std::exp;
  45. using std::sin;
  46. using std::cos;
  47. using std::string;
  48. using boost::math::quadrature::gauss;
  49. using boost::math::constants::pi;
  50. using boost::math::constants::half_pi;
  51. using boost::math::constants::two_div_pi;
  52. using boost::math::constants::two_pi;
  53. using boost::math::constants::half;
  54. using boost::math::constants::third;
  55. using boost::math::constants::half;
  56. using boost::math::constants::third;
  57. using boost::math::constants::catalan;
  58. using boost::math::constants::ln_two;
  59. using boost::math::constants::root_two;
  60. using boost::math::constants::root_two_pi;
  61. using boost::math::constants::root_pi;
  62. using boost::multiprecision::cpp_bin_float_quad;
  63. //
  64. // Error rates depend only on the number of points in the approximation, not the type being tested,
  65. // define all our expected errors here:
  66. //
  67. enum
  68. {
  69. test_ca_error_id,
  70. test_ca_error_id_2,
  71. test_three_quad_error_id,
  72. test_three_quad_error_id_2,
  73. test_integration_over_real_line_error_id,
  74. test_right_limit_infinite_error_id,
  75. test_left_limit_infinite_error_id
  76. };
  77. template <unsigned Points>
  78. double expected_error(unsigned)
  79. {
  80. return 0; // placeholder, all tests will fail
  81. }
  82. template <>
  83. double expected_error<7>(unsigned id)
  84. {
  85. switch (id)
  86. {
  87. case test_ca_error_id:
  88. return 1e-7;
  89. case test_ca_error_id_2:
  90. return 2e-5;
  91. case test_three_quad_error_id:
  92. return 1e-8;
  93. case test_three_quad_error_id_2:
  94. return 3.5e-3;
  95. case test_integration_over_real_line_error_id:
  96. return 6e-3;
  97. case test_right_limit_infinite_error_id:
  98. case test_left_limit_infinite_error_id:
  99. return 1e-5;
  100. }
  101. return 0; // placeholder, all tests will fail
  102. }
  103. template <>
  104. double expected_error<9>(unsigned id)
  105. {
  106. switch (id)
  107. {
  108. case test_ca_error_id:
  109. return 1e-7;
  110. case test_ca_error_id_2:
  111. return 2e-5;
  112. case test_three_quad_error_id:
  113. return 1e-8;
  114. case test_three_quad_error_id_2:
  115. return 3.5e-3;
  116. case test_integration_over_real_line_error_id:
  117. return 6e-3;
  118. case test_right_limit_infinite_error_id:
  119. case test_left_limit_infinite_error_id:
  120. return 1e-5;
  121. }
  122. return 0; // placeholder, all tests will fail
  123. }
  124. template <>
  125. double expected_error<10>(unsigned id)
  126. {
  127. switch (id)
  128. {
  129. case test_ca_error_id:
  130. return 1e-12;
  131. case test_ca_error_id_2:
  132. return 3e-6;
  133. case test_three_quad_error_id:
  134. return 2e-13;
  135. case test_three_quad_error_id_2:
  136. return 2e-3;
  137. case test_integration_over_real_line_error_id:
  138. return 6e-3; // doesn't get any better with more points!
  139. case test_right_limit_infinite_error_id:
  140. case test_left_limit_infinite_error_id:
  141. return 5e-8;
  142. }
  143. return 0; // placeholder, all tests will fail
  144. }
  145. template <>
  146. double expected_error<15>(unsigned id)
  147. {
  148. switch (id)
  149. {
  150. case test_ca_error_id:
  151. return 6e-20;
  152. case test_ca_error_id_2:
  153. return 3e-7;
  154. case test_three_quad_error_id:
  155. return 1e-19;
  156. case test_three_quad_error_id_2:
  157. return 6e-4;
  158. case test_integration_over_real_line_error_id:
  159. return 6e-3; // doesn't get any better with more points!
  160. case test_right_limit_infinite_error_id:
  161. case test_left_limit_infinite_error_id:
  162. return 5e-11;
  163. }
  164. return 0; // placeholder, all tests will fail
  165. }
  166. template <>
  167. double expected_error<20>(unsigned id)
  168. {
  169. switch (id)
  170. {
  171. case test_ca_error_id:
  172. return 1e-26;
  173. case test_ca_error_id_2:
  174. return 1e-7;
  175. case test_three_quad_error_id:
  176. return 3e-27;
  177. case test_three_quad_error_id_2:
  178. return 3e-4;
  179. case test_integration_over_real_line_error_id:
  180. return 5e-5; // doesn't get any better with more points!
  181. case test_right_limit_infinite_error_id:
  182. case test_left_limit_infinite_error_id:
  183. return 1e-15;
  184. }
  185. return 0; // placeholder, all tests will fail
  186. }
  187. template <>
  188. double expected_error<25>(unsigned id)
  189. {
  190. switch (id)
  191. {
  192. case test_ca_error_id:
  193. return 5e-33;
  194. case test_ca_error_id_2:
  195. return 1e-8;
  196. case test_three_quad_error_id:
  197. return 1e-32;
  198. case test_three_quad_error_id_2:
  199. return 3e-4;
  200. case test_integration_over_real_line_error_id:
  201. return 1e-14;
  202. case test_right_limit_infinite_error_id:
  203. case test_left_limit_infinite_error_id:
  204. return 3e-19;
  205. }
  206. return 0; // placeholder, all tests will fail
  207. }
  208. template <>
  209. double expected_error<30>(unsigned id)
  210. {
  211. switch (id)
  212. {
  213. case test_ca_error_id:
  214. return 2e-34;
  215. case test_ca_error_id_2:
  216. return 5e-9;
  217. case test_three_quad_error_id:
  218. return 4e-34;
  219. case test_three_quad_error_id_2:
  220. return 1e-4;
  221. case test_integration_over_real_line_error_id:
  222. return 1e-16;
  223. case test_right_limit_infinite_error_id:
  224. case test_left_limit_infinite_error_id:
  225. return 3e-23;
  226. }
  227. return 0; // placeholder, all tests will fail
  228. }
  229. template<class Real, unsigned Points>
  230. void test_linear()
  231. {
  232. std::cout << "Testing linear functions are integrated properly by gauss on type " << boost::typeindex::type_id<Real>().pretty_name() << "\n";
  233. Real tol = boost::math::tools::epsilon<Real>() * 10;
  234. auto f = [](const Real& x)
  235. {
  236. return 5*x + 7;
  237. };
  238. Real L1;
  239. Real Q = gauss<Real, Points>::integrate(f, (Real) 0, (Real) 1, &L1);
  240. BOOST_CHECK_CLOSE_FRACTION(Q, 9.5, tol);
  241. BOOST_CHECK_CLOSE_FRACTION(L1, 9.5, tol);
  242. }
  243. template<class Real, unsigned Points>
  244. void test_quadratic()
  245. {
  246. std::cout << "Testing quadratic functions are integrated properly by Gaussian quadrature on type " << boost::typeindex::type_id<Real>().pretty_name() << "\n";
  247. Real tol = boost::math::tools::epsilon<Real>() * 10;
  248. auto f = [](const Real& x) { return 5*x*x + 7*x + 12; };
  249. Real L1;
  250. Real Q = gauss<Real, Points>::integrate(f, 0, 1, &L1);
  251. BOOST_CHECK_CLOSE_FRACTION(Q, (Real) 17 + half<Real>()*third<Real>(), tol);
  252. BOOST_CHECK_CLOSE_FRACTION(L1, (Real) 17 + half<Real>()*third<Real>(), tol);
  253. }
  254. // Examples taken from
  255. //http://crd-legacy.lbl.gov/~dhbailey/dhbpapers/quadrature.pdf
  256. template<class Real, unsigned Points>
  257. void test_ca()
  258. {
  259. std::cout << "Testing integration of C(a) on type " << boost::typeindex::type_id<Real>().pretty_name() << "\n";
  260. Real tol = expected_error<Points>(test_ca_error_id);
  261. Real L1;
  262. auto f1 = [](const Real& x) { return atan(x)/(x*(x*x + 1)) ; };
  263. Real Q = gauss<Real, Points>::integrate(f1, 0, 1, &L1);
  264. Real Q_expected = pi<Real>()*ln_two<Real>()/8 + catalan<Real>()*half<Real>();
  265. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, tol);
  266. BOOST_CHECK_CLOSE_FRACTION(L1, Q_expected, tol);
  267. auto f2 = [](Real x)->Real { Real t0 = x*x + 1; Real t1 = sqrt(t0); return atan(t1)/(t0*t1); };
  268. Q = gauss<Real, Points>::integrate(f2, 0 , 1, &L1);
  269. Q_expected = pi<Real>()/4 - pi<Real>()/root_two<Real>() + 3*atan(root_two<Real>())/root_two<Real>();
  270. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, tol);
  271. BOOST_CHECK_CLOSE_FRACTION(L1, Q_expected, tol);
  272. tol = expected_error<Points>(test_ca_error_id_2);
  273. auto f5 = [](Real t)->Real { return t*t*log(t)/((t*t - 1)*(t*t*t*t + 1)); };
  274. Q = gauss<Real, Points>::integrate(f5, 0 , 1);
  275. Q_expected = pi<Real>()*pi<Real>()*(2 - root_two<Real>())/32;
  276. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, tol);
  277. }
  278. template<class Real, unsigned Points>
  279. void test_three_quadrature_schemes_examples()
  280. {
  281. std::cout << "Testing integral in 'A Comparison of Three High Precision Quadrature Schemes' on type " << boost::typeindex::type_id<Real>().pretty_name() << "\n";
  282. Real tol = expected_error<Points>(test_three_quad_error_id);
  283. Real Q;
  284. Real Q_expected;
  285. // Example 1:
  286. auto f1 = [](const Real& t) { return t*boost::math::log1p(t); };
  287. Q = gauss<Real, Points>::integrate(f1, 0 , 1);
  288. Q_expected = half<Real>()*half<Real>();
  289. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, tol);
  290. // Example 2:
  291. auto f2 = [](const Real& t) { return t*t*atan(t); };
  292. Q = gauss<Real, Points>::integrate(f2, 0 , 1);
  293. Q_expected = (pi<Real>() -2 + 2*ln_two<Real>())/12;
  294. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, 2 * tol);
  295. // Example 3:
  296. auto f3 = [](const Real& t) { return exp(t)*cos(t); };
  297. Q = gauss<Real, Points>::integrate(f3, 0, half_pi<Real>());
  298. Q_expected = boost::math::expm1(half_pi<Real>())*half<Real>();
  299. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, tol);
  300. // Example 4:
  301. auto f4 = [](Real x)->Real { Real t0 = sqrt(x*x + 2); return atan(t0)/(t0*(x*x+1)); };
  302. Q = gauss<Real, Points>::integrate(f4, 0 , 1);
  303. Q_expected = 5*pi<Real>()*pi<Real>()/96;
  304. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, tol);
  305. tol = expected_error<Points>(test_three_quad_error_id_2);
  306. // Example 5:
  307. auto f5 = [](const Real& t) { return sqrt(t)*log(t); };
  308. Q = gauss<Real, Points>::integrate(f5, 0 , 1);
  309. Q_expected = -4/ (Real) 9;
  310. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, tol);
  311. // Example 6:
  312. auto f6 = [](const Real& t) { return sqrt(1 - t*t); };
  313. Q = gauss<Real, Points>::integrate(f6, 0 , 1);
  314. Q_expected = pi<Real>()/4;
  315. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, tol);
  316. }
  317. template<class Real, unsigned Points>
  318. void test_integration_over_real_line()
  319. {
  320. std::cout << "Testing integrals over entire real line in 'A Comparison of Three High Precision Quadrature Schemes' on type " << boost::typeindex::type_id<Real>().pretty_name() << "\n";
  321. Real tol = expected_error<Points>(test_integration_over_real_line_error_id);
  322. Real Q;
  323. Real Q_expected;
  324. Real L1;
  325. auto f1 = [](const Real& t) { return 1/(1+t*t);};
  326. Q = gauss<Real, Points>::integrate(f1, -boost::math::tools::max_value<Real>(), boost::math::tools::max_value<Real>(), &L1);
  327. Q_expected = pi<Real>();
  328. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, tol);
  329. BOOST_CHECK_CLOSE_FRACTION(L1, Q_expected, tol);
  330. }
  331. template<class Real, unsigned Points>
  332. void test_right_limit_infinite()
  333. {
  334. std::cout << "Testing right limit infinite for Gaussian quadrature in 'A Comparison of Three High Precision Quadrature Schemes' on type " << boost::typeindex::type_id<Real>().pretty_name() << "\n";
  335. Real tol = expected_error<Points>(test_right_limit_infinite_error_id);
  336. Real Q;
  337. Real Q_expected;
  338. Real L1;
  339. // Example 11:
  340. auto f1 = [](const Real& t) { return 1/(1+t*t);};
  341. Q = gauss<Real, Points>::integrate(f1, 0, boost::math::tools::max_value<Real>(), &L1);
  342. Q_expected = half_pi<Real>();
  343. BOOST_CHECK_CLOSE(Q, Q_expected, 100*tol);
  344. auto f4 = [](const Real& t) { return 1/(1+t*t); };
  345. Q = gauss<Real, Points>::integrate(f4, 1, boost::math::tools::max_value<Real>(), &L1);
  346. Q_expected = pi<Real>()/4;
  347. BOOST_CHECK_CLOSE(Q, Q_expected, 100*tol);
  348. }
  349. template<class Real, unsigned Points>
  350. void test_left_limit_infinite()
  351. {
  352. std::cout << "Testing left limit infinite for Gaussian quadrature in 'A Comparison of Three High Precision Quadrature Schemes' on type " << boost::typeindex::type_id<Real>().pretty_name() << "\n";
  353. Real tol = expected_error<Points>(test_left_limit_infinite_error_id);
  354. Real Q;
  355. Real Q_expected;
  356. // Example 11:
  357. auto f1 = [](const Real& t) { return 1/(1+t*t);};
  358. Q = gauss<Real, Points>::integrate(f1, -boost::math::tools::max_value<Real>(), Real(0));
  359. Q_expected = half_pi<Real>();
  360. BOOST_CHECK_CLOSE(Q, Q_expected, 100*tol);
  361. }
  362. template<class Complex>
  363. void test_complex_lambert_w()
  364. {
  365. std::cout << "Testing that complex-valued integrands are integrated correctly by Gaussian quadrature on type " << boost::typeindex::type_id<Complex>().pretty_name() << "\n";
  366. typedef typename Complex::value_type Real;
  367. Real tol = 10e-9;
  368. using boost::math::constants::pi;
  369. Complex z{2, 3};
  370. auto lw = [&z](Real v)->Complex {
  371. using std::cos;
  372. using std::sin;
  373. using std::exp;
  374. Real sinv = sin(v);
  375. Real cosv = cos(v);
  376. Real cotv = cosv/sinv;
  377. Real cscv = 1/sinv;
  378. Real t = (1-v*cotv)*(1-v*cotv) + v*v;
  379. Real x = v*cscv*exp(-v*cotv);
  380. Complex den = z + x;
  381. Complex num = t*(z/pi<Real>());
  382. Complex res = num/den;
  383. return res;
  384. };
  385. //N[ProductLog[2+3*I], 150]
  386. Complex Q = gauss<Real, 30>::integrate(lw, (Real) 0, pi<Real>());
  387. BOOST_CHECK_CLOSE_FRACTION(Q.real(), boost::lexical_cast<Real>("1.09007653448579084630177782678166964987102108635357778056449870727913321296238687023915522935120701763447787503167111962008709116746523970476893277703"), tol);
  388. BOOST_CHECK_CLOSE_FRACTION(Q.imag(), boost::lexical_cast<Real>("0.530139720774838801426860213574121741928705631382703178297940568794784362495390544411799468140433404536019992695815009036975117285537382995180319280835"), tol);
  389. }
  390. BOOST_AUTO_TEST_CASE(gauss_quadrature_test)
  391. {
  392. #ifdef TEST1
  393. test_linear<double, 7>();
  394. test_quadratic<double, 7>();
  395. test_ca<double, 7>();
  396. test_three_quadrature_schemes_examples<double, 7>();
  397. test_integration_over_real_line<double, 7>();
  398. test_right_limit_infinite<double, 7>();
  399. test_left_limit_infinite<double, 7>();
  400. test_linear<double, 9>();
  401. test_quadratic<double, 9>();
  402. test_ca<double, 9>();
  403. test_three_quadrature_schemes_examples<double, 9>();
  404. test_integration_over_real_line<double, 9>();
  405. test_right_limit_infinite<double, 9>();
  406. test_left_limit_infinite<double, 9>();
  407. test_linear<cpp_bin_float_quad, 10>();
  408. test_quadratic<cpp_bin_float_quad, 10>();
  409. test_ca<cpp_bin_float_quad, 10>();
  410. test_three_quadrature_schemes_examples<cpp_bin_float_quad, 10>();
  411. test_integration_over_real_line<cpp_bin_float_quad, 10>();
  412. test_right_limit_infinite<cpp_bin_float_quad, 10>();
  413. test_left_limit_infinite<cpp_bin_float_quad, 10>();
  414. #endif
  415. #ifdef TEST2
  416. test_linear<cpp_bin_float_quad, 15>();
  417. test_quadratic<cpp_bin_float_quad, 15>();
  418. test_ca<cpp_bin_float_quad, 15>();
  419. test_three_quadrature_schemes_examples<cpp_bin_float_quad, 15>();
  420. test_integration_over_real_line<cpp_bin_float_quad, 15>();
  421. test_right_limit_infinite<cpp_bin_float_quad, 15>();
  422. test_left_limit_infinite<cpp_bin_float_quad, 15>();
  423. test_linear<cpp_bin_float_quad, 20>();
  424. test_quadratic<cpp_bin_float_quad, 20>();
  425. test_ca<cpp_bin_float_quad, 20>();
  426. test_three_quadrature_schemes_examples<cpp_bin_float_quad, 20>();
  427. test_integration_over_real_line<cpp_bin_float_quad, 20>();
  428. test_right_limit_infinite<cpp_bin_float_quad, 20>();
  429. test_left_limit_infinite<cpp_bin_float_quad, 20>();
  430. test_linear<cpp_bin_float_quad, 25>();
  431. test_quadratic<cpp_bin_float_quad, 25>();
  432. test_ca<cpp_bin_float_quad, 25>();
  433. test_three_quadrature_schemes_examples<cpp_bin_float_quad, 25>();
  434. test_integration_over_real_line<cpp_bin_float_quad, 25>();
  435. test_right_limit_infinite<cpp_bin_float_quad, 25>();
  436. test_left_limit_infinite<cpp_bin_float_quad, 25>();
  437. test_linear<cpp_bin_float_quad, 30>();
  438. test_quadratic<cpp_bin_float_quad, 30>();
  439. test_ca<cpp_bin_float_quad, 30>();
  440. test_three_quadrature_schemes_examples<cpp_bin_float_quad, 30>();
  441. test_integration_over_real_line<cpp_bin_float_quad, 30>();
  442. test_right_limit_infinite<cpp_bin_float_quad, 30>();
  443. test_left_limit_infinite<cpp_bin_float_quad, 30>();
  444. #endif
  445. #ifdef TEST3
  446. test_left_limit_infinite<cpp_bin_float_quad, 30>();
  447. test_complex_lambert_w<std::complex<double>>();
  448. test_complex_lambert_w<std::complex<long double>>();
  449. #ifdef BOOST_HAS_FLOAT128
  450. test_left_limit_infinite<boost::multiprecision::float128, 30>();
  451. test_complex_lambert_w<boost::multiprecision::complex128>();
  452. #endif
  453. test_complex_lambert_w<boost::multiprecision::cpp_complex_quad>();
  454. #endif
  455. }
  456. #else
  457. int main() { return 0; }
  458. #endif