lanczos_smoothing_test.cpp 21 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708
  1. /*
  2. * Copyright Nick Thompson, 2019
  3. * Use, modification and distribution are subject to the
  4. * Boost Software License, Version 1.0. (See accompanying file
  5. * LICENSE_1_0.txt or copy at http://www.boost.org/LICENSE_1_0.txt)
  6. */
  7. #define BOOST_TEST_MODULE lanczos_smoothing_test
  8. #include <random>
  9. #include <array>
  10. #include <boost/range.hpp>
  11. #include <boost/numeric/ublas/vector.hpp>
  12. #include <boost/math/constants/constants.hpp>
  13. #include <boost/test/included/unit_test.hpp>
  14. #include <boost/test/tools/floating_point_comparison.hpp>
  15. #include <boost/math/differentiation/lanczos_smoothing.hpp>
  16. #include <boost/multiprecision/cpp_bin_float.hpp>
  17. #include <boost/math/special_functions/next.hpp> // for float_distance
  18. #include <boost/math/tools/condition_numbers.hpp>
  19. using std::abs;
  20. using std::pow;
  21. using std::sqrt;
  22. using std::sin;
  23. using boost::math::constants::two_pi;
  24. using boost::multiprecision::cpp_bin_float_50;
  25. using boost::multiprecision::cpp_bin_float_100;
  26. using boost::math::differentiation::discrete_lanczos_derivative;
  27. using boost::math::differentiation::detail::discrete_legendre;
  28. using boost::math::differentiation::detail::interior_velocity_filter;
  29. using boost::math::differentiation::detail::boundary_velocity_filter;
  30. using boost::math::tools::summation_condition_number;
  31. template<class Real>
  32. void test_dlp_norms()
  33. {
  34. std::cout << "Testing Discrete Legendre Polynomial norms on type " << typeid(Real).name() << "\n";
  35. Real tol = std::numeric_limits<Real>::epsilon();
  36. auto dlp = discrete_legendre<Real>(1, Real(0));
  37. BOOST_CHECK_CLOSE_FRACTION(dlp.norm_sq(0), 3, tol);
  38. BOOST_CHECK_CLOSE_FRACTION(dlp.norm_sq(1), 2, tol);
  39. dlp = discrete_legendre<Real>(2, Real(0));
  40. BOOST_CHECK_CLOSE_FRACTION(dlp.norm_sq(0), Real(5)/Real(2), tol);
  41. BOOST_CHECK_CLOSE_FRACTION(dlp.norm_sq(1), Real(5)/Real(4), tol);
  42. BOOST_CHECK_CLOSE_FRACTION(dlp.norm_sq(2), Real(3*3*7)/Real(pow(2,6)), 2*tol);
  43. dlp = discrete_legendre<Real>(200, Real(0));
  44. for(size_t r = 0; r < 10; ++r)
  45. {
  46. Real calc = dlp.norm_sq(r);
  47. Real expected = Real(2)/Real(2*r+1);
  48. // As long as r << n, ||q_r||^2 -> 2/(2r+1) as n->infty
  49. BOOST_CHECK_CLOSE_FRACTION(calc, expected, 0.05);
  50. }
  51. }
  52. template<class Real>
  53. void test_dlp_evaluation()
  54. {
  55. std::cout << "Testing evaluation of Discrete Legendre polynomials on type " << typeid(Real).name() << "\n";
  56. Real tol = std::numeric_limits<Real>::epsilon();
  57. size_t n = 25;
  58. Real x = 0.72;
  59. auto dlp = discrete_legendre<Real>(n, x);
  60. Real q0 = dlp(x, 0);
  61. BOOST_TEST(q0 == 1);
  62. Real q1 = dlp(x, 1);
  63. BOOST_TEST(q1 == x);
  64. Real q2 = dlp(x, 2);
  65. int N = 2*n+1;
  66. Real expected = 0.5*(3*x*x - Real(N*N - 1)/Real(4*n*n));
  67. BOOST_CHECK_CLOSE_FRACTION(q2, expected, tol);
  68. Real q3 = dlp(x, 3);
  69. expected = (x/3)*(5*expected - (Real(N*N - 4))/(2*n*n));
  70. BOOST_CHECK_CLOSE_FRACTION(q3, expected, 2*tol);
  71. // q_r(x) is even for even r, and odd for odd r:
  72. for (size_t n = 8; n < 22; ++n)
  73. {
  74. dlp = discrete_legendre<Real>(n, x);
  75. for(size_t r = 2; r <= n; ++r)
  76. {
  77. if (r & 1)
  78. {
  79. Real q1 = dlp(x, r);
  80. Real q2 = -dlp(-x, r);
  81. BOOST_CHECK_CLOSE_FRACTION(q1, q2, tol);
  82. }
  83. else
  84. {
  85. Real q1 = dlp(x, r);
  86. Real q2 = dlp(-x, r);
  87. BOOST_CHECK_CLOSE_FRACTION(q1, q2, tol);
  88. }
  89. Real l2_sq = 0;
  90. for (int j = -(int)n; j <= (int) n; ++j)
  91. {
  92. Real y = Real(j)/Real(n);
  93. Real term = dlp(y, r);
  94. l2_sq += term*term;
  95. }
  96. l2_sq /= n;
  97. Real l2_sq_expected = dlp.norm_sq(r);
  98. BOOST_CHECK_CLOSE_FRACTION(l2_sq, l2_sq_expected, 20*tol);
  99. }
  100. }
  101. }
  102. template<class Real>
  103. void test_dlp_next()
  104. {
  105. std::cout << "Testing Discrete Legendre polynomial 'next' function on type " << typeid(Real).name() << "\n";
  106. Real tol = std::numeric_limits<Real>::epsilon();
  107. for(size_t n = 2; n < 20; ++n)
  108. {
  109. for(Real x = -1; x <= 1; x += 0.1)
  110. {
  111. auto dlp = discrete_legendre<Real>(n, x);
  112. for (size_t k = 2; k < n; ++k)
  113. {
  114. BOOST_CHECK_CLOSE(dlp.next(), dlp(x, k), tol);
  115. }
  116. dlp = discrete_legendre<Real>(n, x);
  117. for (size_t k = 2; k < n; ++k)
  118. {
  119. BOOST_CHECK_CLOSE(dlp.next_prime(), dlp.prime(x, k), tol);
  120. }
  121. }
  122. }
  123. }
  124. template<class Real>
  125. void test_dlp_derivatives()
  126. {
  127. std::cout << "Testing Discrete Legendre polynomial derivatives on type " << typeid(Real).name() << "\n";
  128. Real tol = 10*std::numeric_limits<Real>::epsilon();
  129. int n = 25;
  130. Real x = 0.72;
  131. auto dlp = discrete_legendre<Real>(n, x);
  132. Real q0p = dlp.prime(x, 0);
  133. BOOST_TEST(q0p == 0);
  134. Real q1p = dlp.prime(x, 1);
  135. BOOST_TEST(q1p == 1);
  136. Real q2p = dlp.prime(x, 2);
  137. Real expected = 3*x;
  138. BOOST_CHECK_CLOSE_FRACTION(q2p, expected, tol);
  139. }
  140. template<class Real>
  141. void test_dlp_second_derivative()
  142. {
  143. std::cout << "Testing Discrete Legendre polynomial derivatives on type " << typeid(Real).name() << "\n";
  144. int n = 25;
  145. Real x = Real(1)/Real(3);
  146. auto dlp = discrete_legendre<Real>(n, x);
  147. Real q2pp = dlp.next_dbl_prime();
  148. BOOST_TEST(q2pp == 3);
  149. }
  150. template<class Real>
  151. void test_interior_velocity_filter()
  152. {
  153. using boost::math::constants::half;
  154. std::cout << "Testing interior filter on type " << typeid(Real).name() << "\n";
  155. Real tol = std::numeric_limits<Real>::epsilon();
  156. for(int n = 1; n < 10; ++n)
  157. {
  158. for (int p = 1; p < n; p += 2)
  159. {
  160. auto f = interior_velocity_filter<Real>(n,p);
  161. // Since we only store half the filter coefficients,
  162. // we need to reindex the moment sums:
  163. auto cond = summation_condition_number<Real>(0);
  164. for (size_t j = 0; j < f.size(); ++j)
  165. {
  166. cond += j*f[j];
  167. }
  168. BOOST_CHECK_CLOSE_FRACTION(cond.sum(), half<Real>(), 2*cond()*tol);
  169. for (int l = 3; l <= p; l += 2)
  170. {
  171. cond = summation_condition_number<Real>(0);
  172. for (size_t j = 0; j < f.size() - 1; ++j)
  173. {
  174. cond += pow(Real(j), l)*f[j];
  175. }
  176. Real expected = -pow(Real(f.size() - 1), l)*f[f.size()-1];
  177. BOOST_CHECK_CLOSE_FRACTION(expected, cond.sum(), 7*cond()*tol);
  178. }
  179. //std::cout << "(n,p) = (" << n << "," << p << ") = {";
  180. //for (auto & x : f)
  181. //{
  182. // std::cout << x << ", ";
  183. //}
  184. //std::cout << "}\n";
  185. }
  186. }
  187. }
  188. template<class Real>
  189. void test_interior_lanczos()
  190. {
  191. std::cout << "Testing interior Lanczos on type " << typeid(Real).name() << "\n";
  192. Real tol = std::numeric_limits<Real>::epsilon();
  193. std::vector<Real> v(500);
  194. std::fill(v.begin(), v.end(), 7);
  195. for (size_t n = 1; n < 10; ++n)
  196. {
  197. for (size_t p = 2; p < 2*n; p += 2)
  198. {
  199. auto dld = discrete_lanczos_derivative(Real(0.1), n, p);
  200. for (size_t m = n; m < v.size() - n; ++m)
  201. {
  202. Real dvdt = dld(v, m);
  203. BOOST_CHECK_SMALL(dvdt, tol);
  204. }
  205. auto dvdt = dld(v);
  206. for (size_t m = n; m < v.size() - n; ++m)
  207. {
  208. BOOST_CHECK_SMALL(dvdt[m], tol);
  209. }
  210. }
  211. }
  212. for(size_t i = 0; i < v.size(); ++i)
  213. {
  214. v[i] = 7*i+8;
  215. }
  216. for (size_t n = 1; n < 10; ++n)
  217. {
  218. for (size_t p = 2; p < 2*n; p += 2)
  219. {
  220. auto dld = discrete_lanczos_derivative(Real(1), n, p);
  221. for (size_t m = n; m < v.size() - n; ++m)
  222. {
  223. Real dvdt = dld(v, m);
  224. BOOST_CHECK_CLOSE_FRACTION(dvdt, 7, 2000*tol);
  225. }
  226. auto dvdt = dld(v);
  227. for (size_t m = n; m < v.size() - n; ++m)
  228. {
  229. BOOST_CHECK_CLOSE_FRACTION(dvdt[m], 7, 2000*tol);
  230. }
  231. }
  232. }
  233. //std::random_device rd{};
  234. //auto seed = rd();
  235. //std::cout << "Seed = " << seed << "\n";
  236. std::mt19937 gen(4172378669);
  237. std::normal_distribution<> dis{0, 0.01};
  238. for (size_t i = 0; i < v.size(); ++i)
  239. {
  240. v[i] = 7*i+8 + dis(gen);
  241. }
  242. for (size_t n = 1; n < 10; ++n)
  243. {
  244. for (size_t p = 2; p < 2*n; p += 2)
  245. {
  246. auto dld = discrete_lanczos_derivative(Real(1), n, p);
  247. for (size_t m = n; m < v.size() - n; ++m)
  248. {
  249. BOOST_CHECK_CLOSE_FRACTION(dld(v, m), Real(7), Real(0.0042));
  250. }
  251. }
  252. }
  253. for (size_t i = 0; i < v.size(); ++i)
  254. {
  255. v[i] = 15*i*i + 7*i+8 + dis(gen);
  256. }
  257. for (size_t n = 1; n < 10; ++n)
  258. {
  259. for (size_t p = 2; p < 2*n; p += 2)
  260. {
  261. auto dld = discrete_lanczos_derivative(Real(1), n, p);
  262. for (size_t m = n; m < v.size() - n; ++m)
  263. {
  264. BOOST_CHECK_CLOSE_FRACTION(dld(v,m), Real(30*m + 7), Real(0.00008));
  265. }
  266. }
  267. }
  268. std::normal_distribution<> dis1{0, 0.0001};
  269. Real omega = Real(1)/Real(16);
  270. for (size_t i = 0; i < v.size(); ++i)
  271. {
  272. v[i] = sin(i*omega) + dis1(gen);
  273. }
  274. for (size_t n = 10; n < 20; ++n)
  275. {
  276. for (size_t p = 3; p < 100 && p < n/2; p += 2)
  277. {
  278. auto dld = discrete_lanczos_derivative(Real(1), n, p);
  279. for (size_t m = n; m < v.size() - n && m < n + 10; ++m)
  280. {
  281. BOOST_CHECK_CLOSE_FRACTION(dld(v,m), omega*cos(omega*m), Real(0.03));
  282. }
  283. }
  284. }
  285. }
  286. template<class Real>
  287. void test_boundary_velocity_filters()
  288. {
  289. std::cout << "Testing boundary filters on type " << typeid(Real).name() << "\n";
  290. Real tol = std::numeric_limits<Real>::epsilon();
  291. for(int n = 1; n < 5; ++n)
  292. {
  293. for (int p = 1; p < 2*n+1; ++p)
  294. {
  295. for (int s = -n; s <= n; ++s)
  296. {
  297. auto f = boundary_velocity_filter<Real>(n, p, s);
  298. // Sum is zero:
  299. auto cond = summation_condition_number<Real>(0);
  300. for (size_t i = 0; i < f.size() - 1; ++i)
  301. {
  302. cond += f[i];
  303. }
  304. BOOST_CHECK_CLOSE_FRACTION(cond.sum(), -f[f.size()-1], 6*cond()*tol);
  305. cond = summation_condition_number<Real>(0);
  306. for (size_t k = 0; k < f.size(); ++k)
  307. {
  308. Real j = Real(k) - Real(n);
  309. // note the shifted index here:
  310. cond += (j-s)*f[k];
  311. }
  312. BOOST_CHECK_CLOSE_FRACTION(cond.sum(), 1, 6*cond()*tol);
  313. for (int l = 2; l <= p; ++l)
  314. {
  315. cond = summation_condition_number<Real>(0);
  316. for (size_t k = 0; k < f.size() - 1; ++k)
  317. {
  318. Real j = Real(k) - Real(n);
  319. // The condition number of this sum is infinite!
  320. // No need to get to worked up about the tolerance.
  321. cond += pow(j-s, l)*f[k];
  322. }
  323. Real expected = -pow(Real(f.size()-1) - Real(n) - Real(s), l)*f[f.size()-1];
  324. if (expected == 0)
  325. {
  326. BOOST_CHECK_SMALL(cond.sum(), cond()*tol);
  327. }
  328. else
  329. {
  330. BOOST_CHECK_CLOSE_FRACTION(expected, cond.sum(), 200*cond()*tol);
  331. }
  332. }
  333. //std::cout << "(n,p,s) = ("<< n << ", " << p << "," << s << ") = {";
  334. //for (auto & x : f)
  335. //{
  336. // std::cout << x << ", ";
  337. //}
  338. //std::cout << "}\n";*/
  339. }
  340. }
  341. }
  342. }
  343. template<class Real>
  344. void test_boundary_lanczos()
  345. {
  346. std::cout << "Testing Lanczos boundary on type " << typeid(Real).name() << "\n";
  347. Real tol = std::numeric_limits<Real>::epsilon();
  348. std::vector<Real> v(500, 7);
  349. for (size_t n = 1; n < 10; ++n)
  350. {
  351. for (size_t p = 2; p < 2*n; ++p)
  352. {
  353. auto lsd = discrete_lanczos_derivative(Real(0.0125), n, p);
  354. for (size_t m = 0; m < n; ++m)
  355. {
  356. Real dvdt = lsd(v,m);
  357. BOOST_CHECK_SMALL(dvdt, 4*sqrt(tol));
  358. }
  359. for (size_t m = v.size() - n; m < v.size(); ++m)
  360. {
  361. Real dvdt = lsd(v,m);
  362. BOOST_CHECK_SMALL(dvdt, 4*sqrt(tol));
  363. }
  364. }
  365. }
  366. for(size_t i = 0; i < v.size(); ++i)
  367. {
  368. v[i] = 7*i+8;
  369. }
  370. for (size_t n = 3; n < 10; ++n)
  371. {
  372. for (size_t p = 2; p < 2*n; ++p)
  373. {
  374. auto lsd = discrete_lanczos_derivative(Real(1), n, p);
  375. for (size_t m = 0; m < n; ++m)
  376. {
  377. Real dvdt = lsd(v,m);
  378. BOOST_CHECK_CLOSE_FRACTION(dvdt, 7, sqrt(tol));
  379. }
  380. for (size_t m = v.size() - n; m < v.size(); ++m)
  381. {
  382. Real dvdt = lsd(v,m);
  383. BOOST_CHECK_CLOSE_FRACTION(dvdt, 7, 4*sqrt(tol));
  384. }
  385. }
  386. }
  387. for (size_t i = 0; i < v.size(); ++i)
  388. {
  389. v[i] = 15*i*i + 7*i+8;
  390. }
  391. for (size_t n = 1; n < 10; ++n)
  392. {
  393. for (size_t p = 2; p < 2*n; ++p)
  394. {
  395. auto lsd = discrete_lanczos_derivative(Real(1), n, p);
  396. for (size_t m = 0; m < v.size(); ++m)
  397. {
  398. BOOST_CHECK_CLOSE_FRACTION(lsd(v,m), 30*m+7, 30*sqrt(tol));
  399. }
  400. }
  401. }
  402. // Demonstrate that the boundary filters are also denoising:
  403. //std::random_device rd{};
  404. //auto seed = rd();
  405. //std::cout << "seed = " << seed << "\n";
  406. std::mt19937 gen(311354333);
  407. std::normal_distribution<> dis{0, 0.01};
  408. for (size_t i = 0; i < v.size(); ++i)
  409. {
  410. v[i] += dis(gen);
  411. }
  412. for (size_t n = 1; n < 10; ++n)
  413. {
  414. for (size_t p = 2; p < n; ++p)
  415. {
  416. auto lsd = discrete_lanczos_derivative(Real(1), n, p);
  417. for (size_t m = 0; m < v.size(); ++m)
  418. {
  419. BOOST_CHECK_CLOSE_FRACTION(lsd(v,m), 30*m+7, 0.005);
  420. }
  421. auto dvdt = lsd(v);
  422. for (size_t m = 0; m < v.size(); ++m)
  423. {
  424. BOOST_CHECK_CLOSE_FRACTION(dvdt[m], 30*m+7, 0.005);
  425. }
  426. }
  427. }
  428. }
  429. template<class Real>
  430. void test_acceleration_filters()
  431. {
  432. Real eps = std::numeric_limits<Real>::epsilon();
  433. for (size_t n = 1; n < 5; ++n)
  434. {
  435. for(size_t p = 3; p <= 2*n; ++p)
  436. {
  437. for(int64_t s = -int64_t(n); s <= 0; ++s)
  438. {
  439. auto g = boost::math::differentiation::detail::acceleration_filter<long double>(n,p,s);
  440. std::vector<Real> f(g.size());
  441. for (size_t i = 0; i < g.size(); ++i)
  442. {
  443. f[i] = static_cast<Real>(g[i]);
  444. }
  445. auto cond = summation_condition_number<Real>(0);
  446. for (size_t i = 0; i < f.size() - 1; ++i)
  447. {
  448. cond += f[i];
  449. }
  450. BOOST_CHECK_CLOSE_FRACTION(cond.sum(), -f[f.size()-1], 10*cond()*eps);
  451. cond = summation_condition_number<Real>(0);
  452. for (size_t k = 0; k < f.size() -1; ++k)
  453. {
  454. Real j = Real(k) - Real(n);
  455. cond += (j-s)*f[k];
  456. }
  457. Real expected = -(Real(f.size()-1)- Real(n) - s)*f[f.size()-1];
  458. BOOST_CHECK_CLOSE_FRACTION(cond.sum(), expected, 10*cond()*eps);
  459. cond = summation_condition_number<Real>(0);
  460. for (size_t k = 0; k < f.size(); ++k)
  461. {
  462. Real j = Real(k) - Real(n);
  463. cond += (j-s)*(j-s)*f[k];
  464. }
  465. BOOST_CHECK_CLOSE_FRACTION(cond.sum(), 2, 100*cond()*eps);
  466. // See unlabelled equation in McDevitt, 2012, just after equation 26:
  467. // It appears that there is an off-by-one error in that equation, since p + 1 moments don't vanish, only p.
  468. // This test is itself suspect; the condition number of the moment sum is infinite.
  469. // So the *slightest* error in the filter gets amplified by the test; in terms of the
  470. // behavior of the actual filter, it's not a big deal.
  471. for (size_t l = 3; l <= p; ++l)
  472. {
  473. cond = summation_condition_number<Real>(0);
  474. for (size_t k = 0; k < f.size() - 1; ++k)
  475. {
  476. Real j = Real(k) - Real(n);
  477. cond += pow((j-s), l)*f[k];
  478. }
  479. Real expected = -pow(Real(f.size()- 1 - n -s), l)*f[f.size()-1];
  480. BOOST_CHECK_CLOSE_FRACTION(cond.sum(), expected, 1000*cond()*eps);
  481. }
  482. }
  483. }
  484. }
  485. }
  486. template<class Real>
  487. void test_lanczos_acceleration()
  488. {
  489. Real eps = std::numeric_limits<Real>::epsilon();
  490. std::vector<Real> v(100, 7);
  491. auto lanczos = discrete_lanczos_derivative<Real, 2>(Real(1), 4, 3);
  492. for (size_t i = 0; i < v.size(); ++i)
  493. {
  494. BOOST_CHECK_SMALL(lanczos(v, i), eps);
  495. }
  496. for(size_t i = 0; i < v.size(); ++i)
  497. {
  498. v[i] = 7*i + 6;
  499. }
  500. for (size_t i = 0; i < v.size(); ++i)
  501. {
  502. BOOST_CHECK_SMALL(lanczos(v,i), 200*eps);
  503. }
  504. for(size_t i = 0; i < v.size(); ++i)
  505. {
  506. v[i] = 7*i*i + 9*i + 6;
  507. }
  508. for (size_t i = 0; i < v.size(); ++i)
  509. {
  510. BOOST_CHECK_CLOSE_FRACTION(lanczos(v, i), 14, 1500*eps);
  511. }
  512. // Now add noise, and kick up the smoothing of the Lanzcos derivative (increase n):
  513. //std::random_device rd{};
  514. //auto seed = rd();
  515. //std::cout << "seed = " << seed << "\n";
  516. size_t seed = 2507134629;
  517. std::mt19937 gen(seed);
  518. Real std_dev = 0.1;
  519. std::normal_distribution<Real> dis{0, std_dev};
  520. for (size_t i = 0; i < v.size(); ++i)
  521. {
  522. v[i] += dis(gen);
  523. }
  524. lanczos = discrete_lanczos_derivative<Real, 2>(Real(1), 18, 3);
  525. auto w = lanczos(v);
  526. for (size_t i = 0; i < v.size(); ++i)
  527. {
  528. BOOST_CHECK_CLOSE_FRACTION(w[i], 14, std_dev/200);
  529. }
  530. }
  531. template<class Real>
  532. void test_rescaling()
  533. {
  534. std::cout << "Test rescaling on type " << typeid(Real).name() << "\n";
  535. Real tol = std::numeric_limits<Real>::epsilon();
  536. std::vector<Real> v(500);
  537. for(size_t i = 0; i < v.size(); ++i)
  538. {
  539. v[i] = 7*i*i + 9*i + 6;
  540. }
  541. std::vector<Real> dvdt1(500);
  542. std::vector<Real> dvdt2(500);
  543. auto lanczos1 = discrete_lanczos_derivative(Real(1));
  544. auto lanczos2 = discrete_lanczos_derivative(Real(2));
  545. lanczos1(v, dvdt1);
  546. lanczos2(v, dvdt2);
  547. for(size_t i = 0; i < v.size(); ++i)
  548. {
  549. BOOST_CHECK_CLOSE_FRACTION(dvdt1[i], 2*dvdt2[i], tol);
  550. }
  551. auto lanczos3 = discrete_lanczos_derivative<Real, 2>(Real(1));
  552. auto lanczos4 = discrete_lanczos_derivative<Real, 2>(Real(2));
  553. std::vector<Real> dv2dt21(500);
  554. std::vector<Real> dv2dt22(500);
  555. for(size_t i = 0; i < v.size(); ++i)
  556. {
  557. BOOST_CHECK_CLOSE_FRACTION(dv2dt21[i], 4*dv2dt22[i], tol);
  558. }
  559. }
  560. template<class Real>
  561. void test_data_representations()
  562. {
  563. std::cout << "Test rescaling on type " << typeid(Real).name() << "\n";
  564. Real tol = 150*std::numeric_limits<Real>::epsilon();
  565. std::array<Real, 500> v;
  566. for(size_t i = 0; i < v.size(); ++i)
  567. {
  568. v[i] = 9*i + 6;
  569. }
  570. std::array<Real, 500> dvdt;
  571. auto lanczos = discrete_lanczos_derivative(Real(1));
  572. lanczos(v, dvdt);
  573. for(size_t i = 0; i < v.size(); ++i)
  574. {
  575. BOOST_CHECK_CLOSE_FRACTION(dvdt[i], 9, tol);
  576. }
  577. boost::numeric::ublas::vector<Real> w(500);
  578. boost::numeric::ublas::vector<Real> dwdt(500);
  579. for(size_t i = 0; i < w.size(); ++i)
  580. {
  581. w[i] = 9*i + 6;
  582. }
  583. lanczos(w, dwdt);
  584. for(size_t i = 0; i < v.size(); ++i)
  585. {
  586. BOOST_CHECK_CLOSE_FRACTION(dwdt[i], 9, tol);
  587. }
  588. auto v1 = boost::make_iterator_range(v.begin(), v.end());
  589. auto v2 = boost::make_iterator_range(dvdt.begin(), dvdt.end());
  590. lanczos(v1, v2);
  591. for(size_t i = 0; i < v2.size(); ++i)
  592. {
  593. BOOST_CHECK_CLOSE_FRACTION(v2[i], 9, tol);
  594. }
  595. auto lanczos2 = discrete_lanczos_derivative<Real, 2>(Real(1));
  596. lanczos2(v1, v2);
  597. for(size_t i = 0; i < v2.size(); ++i)
  598. {
  599. BOOST_CHECK_SMALL(v2[i], 10*tol);
  600. }
  601. }
  602. BOOST_AUTO_TEST_CASE(lanczos_smoothing_test)
  603. {
  604. test_dlp_second_derivative<double>();
  605. test_dlp_norms<double>();
  606. test_dlp_evaluation<double>();
  607. test_dlp_derivatives<double>();
  608. test_dlp_next<double>();
  609. // Takes too long!
  610. //test_dlp_norms<cpp_bin_float_50>();
  611. test_boundary_velocity_filters<double>();
  612. test_boundary_velocity_filters<long double>();
  613. // Takes too long!
  614. //test_boundary_velocity_filters<cpp_bin_float_50>();
  615. test_boundary_lanczos<double>();
  616. test_boundary_lanczos<long double>();
  617. // Takes too long!
  618. //test_boundary_lanczos<cpp_bin_float_50>();
  619. test_interior_velocity_filter<double>();
  620. test_interior_velocity_filter<long double>();
  621. test_interior_lanczos<double>();
  622. test_acceleration_filters<double>();
  623. test_lanczos_acceleration<float>();
  624. test_lanczos_acceleration<double>();
  625. test_rescaling<double>();
  626. test_data_representations<double>();
  627. }