tanh_sinh_quadrature_test.cpp 40 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893894895896897898899900901902903904905906907908909910911912913914915916917918919920921922923924925926927928929930931932933934935936937938939940941942943944945946947948949950951952953954955956957958959960961962963964965966967968969970971972973974975976977978979980981982983984985986987988989
  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 <boost/config.hpp>
  8. #include <boost/detail/workaround.hpp>
  9. #if !defined(BOOST_NO_CXX11_DECLTYPE) && !defined(BOOST_NO_CXX11_TRAILING_RESULT_TYPES) && !defined(BOOST_NO_SFINAE_EXPR)
  10. #include <boost/math/concepts/real_concept.hpp>
  11. #include <boost/test/included/unit_test.hpp>
  12. #include <boost/test/tools/floating_point_comparison.hpp>
  13. #include <boost/math/quadrature/tanh_sinh.hpp>
  14. #include <boost/math/special_functions/sinc.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/beta.hpp>
  20. #include <boost/math/special_functions/ellint_rc.hpp>
  21. #include <boost/math/special_functions/ellint_rj.hpp>
  22. #ifdef BOOST_HAS_FLOAT128
  23. #include <boost/multiprecision/float128.hpp>
  24. #endif
  25. #ifdef _MSC_VER
  26. #pragma warning(disable:4127) // Conditional expression is constant
  27. #endif
  28. #if !defined(TEST1) && !defined(TEST2) && !defined(TEST3) && !defined(TEST4) && !defined(TEST5) && !defined(TEST6) && !defined(TEST7) && !defined(TEST8)\
  29. && !defined(TEST1A) && !defined(TEST1B) && !defined(TEST2A) && !defined(TEST3A) && !defined(TEST6A) && !defined(TEST9)
  30. # define TEST1
  31. # define TEST2
  32. # define TEST3
  33. # define TEST4
  34. # define TEST5
  35. # define TEST6
  36. # define TEST7
  37. # define TEST8
  38. # define TEST1A
  39. # define TEST1B
  40. # define TEST2A
  41. # define TEST3A
  42. # define TEST6A
  43. # define TEST9
  44. #endif
  45. using std::expm1;
  46. using std::atan;
  47. using std::tan;
  48. using std::log;
  49. using std::log1p;
  50. using std::asinh;
  51. using std::atanh;
  52. using std::sqrt;
  53. using std::isnormal;
  54. using std::abs;
  55. using std::sinh;
  56. using std::tanh;
  57. using std::cosh;
  58. using std::pow;
  59. using std::exp;
  60. using std::sin;
  61. using std::cos;
  62. using std::string;
  63. using boost::multiprecision::cpp_bin_float_50;
  64. using boost::multiprecision::cpp_bin_float_100;
  65. using boost::multiprecision::cpp_dec_float_50;
  66. using boost::multiprecision::cpp_dec_float_100;
  67. using boost::multiprecision::cpp_bin_float_quad;
  68. using boost::math::sinc_pi;
  69. using boost::math::quadrature::tanh_sinh;
  70. using boost::math::quadrature::detail::tanh_sinh_detail;
  71. using boost::math::constants::pi;
  72. using boost::math::constants::half_pi;
  73. using boost::math::constants::two_div_pi;
  74. using boost::math::constants::two_pi;
  75. using boost::math::constants::half;
  76. using boost::math::constants::third;
  77. using boost::math::constants::half;
  78. using boost::math::constants::third;
  79. using boost::math::constants::catalan;
  80. using boost::math::constants::ln_two;
  81. using boost::math::constants::root_two;
  82. using boost::math::constants::root_two_pi;
  83. using boost::math::constants::root_pi;
  84. template <class T>
  85. void print_levels(const T& v, const char* suffix)
  86. {
  87. std::cout << "{\n";
  88. for (unsigned i = 0; i < v.size(); ++i)
  89. {
  90. std::cout << " { ";
  91. for (unsigned j = 0; j < v[i].size(); ++j)
  92. {
  93. std::cout << v[i][j] << suffix << ", ";
  94. }
  95. std::cout << "},\n";
  96. }
  97. std::cout << " };\n";
  98. }
  99. template <class T>
  100. void print_complement_indexes(const T& v)
  101. {
  102. std::cout << "\n {";
  103. for (unsigned i = 0; i < v.size(); ++i)
  104. {
  105. unsigned index = 0;
  106. while (v[i][index] >= 0)
  107. ++index;
  108. std::cout << index << ", ";
  109. }
  110. std::cout << "};\n";
  111. }
  112. template <class T>
  113. void print_levels(const std::pair<T, T>& p, const char* suffix = "")
  114. {
  115. std::cout << " static const std::vector<std::vector<Real> > abscissa = ";
  116. print_levels(p.first, suffix);
  117. std::cout << " static const std::vector<std::vector<Real> > weights = ";
  118. print_levels(p.second, suffix);
  119. std::cout << " static const std::vector<unsigned > indexes = ";
  120. print_complement_indexes(p.first);
  121. }
  122. template <class Real>
  123. std::pair<std::vector<std::vector<Real>>, std::vector<std::vector<Real>> > generate_constants(unsigned max_index, unsigned max_rows)
  124. {
  125. using boost::math::constants::half_pi;
  126. using boost::math::constants::two_div_pi;
  127. using boost::math::constants::pi;
  128. auto g = [](Real t) { return tanh(half_pi<Real>()*sinh(t)); };
  129. auto w = [](Real t) { Real cs = cosh(half_pi<Real>() * sinh(t)); return half_pi<Real>() * cosh(t) / (cs * cs); };
  130. auto gc = [](Real t) { Real u2 = half_pi<Real>() * sinh(t); return 1 / (exp(u2) *cosh(u2)); };
  131. auto g_inv = [](float x)->float { return asinh(two_div_pi<float>()*atanh(x)); };
  132. auto gc_inv = [](float x)
  133. {
  134. float l = log(sqrt((2 - x) / x));
  135. return log((sqrt(4 * l * l + pi<float>() * pi<float>()) + 2 * l) / pi<float>());
  136. };
  137. std::vector<std::vector<Real>> abscissa, weights;
  138. std::vector<Real> temp;
  139. float t_crossover = gc_inv(0.5f);
  140. Real h = 1;
  141. for (unsigned i = 0; i < max_index; ++i)
  142. {
  143. temp.push_back((i < t_crossover) ? g(i * h) : -gc(i * h));
  144. }
  145. abscissa.push_back(temp);
  146. temp.clear();
  147. for (unsigned i = 0; i < max_index; ++i)
  148. {
  149. temp.push_back(w(i * h));
  150. }
  151. weights.push_back(temp);
  152. temp.clear();
  153. for (unsigned row = 1; row < max_rows; ++row)
  154. {
  155. h /= 2;
  156. for (Real t = h; t < max_index - 1; t += 2 * h)
  157. temp.push_back((t < t_crossover) ? g(t) : -gc(t));
  158. abscissa.push_back(temp);
  159. temp.clear();
  160. }
  161. h = 1;
  162. for (unsigned row = 1; row < max_rows; ++row)
  163. {
  164. h /= 2;
  165. for (Real t = h; t < max_index - 1; t += 2 * h)
  166. temp.push_back(w(t));
  167. weights.push_back(temp);
  168. temp.clear();
  169. }
  170. return std::make_pair(abscissa, weights);
  171. }
  172. template <class Real>
  173. const tanh_sinh<Real>& get_integrator()
  174. {
  175. static const tanh_sinh<Real> integrator(15);
  176. return integrator;
  177. }
  178. template <class Real>
  179. Real get_convergence_tolerance()
  180. {
  181. return boost::math::tools::root_epsilon<Real>();
  182. }
  183. template<class Real>
  184. void test_linear()
  185. {
  186. std::cout << "Testing linear functions are integrated properly by tanh_sinh on type " << boost::typeindex::type_id<Real>().pretty_name() << "\n";
  187. Real tol = 10*boost::math::tools::epsilon<Real>();
  188. auto integrator = get_integrator<Real>();
  189. auto f = [](const Real& x)
  190. {
  191. return 5*x + 7;
  192. };
  193. Real error;
  194. Real L1;
  195. Real Q = integrator.integrate(f, (Real) 0, (Real) 1, get_convergence_tolerance<Real>(), &error, &L1);
  196. BOOST_CHECK_CLOSE_FRACTION(Q, 9.5, tol);
  197. BOOST_CHECK_CLOSE_FRACTION(L1, 9.5, tol);
  198. }
  199. template<class Real>
  200. void test_quadratic()
  201. {
  202. std::cout << "Testing quadratic functions are integrated properly by tanh_sinh on type " << boost::typeindex::type_id<Real>().pretty_name() << "\n";
  203. Real tol = 10*boost::math::tools::epsilon<Real>();
  204. auto integrator = get_integrator<Real>();
  205. auto f = [](const Real& x) { return 5*x*x + 7*x + 12; };
  206. Real error;
  207. Real L1;
  208. Real Q = integrator.integrate(f, (Real) 0, (Real) 1, get_convergence_tolerance<Real>(), &error, &L1);
  209. BOOST_CHECK_CLOSE_FRACTION(Q, (Real) 17 + half<Real>()*third<Real>(), tol);
  210. BOOST_CHECK_CLOSE_FRACTION(L1, (Real) 17 + half<Real>()*third<Real>(), tol);
  211. }
  212. template<class Real>
  213. void test_singular()
  214. {
  215. std::cout << "Testing integration of endpoint singularities on type " << boost::typeindex::type_id<Real>().pretty_name() << "\n";
  216. Real tol = 10*boost::math::tools::epsilon<Real>();
  217. Real error;
  218. Real L1;
  219. auto integrator = get_integrator<Real>();
  220. auto f = [](const Real& x) { return log(x)*log(1-x); };
  221. Real Q = integrator.integrate(f, (Real) 0, (Real) 1, get_convergence_tolerance<Real>(), &error, &L1);
  222. Real Q_expected = 2 - pi<Real>()*pi<Real>()*half<Real>()*third<Real>();
  223. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, tol);
  224. BOOST_CHECK_CLOSE_FRACTION(L1, Q_expected, tol);
  225. }
  226. // Examples taken from
  227. //http://crd-legacy.lbl.gov/~dhbailey/dhbpapers/quadrature.pdf
  228. template<class Real>
  229. void test_ca()
  230. {
  231. std::cout << "Testing integration of C(a) on type " << boost::typeindex::type_id<Real>().pretty_name() << "\n";
  232. Real tol = 10 * boost::math::tools::epsilon<Real>();
  233. Real error;
  234. Real L1;
  235. auto integrator = get_integrator<Real>();
  236. auto f1 = [](const Real& x) { return atan(x)/(x*(x*x + 1)) ; };
  237. Real Q = integrator.integrate(f1, (Real) 0, (Real) 1, get_convergence_tolerance<Real>(), &error, &L1);
  238. Real Q_expected = pi<Real>()*ln_two<Real>()/8 + catalan<Real>()*half<Real>();
  239. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, tol);
  240. BOOST_CHECK_CLOSE_FRACTION(L1, Q_expected, tol);
  241. auto f2 = [](Real x)->Real { Real t0 = x*x + 1; Real t1 = sqrt(t0); return atan(t1)/(t0*t1); };
  242. Q = integrator.integrate(f2, (Real) 0 , (Real) 1, get_convergence_tolerance<Real>(), &error, &L1);
  243. Q_expected = pi<Real>()/4 - pi<Real>()/root_two<Real>() + 3*atan(root_two<Real>())/root_two<Real>();
  244. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, tol);
  245. BOOST_CHECK_CLOSE_FRACTION(L1, Q_expected, tol);
  246. auto f5 = [](Real t)->Real { return t*t*log(t)/((t*t - 1)*(t*t*t*t + 1)); };
  247. Q = integrator.integrate(f5, (Real) 0 , (Real) 1);
  248. Q_expected = pi<Real>()*pi<Real>()*(2 - root_two<Real>())/32;
  249. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, tol);
  250. // Oh it suffers on this one.
  251. auto f6 = [](Real t)->Real { Real x = log(t); return x*x; };
  252. Q = integrator.integrate(f6, (Real) 0 , (Real) 1, get_convergence_tolerance<Real>(), &error, &L1);
  253. Q_expected = 2;
  254. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, 50*tol);
  255. BOOST_CHECK_CLOSE_FRACTION(L1, Q_expected, 50*tol);
  256. // Although it doesn't get to the requested tolerance on this integral, the error bounds can be queried and are reasonable:
  257. tol = sqrt(boost::math::tools::epsilon<Real>());
  258. auto f7 = [](const Real& t) { return sqrt(tan(t)); };
  259. Q = integrator.integrate(f7, (Real) 0 , (Real) half_pi<Real>(), get_convergence_tolerance<Real>(), &error, &L1);
  260. Q_expected = pi<Real>()/root_two<Real>();
  261. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, tol);
  262. BOOST_CHECK_CLOSE_FRACTION(L1, Q_expected, tol);
  263. auto f8 = [](const Real& t) { return log(cos(t)); };
  264. Q = integrator.integrate(f8, (Real) 0 , half_pi<Real>(), get_convergence_tolerance<Real>(), &error, &L1);
  265. Q_expected = -pi<Real>()*ln_two<Real>()*half<Real>();
  266. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, tol);
  267. BOOST_CHECK_CLOSE_FRACTION(L1, -Q_expected, tol);
  268. }
  269. template<class Real>
  270. void test_three_quadrature_schemes_examples()
  271. {
  272. std::cout << "Testing integral in 'A Comparison of Three High Precision Quadrature Schemes' on type " << boost::typeindex::type_id<Real>().pretty_name() << "\n";
  273. Real tol = 10 * boost::math::tools::epsilon<Real>();
  274. Real Q;
  275. Real Q_expected;
  276. auto integrator = get_integrator<Real>();
  277. // Example 1:
  278. auto f1 = [](const Real& t) { return t*boost::math::log1p(t); };
  279. Q = integrator.integrate(f1, (Real) 0 , (Real) 1);
  280. Q_expected = half<Real>()*half<Real>();
  281. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, tol);
  282. // Example 2:
  283. auto f2 = [](const Real& t) { return t*t*atan(t); };
  284. Q = integrator.integrate(f2, (Real) 0 , (Real) 1);
  285. Q_expected = (pi<Real>() -2 + 2*ln_two<Real>())/12;
  286. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, 2 * tol);
  287. // Example 3:
  288. auto f3 = [](const Real& t) { return exp(t)*cos(t); };
  289. Q = integrator.integrate(f3, (Real) 0, half_pi<Real>());
  290. Q_expected = boost::math::expm1(half_pi<Real>())*half<Real>();
  291. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, tol);
  292. // Example 4:
  293. auto f4 = [](Real x)->Real { Real t0 = sqrt(x*x + 2); return atan(t0)/(t0*(x*x+1)); };
  294. Q = integrator.integrate(f4, (Real) 0 , (Real) 1);
  295. Q_expected = 5*pi<Real>()*pi<Real>()/96;
  296. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, tol);
  297. // Example 5:
  298. auto f5 = [](const Real& t) { return sqrt(t)*log(t); };
  299. Q = integrator.integrate(f5, (Real) 0 , (Real) 1);
  300. Q_expected = -4/ (Real) 9;
  301. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, tol);
  302. // Example 6:
  303. auto f6 = [](const Real& t) { return sqrt(1 - t*t); };
  304. Q = integrator.integrate(f6, (Real) 0 , (Real) 1);
  305. Q_expected = pi<Real>()/4;
  306. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, tol);
  307. }
  308. template<class Real>
  309. void test_integration_over_real_line()
  310. {
  311. 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";
  312. Real tol = 10 * boost::math::tools::epsilon<Real>();
  313. Real Q;
  314. Real Q_expected;
  315. Real error;
  316. Real L1;
  317. auto integrator = get_integrator<Real>();
  318. auto f1 = [](const Real& t) { return 1/(1+t*t);};
  319. Q = integrator.integrate(f1, std::numeric_limits<Real>::has_infinity ? -std::numeric_limits<Real>::infinity() : -boost::math::tools::max_value<Real>(), std::numeric_limits<Real>::has_infinity ? std::numeric_limits<Real>::infinity() : boost::math::tools::max_value<Real>(), get_convergence_tolerance<Real>(), &error, &L1);
  320. Q_expected = pi<Real>();
  321. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, tol);
  322. BOOST_CHECK_CLOSE_FRACTION(L1, Q_expected, tol);
  323. auto f2 = [](const Real& t) { return exp(-t*t*half<Real>()); };
  324. Q = integrator.integrate(f2, std::numeric_limits<Real>::has_infinity ? -std::numeric_limits<Real>::infinity() : -boost::math::tools::max_value<Real>(), std::numeric_limits<Real>::has_infinity ? std::numeric_limits<Real>::infinity() : boost::math::tools::max_value<Real>(), get_convergence_tolerance<Real>(), &error, &L1);
  325. Q_expected = root_two_pi<Real>();
  326. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, tol * 2);
  327. BOOST_CHECK_CLOSE_FRACTION(L1, Q_expected, tol * 2);
  328. // This test shows how oscillatory integrals are approximated very poorly by this method:
  329. //std::cout << "Testing sinc integral: \n";
  330. //Q = integrator.integrate(boost::math::sinc_pi<Real>, -std::numeric_limits<Real>::infinity(), std::numeric_limits<Real>::infinity(), &error, &L1);
  331. //std::cout << "Error estimate of sinc integral: " << error << std::endl;
  332. //std::cout << "L1 norm of sinc integral " << L1 << std::endl;
  333. //Q_expected = pi<Real>();
  334. //BOOST_CHECK_CLOSE(Q, Q_expected, 100*tol);
  335. auto f4 = [](const Real& t) { return 1/cosh(t);};
  336. Q = integrator.integrate(f4, std::numeric_limits<Real>::has_infinity ? -std::numeric_limits<Real>::infinity() : -boost::math::tools::max_value<Real>(), std::numeric_limits<Real>::has_infinity ? std::numeric_limits<Real>::infinity() : boost::math::tools::max_value<Real>(), get_convergence_tolerance<Real>(), &error, &L1);
  337. Q_expected = pi<Real>();
  338. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, tol);
  339. BOOST_CHECK_CLOSE_FRACTION(L1, Q_expected, tol);
  340. }
  341. template<class Real>
  342. void test_right_limit_infinite()
  343. {
  344. 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";
  345. Real tol = 10 * boost::math::tools::epsilon<Real>();
  346. Real Q;
  347. Real Q_expected;
  348. Real error;
  349. Real L1;
  350. auto integrator = get_integrator<Real>();
  351. // Example 11:
  352. auto f1 = [](const Real& t) { return 1/(1+t*t);};
  353. Q = integrator.integrate(f1, 0, std::numeric_limits<Real>::has_infinity ? std::numeric_limits<Real>::infinity() : boost::math::tools::max_value<Real>(), get_convergence_tolerance<Real>(), &error, &L1);
  354. Q_expected = half_pi<Real>();
  355. BOOST_CHECK_CLOSE(Q, Q_expected, 100*tol);
  356. // Example 12
  357. auto f2 = [](const Real& t) { return exp(-t)/sqrt(t); };
  358. Q = integrator.integrate(f2, 0, std::numeric_limits<Real>::has_infinity ? std::numeric_limits<Real>::infinity() : boost::math::tools::max_value<Real>(), get_convergence_tolerance<Real>(), &error, &L1);
  359. Q_expected = root_pi<Real>();
  360. BOOST_CHECK_CLOSE(Q, Q_expected, 1000*tol);
  361. auto f3 = [](const Real& t) { return exp(-t)*cos(t); };
  362. Q = integrator.integrate(f3, 0, std::numeric_limits<Real>::has_infinity ? std::numeric_limits<Real>::infinity() : boost::math::tools::max_value<Real>(), get_convergence_tolerance<Real>(), &error, &L1);
  363. Q_expected = half<Real>();
  364. BOOST_CHECK_CLOSE(Q, Q_expected, 100*tol);
  365. auto f4 = [](const Real& t) { return 1/(1+t*t); };
  366. 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);
  367. Q_expected = pi<Real>()/4;
  368. BOOST_CHECK_CLOSE(Q, Q_expected, 100*tol);
  369. }
  370. template<class Real>
  371. void test_left_limit_infinite()
  372. {
  373. std::cout << "Testing left limit infinite for tanh_sinh in 'A Comparison of Three High Precision Quadrature Schemes' on type " << boost::typeindex::type_id<Real>().pretty_name() << "\n";
  374. Real tol = 10 * boost::math::tools::epsilon<Real>();
  375. Real Q;
  376. Real Q_expected;
  377. auto integrator = get_integrator<Real>();
  378. // Example 11:
  379. auto f1 = [](const Real& t) { return 1/(1+t*t);};
  380. Q = integrator.integrate(f1, std::numeric_limits<Real>::has_infinity ? -std::numeric_limits<Real>::infinity() : -boost::math::tools::max_value<Real>(), Real(0));
  381. Q_expected = half_pi<Real>();
  382. BOOST_CHECK_CLOSE(Q, Q_expected, 100*tol);
  383. }
  384. // A horrible function taken from
  385. // http://www.chebfun.org/examples/quad/GaussClenCurt.html
  386. template<class Real>
  387. void test_horrible()
  388. {
  389. std::cout << "Testing Trefenthen's horrible integral on type " << boost::typeindex::type_id<Real>().pretty_name() << "\n";
  390. // We only know the integral to double precision, so requesting a higher tolerance doesn't make sense.
  391. Real tol = 10 * std::numeric_limits<float>::epsilon();
  392. Real Q;
  393. Real Q_expected;
  394. Real error;
  395. Real L1;
  396. auto integrator = get_integrator<Real>();
  397. auto f = [](Real x)->Real { return x*sin(2*exp(2*sin(2*exp(2*x) ) ) ); };
  398. Q = integrator.integrate(f, (Real) -1, (Real) 1, get_convergence_tolerance<Real>(), &error, &L1);
  399. // NIntegrate[x*Sin[2*Exp[2*Sin[2*Exp[2*x]]]], {x, -1, 1}, WorkingPrecision -> 130, MaxRecursion -> 100]
  400. Q_expected = boost::lexical_cast<Real>("0.33673283478172753598559003181355241139806404130031017259552729882281");
  401. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, tol);
  402. // Over again without specifying the bounds:
  403. Q = integrator.integrate(f, get_convergence_tolerance<Real>(), &error, &L1);
  404. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, tol);
  405. }
  406. // Some examples of tough integrals from NR, section 4.5.4:
  407. template<class Real>
  408. void test_nr_examples()
  409. {
  410. std::cout << "Testing singular integrals from NR 4.5.4 on type " << boost::typeindex::type_id<Real>().pretty_name() << "\n";
  411. Real tol = 10 * boost::math::tools::epsilon<Real>();
  412. Real Q;
  413. Real Q_expected;
  414. Real error;
  415. Real L1;
  416. auto integrator = get_integrator<Real>();
  417. auto f1 = [](Real x)->Real
  418. {
  419. return (sin(x * half<Real>()) * exp(-x) / x) / sqrt(x);
  420. };
  421. Q = integrator.integrate(f1, 0, std::numeric_limits<Real>::has_infinity ? std::numeric_limits<Real>::infinity() : boost::math::tools::max_value<Real>(), get_convergence_tolerance<Real>(), &error, &L1);
  422. Q_expected = sqrt(pi<Real>()*(sqrt((Real) 5) - 2));
  423. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, 25*tol);
  424. auto f2 = [](Real x)->Real { return pow(x, -(Real) 2/(Real) 7)*exp(-x*x); };
  425. Q = integrator.integrate(f2, 0, std::numeric_limits<Real>::has_infinity ? std::numeric_limits<Real>::infinity() : boost::math::tools::max_value<Real>());
  426. Q_expected = half<Real>()*boost::math::tgamma((Real) 5/ (Real) 14);
  427. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, tol * 6);
  428. }
  429. // Test integrand known to fool some termination schemes:
  430. template<class Real>
  431. void test_early_termination()
  432. {
  433. std::cout << "Testing Clenshaw & Curtis's example of integrand which fools termination schemes on type " << boost::typeindex::type_id<Real>().pretty_name() << "\n";
  434. Real tol = 10 * boost::math::tools::epsilon<Real>();
  435. Real Q;
  436. Real Q_expected;
  437. Real error;
  438. Real L1;
  439. auto integrator = get_integrator<Real>();
  440. auto f1 = [](Real x)->Real { return 23*cosh(x)/ (Real) 25 - cos(x) ; };
  441. Q = integrator.integrate(f1, (Real) -1, (Real) 1, get_convergence_tolerance<Real>(), &error, &L1);
  442. Q_expected = 46*sinh((Real) 1)/(Real) 25 - 2*sin((Real) 1);
  443. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, tol);
  444. // Over again with no bounds:
  445. Q = integrator.integrate(f1);
  446. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, tol);
  447. }
  448. // Test some definite integrals from the CRC handbook, 32nd edition:
  449. template<class Real>
  450. void test_crc()
  451. {
  452. std::cout << "Testing CRC formulas on type " << boost::typeindex::type_id<Real>().pretty_name() << "\n";
  453. Real tol = 10 * boost::math::tools::epsilon<Real>();
  454. Real Q;
  455. Real Q_expected;
  456. Real error;
  457. Real L1;
  458. auto integrator = get_integrator<Real>();
  459. // CRC Definite integral 585
  460. auto f1 = [](Real x)->Real { Real t = log(1/x); return x*x*t*t*t; };
  461. Q = integrator.integrate(f1, (Real) 0, (Real) 1, get_convergence_tolerance<Real>(), &error, &L1);
  462. Q_expected = (Real) 2/ (Real) 27;
  463. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, tol);
  464. // CRC 636:
  465. auto f2 = [](Real x)->Real { return sqrt(cos(x)); };
  466. Q = integrator.integrate(f2, (Real) 0, (Real) half_pi<Real>(), get_convergence_tolerance<Real>(), &error, &L1);
  467. //Q_expected = pow(two_pi<Real>(), 3*half<Real>())/pow(boost::math::tgamma((Real) 1/ (Real) 4), 2);
  468. Q_expected = boost::lexical_cast<Real>("1.198140234735592207439922492280323878227212663215651558263674952946405214143915670835885556489793389375907225");
  469. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, tol);
  470. // CRC Section 5.5, integral 585:
  471. for (int n = 0; n < 3; ++n) {
  472. for (int m = 0; m < 3; ++m) {
  473. auto f = [&](Real x)->Real { return pow(x, Real(m))*pow(log(1/x), Real(n)); };
  474. Q = integrator.integrate(f, (Real) 0, (Real) 1, get_convergence_tolerance<Real>(), &error, &L1);
  475. // Calculation of the tgamma function is not exact, giving spurious failures.
  476. // Casting to cpp_bin_float_100 beforehand fixes most of them.
  477. cpp_bin_float_100 np1 = n + 1;
  478. cpp_bin_float_100 mp1 = m + 1;
  479. Q_expected = boost::lexical_cast<Real>((tgamma(np1)/pow(mp1, np1)).str());
  480. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, tol);
  481. }
  482. }
  483. // CRC Section 5.5, integral 591
  484. // The parameter p allows us to control the strength of the singularity.
  485. // Rapid convergence is not guaranteed for this function, as the branch cut makes it non-analytic on a disk.
  486. // This converges only when our test type has an extended exponent range as all the area of the integral
  487. // occurs so close to 0 (or 1) that we need abscissa values exceptionally small to find it.
  488. // "There's a lot of room at the bottom".
  489. // We also use a 2 argument functor so that 1-x is evaluated accurately:
  490. if (std::numeric_limits<Real>::max_exponent > std::numeric_limits<double>::max_exponent)
  491. {
  492. for (Real p = Real (-0.99); p < 1; p += Real(0.1)) {
  493. auto f = [&](Real x, Real cx)->Real
  494. {
  495. //return pow(x, p) / pow(1 - x, p);
  496. return cx < 0 ? exp(p * (log(x) - boost::math::log1p(-x))) : pow(x, p) / pow(cx, p);
  497. };
  498. Q = integrator.integrate(f, (Real)0, (Real)1, get_convergence_tolerance<Real>(), &error, &L1);
  499. Q_expected = 1 / sinc_pi(p*pi<Real>());
  500. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, 10 * tol);
  501. }
  502. }
  503. // There is an alternative way to evaluate the above integral: by noticing that all the area of the integral
  504. // is near zero for p < 0 and near 1 for p > 0 we can substitute exp(-x) for x and remap the integral to the
  505. // domain (0, INF). Internally we need to expand out the terms and evaluate using logs to avoid spurous overflow,
  506. // this gives us
  507. // for p > 0:
  508. for (Real p = Real(0.99); p > 0; p -= Real(0.1)) {
  509. auto f = [&](Real x)->Real
  510. {
  511. return exp(-x * (1 - p) + p * log(-boost::math::expm1(-x)));
  512. };
  513. Q = integrator.integrate(f, 0, boost::math::tools::max_value<Real>(), get_convergence_tolerance<Real>(), &error, &L1);
  514. Q_expected = 1 / sinc_pi(p*pi<Real>());
  515. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, 10 * tol);
  516. }
  517. // and for p < 1:
  518. for (Real p = Real (-0.99); p < 0; p += Real(0.1)) {
  519. auto f = [&](Real x)->Real
  520. {
  521. return exp(-p * log(-boost::math::expm1(-x)) - (1 + p) * x);
  522. };
  523. Q = integrator.integrate(f, 0, boost::math::tools::max_value<Real>(), get_convergence_tolerance<Real>(), &error, &L1);
  524. Q_expected = 1 / sinc_pi(p*pi<Real>());
  525. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, 10 * tol);
  526. }
  527. // CRC Section 5.5, integral 635
  528. for (int m = 0; m < 10; ++m) {
  529. auto f = [&](Real x)->Real { return 1/(1 + pow(tan(x), m)); };
  530. Q = integrator.integrate(f, (Real) 0, half_pi<Real>(), get_convergence_tolerance<Real>(), &error, &L1);
  531. Q_expected = half_pi<Real>()/2;
  532. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, tol);
  533. }
  534. // CRC Section 5.5, integral 637:
  535. //
  536. // When h gets very close to 1, the strength of the singularity gradually increases until we
  537. // no longer have enough exponent range to evaluate the integral....
  538. // We also have to use the 2-argument functor version of the integrator to avoid
  539. // cancellation error, since the singularity is near PI/2.
  540. //
  541. Real limit = std::numeric_limits<Real>::max_exponent > std::numeric_limits<double>::max_exponent
  542. ? .95f : .85f;
  543. for (Real h = Real(0.01); h < limit; h += Real(0.1)) {
  544. auto f = [&](Real x, Real xc)->Real { return xc > 0 ? pow(1/tan(xc), h) : pow(tan(x), h); };
  545. Q = integrator.integrate(f, (Real) 0, half_pi<Real>(), get_convergence_tolerance<Real>(), &error, &L1);
  546. Q_expected = half_pi<Real>()/cos(h*half_pi<Real>());
  547. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, tol);
  548. }
  549. // CRC Section 5.5, integral 637:
  550. //
  551. // Over again, but with argument transformation, we want:
  552. //
  553. // Integral of tan(x)^h over (0, PI/2)
  554. //
  555. // Note that the bulk of the area is next to the singularity at PI/2,
  556. // so we'll start by replacing x by PI/2 - x, and that tan(PI/2 - x) == 1/tan(x)
  557. // so we now have:
  558. //
  559. // Integral of 1/tan(x)^h over (0, PI/2)
  560. //
  561. // Which is almost the ideal form, except that when h is very close to 1
  562. // we run out of exponent range in evaluating the integral arbitrarily close to 0.
  563. // So now we substitute exp(-x) for x: this stretches out the range of the
  564. // integral to (-log(PI/2), +INF) with the singularity at +INF giving:
  565. //
  566. // Integral of exp(-x)/tan(exp(-x))^h over (-log(PI/2), +INF)
  567. //
  568. // We just need a way to evaluate the function without spurious under/overflow
  569. // in the exp terms. Note that for small x: tan(x) ~= x, so making this
  570. // substitution and evaluating by logs we have:
  571. //
  572. // exp(-x)/tan(exp(-x))^h ~= exp((h - 1) * x) for x > -log(epsion);
  573. //
  574. // Here's how that looks in code:
  575. //
  576. for (Real i = 80; i < 100; ++i) {
  577. Real h = i / 100;
  578. auto f = [&](Real x)->Real
  579. {
  580. if (x > -log(boost::math::tools::epsilon<Real>()))
  581. return exp((h - 1) * x);
  582. else
  583. {
  584. Real et = exp(-x);
  585. // Need to deal with numeric instability causing et to be greater than PI/2:
  586. return et >= boost::math::constants::half_pi<Real>() ? 0 : et * pow(1 / tan(et), h);
  587. }
  588. };
  589. Q = integrator.integrate(f, -log(half_pi<Real>()), boost::math::tools::max_value<Real>(), get_convergence_tolerance<Real>(), &error, &L1);
  590. Q_expected = half_pi<Real>() / cos(h*half_pi<Real>());
  591. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, 5 * tol);
  592. }
  593. // CRC Section 5.5, integral 670:
  594. auto f3 = [](Real x)->Real { return sqrt(log(1/x)); };
  595. Q = integrator.integrate(f3, (Real) 0, (Real) 1, get_convergence_tolerance<Real>(), &error, &L1);
  596. Q_expected = root_pi<Real>()/2;
  597. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, tol);
  598. }
  599. template <class Real>
  600. void test_sf()
  601. {
  602. using std::sqrt;
  603. // Test some special functions that we already know how to evaluate:
  604. std::cout << "Testing special functions on type " << boost::typeindex::type_id<Real>().pretty_name() << "\n";
  605. Real tol = 10 * boost::math::tools::epsilon<Real>();
  606. auto integrator = get_integrator<Real>();
  607. // incomplete beta:
  608. if (std::numeric_limits<Real>::digits10 < 37) // Otherwise too slow
  609. {
  610. Real a(100), b(15);
  611. auto f = [&](Real x)->Real { return boost::math::ibeta_derivative(a, b, x); };
  612. BOOST_CHECK_CLOSE_FRACTION(integrator.integrate(f, 0, Real(0.25)), boost::math::ibeta(100, 15, Real(0.25)), tol * 10);
  613. // Check some really extreme versions:
  614. a = 1000;
  615. b = 500;
  616. BOOST_CHECK_CLOSE_FRACTION(integrator.integrate(f, 0, 1), Real(1), tol * 15);
  617. //
  618. // This is as extreme as we can get in this domain: otherwise the function has all it's
  619. // area so close to zero we never get in there no matter how many levels we go down:
  620. //
  621. a = Real(1) / 15;
  622. b = 30;
  623. BOOST_CHECK_CLOSE_FRACTION(integrator.integrate(f, 0, 1), Real(1), tol * 25);
  624. }
  625. Real x = 2, y = 3, z = 0.5, p = 0.25;
  626. // Elliptic integral RC:
  627. Real Q = integrator.integrate([&](const Real& t)->Real { return 1 / (sqrt(t + x) * (t + y)); }, 0, std::numeric_limits<Real>::has_infinity ? std::numeric_limits<Real>::infinity() : boost::math::tools::max_value<Real>()) / 2;
  628. BOOST_CHECK_CLOSE_FRACTION(Q, boost::math::ellint_rc(x, y), tol);
  629. // Elliptic Integral RJ:
  630. BOOST_CHECK_CLOSE_FRACTION(Real((Real(3) / 2) * integrator.integrate([&](Real t)->Real { return 1 / (sqrt((t + x) * (t + y) * (t + z)) * (t + p)); }, 0, std::numeric_limits<Real>::has_infinity ? std::numeric_limits<Real>::infinity() : boost::math::tools::max_value<Real>())), boost::math::ellint_rj(x, y, z, p), tol);
  631. z = 5.5;
  632. if (std::numeric_limits<Real>::digits10 > 40)
  633. tol *= 200;
  634. else if (!std::numeric_limits<Real>::is_specialized)
  635. tol *= 3;
  636. // tgamma expressed as an integral:
  637. BOOST_CHECK_CLOSE_FRACTION(integrator.integrate([&](Real t)->Real { using std::pow; using std::exp; return t > 10000 ? Real(0) : Real(pow(t, z - 1) * exp(-t)); }, 0, std::numeric_limits<Real>::has_infinity ? std::numeric_limits<Real>::infinity() : boost::math::tools::max_value<Real>()),
  638. boost::math::tgamma(z), tol);
  639. BOOST_CHECK_CLOSE_FRACTION(integrator.integrate([](const Real& t)->Real { using std::exp; return exp(-t*t); }, std::numeric_limits<Real>::has_infinity ? -std::numeric_limits<Real>::infinity() : -boost::math::tools::max_value<Real>(), std::numeric_limits<Real>::has_infinity ? std::numeric_limits<Real>::infinity() : boost::math::tools::max_value<Real>()),
  640. boost::math::constants::root_pi<Real>(), tol);
  641. }
  642. template <class Real>
  643. void test_2_arg()
  644. {
  645. BOOST_MATH_STD_USING
  646. std::cout << "Testing 2 argument functors on type " << boost::typeindex::type_id<Real>().pretty_name() << "\n";
  647. Real tol = 10 * boost::math::tools::epsilon<Real>();
  648. auto integrator = get_integrator<Real>();
  649. //
  650. // There are a whole family of integrals of the general form
  651. // x(1-x)^-N ; N < 1
  652. // which have all the interesting behaviour near the 2 singularities
  653. // and all converge, see: http://www.wolframalpha.com/input/?i=integrate+(x+*+(1-x))%5E-1%2FN+from+0+to+1
  654. //
  655. Real Q = integrator.integrate([&](const Real& t, const Real & tc)->Real
  656. {
  657. return tc < 0 ? 1 / sqrt(t * (1-t)) : 1 / sqrt(t * tc);
  658. }, 0, 1);
  659. BOOST_CHECK_CLOSE_FRACTION(Q, boost::math::constants::pi<Real>(), tol);
  660. Q = integrator.integrate([&](const Real& t, const Real & tc)->Real
  661. {
  662. return tc < 0 ? 1 / boost::math::cbrt(t * (1-t)) : 1 / boost::math::cbrt(t * tc);
  663. }, 0, 1);
  664. BOOST_CHECK_CLOSE_FRACTION(Q, boost::math::pow<2>(boost::math::tgamma(Real(2) / 3)) / boost::math::tgamma(Real(4) / 3), tol * 3);
  665. //
  666. // We can do the same thing with ((1+x)(1-x))^-N ; N < 1
  667. //
  668. Q = integrator.integrate([&](const Real& t, const Real & tc)->Real
  669. {
  670. return t < 0 ? 1 / sqrt(-tc * (1-t)) : 1 / sqrt((t + 1) * tc);
  671. }, -1, 1);
  672. BOOST_CHECK_CLOSE_FRACTION(Q, boost::math::constants::pi<Real>(), tol);
  673. Q = integrator.integrate([&](const Real& t, const Real & tc)->Real
  674. {
  675. return t < 0 ? 1 / sqrt(-tc * (1-t)) : 1 / sqrt((t + 1) * tc);
  676. });
  677. BOOST_CHECK_CLOSE_FRACTION(Q, boost::math::constants::pi<Real>(), tol);
  678. Q = integrator.integrate([&](const Real& t, const Real & tc)->Real
  679. {
  680. return t < 0 ? 1 / boost::math::cbrt(-tc * (1-t)) : 1 / boost::math::cbrt((t + 1) * tc);
  681. }, -1, 1);
  682. BOOST_CHECK_CLOSE_FRACTION(Q, sqrt(boost::math::constants::pi<Real>()) * boost::math::tgamma(Real(2) / 3) / boost::math::tgamma(Real(7) / 6), tol * 10);
  683. Q = integrator.integrate([&](const Real& t, const Real & tc)->Real
  684. {
  685. return t < 0 ? 1 / boost::math::cbrt(-tc * (1-t)) : 1 / boost::math::cbrt((t + 1) * tc);
  686. });
  687. BOOST_CHECK_CLOSE_FRACTION(Q, sqrt(boost::math::constants::pi<Real>()) * boost::math::tgamma(Real(2) / 3) / boost::math::tgamma(Real(7) / 6), tol * 10);
  688. //
  689. // These are taken from above, and do not get to full precision via the single arg functor:
  690. //
  691. auto f7 = [](const Real& t, const Real& tc) { return t < 1 ? sqrt(tan(t)) : sqrt(1/tan(tc)); };
  692. Real error, L1;
  693. Q = integrator.integrate(f7, (Real)0, (Real)half_pi<Real>(), get_convergence_tolerance<Real>(), &error, &L1);
  694. Real Q_expected = pi<Real>() / root_two<Real>();
  695. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, tol);
  696. BOOST_CHECK_CLOSE_FRACTION(L1, Q_expected, tol);
  697. auto f8 = [](const Real& t, const Real& tc) { return t < 1 ? log(cos(t)) : log(sin(tc)); };
  698. Q = integrator.integrate(f8, (Real)0, half_pi<Real>(), get_convergence_tolerance<Real>(), &error, &L1);
  699. Q_expected = -pi<Real>()*ln_two<Real>()*half<Real>();
  700. BOOST_CHECK_CLOSE_FRACTION(Q, Q_expected, tol);
  701. BOOST_CHECK_CLOSE_FRACTION(L1, -Q_expected, tol);
  702. }
  703. template <class Complex>
  704. void test_complex()
  705. {
  706. typedef typename Complex::value_type value_type;
  707. //
  708. // Integral version of the confluent hypergeometric function:
  709. // https://dlmf.nist.gov/13.4#E1
  710. //
  711. value_type tol = 10 * boost::math::tools::epsilon<value_type>();
  712. Complex a(2, 3), b(3, 4), z(0.5, -2);
  713. auto f = [&](value_type t)
  714. {
  715. return exp(z * t) * pow(t, a - value_type(1)) * pow(value_type(1) - t, b - a - value_type(1));
  716. };
  717. auto integrator = get_integrator<value_type>();
  718. auto Q = integrator.integrate(f, value_type(0), value_type(1), get_convergence_tolerance<value_type>());
  719. //
  720. // Expected result computed from http://www.wolframalpha.com/input/?i=1F1%5B(2%2B3i),+(3%2B4i);+(0.5-2i)%5D+*+gamma(2%2B3i)+*+gamma(1%2Bi)+%2F+gamma(3%2B4i)
  721. //
  722. Complex expected(boost::lexical_cast<value_type>("-0.2911081612888249710582867318081776512805281815037891183828405999609246645054069649838607112484426042883371996"),
  723. boost::lexical_cast<value_type>("0.4507983563969959578849120188097153649211346293694903758252662015991543519595834937475296809912196906074655385"));
  724. value_type error = abs(expected - Q);
  725. BOOST_CHECK_LE(error, tol);
  726. //
  727. // Sin Integral https://dlmf.nist.gov/6.2#E9
  728. //
  729. auto f2 = [z](value_type t)
  730. {
  731. return -exp(-z * cos(t)) * cos(z * sin(t));
  732. };
  733. Q = integrator.integrate(f2, value_type(0), boost::math::constants::half_pi<value_type>(), get_convergence_tolerance<value_type>());
  734. expected = Complex(boost::lexical_cast<value_type>("0.8893822921008980697856313681734926564752476188106405688951257340480164694708337246829840859633322683740376134733"),
  735. -boost::lexical_cast<value_type>("2.381380802906111364088958767973164614925936185337231718483495612539455538280372745733208000514737758457795502168"));
  736. expected -= boost::math::constants::half_pi<value_type>();
  737. error = abs(expected - Q);
  738. BOOST_CHECK_LE(error, tol);
  739. }
  740. BOOST_AUTO_TEST_CASE(tanh_sinh_quadrature_test)
  741. {
  742. #ifdef GENERATE_CONSTANTS
  743. //
  744. // Generate pre-computed coefficients:
  745. std::cout << std::setprecision(35);
  746. print_levels(generate_constants<cpp_bin_float_100>(10, 8), "L");
  747. #else
  748. #ifdef TEST1
  749. test_right_limit_infinite<float>();
  750. test_left_limit_infinite<float>();
  751. test_linear<float>();
  752. test_quadratic<float>();
  753. test_singular<float>();
  754. test_ca<float>();
  755. test_three_quadrature_schemes_examples<float>();
  756. test_horrible<float>();
  757. test_integration_over_real_line<float>();
  758. test_nr_examples<float>();
  759. #endif
  760. #ifdef TEST1A
  761. test_early_termination<float>();
  762. test_2_arg<float>();
  763. #endif
  764. #ifdef TEST1B
  765. test_crc<float>();
  766. #endif
  767. #ifdef TEST2
  768. test_right_limit_infinite<double>();
  769. test_left_limit_infinite<double>();
  770. test_linear<double>();
  771. test_quadratic<double>();
  772. test_singular<double>();
  773. test_ca<double>();
  774. test_three_quadrature_schemes_examples<double>();
  775. test_horrible<double>();
  776. test_integration_over_real_line<double>();
  777. test_nr_examples<double>();
  778. test_early_termination<double>();
  779. test_sf<double>();
  780. test_2_arg<double>();
  781. #endif
  782. #ifdef TEST2A
  783. test_crc<double>();
  784. #endif
  785. #ifndef BOOST_MATH_NO_LONG_DOUBLE_MATH_FUNCTIONS
  786. #ifdef TEST3
  787. test_right_limit_infinite<long double>();
  788. test_left_limit_infinite<long double>();
  789. test_linear<long double>();
  790. test_quadratic<long double>();
  791. test_singular<long double>();
  792. test_ca<long double>();
  793. test_three_quadrature_schemes_examples<long double>();
  794. test_horrible<long double>();
  795. test_integration_over_real_line<long double>();
  796. test_nr_examples<long double>();
  797. test_early_termination<long double>();
  798. test_sf<long double>();
  799. test_2_arg<long double>();
  800. #endif
  801. #ifdef TEST3A
  802. test_crc<long double>();
  803. #endif
  804. #endif
  805. #ifdef TEST4
  806. test_right_limit_infinite<cpp_bin_float_quad>();
  807. test_left_limit_infinite<cpp_bin_float_quad>();
  808. test_linear<cpp_bin_float_quad>();
  809. test_quadratic<cpp_bin_float_quad>();
  810. test_singular<cpp_bin_float_quad>();
  811. test_ca<cpp_bin_float_quad>();
  812. test_three_quadrature_schemes_examples<cpp_bin_float_quad>();
  813. test_horrible<cpp_bin_float_quad>();
  814. test_nr_examples<cpp_bin_float_quad>();
  815. test_early_termination<cpp_bin_float_quad>();
  816. test_crc<cpp_bin_float_quad>();
  817. test_sf<cpp_bin_float_quad>();
  818. test_2_arg<cpp_bin_float_quad>();
  819. #endif
  820. #ifdef TEST5
  821. test_sf<cpp_bin_float_50>();
  822. test_sf<cpp_bin_float_100>();
  823. test_sf<boost::multiprecision::number<boost::multiprecision::cpp_bin_float<150> > >();
  824. #endif
  825. #ifndef BOOST_MATH_NO_LONG_DOUBLE_MATH_FUNCTIONS
  826. #ifdef TEST6
  827. test_right_limit_infinite<boost::math::concepts::real_concept>();
  828. test_left_limit_infinite<boost::math::concepts::real_concept>();
  829. test_linear<boost::math::concepts::real_concept>();
  830. test_quadratic<boost::math::concepts::real_concept>();
  831. test_singular<boost::math::concepts::real_concept>();
  832. test_ca<boost::math::concepts::real_concept>();
  833. test_three_quadrature_schemes_examples<boost::math::concepts::real_concept>();
  834. test_horrible<boost::math::concepts::real_concept>();
  835. test_integration_over_real_line<boost::math::concepts::real_concept>();
  836. test_nr_examples<boost::math::concepts::real_concept>();
  837. test_early_termination<boost::math::concepts::real_concept>();
  838. test_sf<boost::math::concepts::real_concept>();
  839. test_2_arg<boost::math::concepts::real_concept>();
  840. #endif
  841. #ifdef TEST6A
  842. test_crc<boost::math::concepts::real_concept>();
  843. #endif
  844. #endif
  845. #ifdef TEST7
  846. test_sf<cpp_dec_float_50>();
  847. #endif
  848. #if defined(TEST8) && defined(BOOST_HAS_FLOAT128)
  849. test_right_limit_infinite<boost::multiprecision::float128>();
  850. test_left_limit_infinite<boost::multiprecision::float128>();
  851. test_linear<boost::multiprecision::float128>();
  852. test_quadratic<boost::multiprecision::float128>();
  853. test_singular<boost::multiprecision::float128>();
  854. test_ca<boost::multiprecision::float128>();
  855. test_three_quadrature_schemes_examples<boost::multiprecision::float128>();
  856. test_horrible<boost::multiprecision::float128>();
  857. test_integration_over_real_line<boost::multiprecision::float128>();
  858. test_nr_examples<boost::multiprecision::float128>();
  859. test_early_termination<boost::multiprecision::float128>();
  860. test_crc<boost::multiprecision::float128>();
  861. test_sf<boost::multiprecision::float128>();
  862. test_2_arg<boost::multiprecision::float128>();
  863. #endif
  864. #ifdef TEST9
  865. test_complex<std::complex<double> >();
  866. test_complex<std::complex<float> >();
  867. #endif
  868. #endif
  869. }
  870. #else
  871. int main() { return 0; }
  872. #endif