gauss_kronrod_quadrature_test.cpp 18 KB

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