lanczos_smoothing.hpp 20 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590
  1. // (C) Copyright Nick Thompson 2019.
  2. // Use, modification and distribution are subject to the
  3. // Boost Software License, Version 1.0. (See accompanying file
  4. // LICENSE_1_0.txt or copy at http://www.boost.org/LICENSE_1_0.txt)
  5. #ifndef BOOST_MATH_DIFFERENTIATION_LANCZOS_SMOOTHING_HPP
  6. #define BOOST_MATH_DIFFERENTIATION_LANCZOS_SMOOTHING_HPP
  7. #include <cmath> // for std::abs
  8. #include <cstddef>
  9. #include <limits> // to nan initialize
  10. #include <vector>
  11. #include <string>
  12. #include <stdexcept>
  13. #include <type_traits>
  14. #include <boost/math/tools/assert.hpp>
  15. #include <boost/math/tools/is_standalone.hpp>
  16. #ifndef BOOST_MATH_STANDALONE
  17. #include <boost/config.hpp>
  18. #ifdef BOOST_NO_CXX17_IF_CONSTEXPR
  19. #error "The header <boost/math/norms.hpp> can only be used in C++17 and later."
  20. #endif
  21. #endif
  22. namespace boost::math::differentiation {
  23. namespace detail {
  24. template <typename Real>
  25. class discrete_legendre {
  26. public:
  27. explicit discrete_legendre(std::size_t n, Real x) : m_n{n}, m_r{2}, m_x{x},
  28. m_qrm2{1}, m_qrm1{x},
  29. m_qrm2p{0}, m_qrm1p{1},
  30. m_qrm2pp{0}, m_qrm1pp{0}
  31. {
  32. using std::abs;
  33. BOOST_MATH_ASSERT_MSG(abs(m_x) <= 1, "Three term recurrence is stable only for |x| <=1.");
  34. // The integer n indexes a family of discrete Legendre polynomials indexed by k <= 2*n
  35. }
  36. Real norm_sq(int r) const
  37. {
  38. Real prod = Real(2) / Real(2 * r + 1);
  39. for (int k = -r; k <= r; ++k) {
  40. prod *= Real(2 * m_n + 1 + k) / Real(2 * m_n);
  41. }
  42. return prod;
  43. }
  44. Real next()
  45. {
  46. Real N = 2 * m_n + 1;
  47. Real num = (m_r - 1) * (N * N - (m_r - 1) * (m_r - 1)) * m_qrm2;
  48. Real tmp = (2 * m_r - 1) * m_x * m_qrm1 - num / Real(4 * m_n * m_n);
  49. m_qrm2 = m_qrm1;
  50. m_qrm1 = tmp / m_r;
  51. ++m_r;
  52. return m_qrm1;
  53. }
  54. Real next_prime()
  55. {
  56. Real N = 2 * m_n + 1;
  57. Real s = (m_r - 1) * (N * N - (m_r - 1) * (m_r - 1)) / Real(4 * m_n * m_n);
  58. Real tmp1 = ((2 * m_r - 1) * m_x * m_qrm1 - s * m_qrm2) / m_r;
  59. Real tmp2 = ((2 * m_r - 1) * (m_qrm1 + m_x * m_qrm1p) - s * m_qrm2p) / m_r;
  60. m_qrm2 = m_qrm1;
  61. m_qrm1 = tmp1;
  62. m_qrm2p = m_qrm1p;
  63. m_qrm1p = tmp2;
  64. ++m_r;
  65. return m_qrm1p;
  66. }
  67. Real next_dbl_prime()
  68. {
  69. Real N = 2*m_n + 1;
  70. Real trm1 = 2*m_r - 1;
  71. Real s = (m_r - 1) * (N * N - (m_r - 1) * (m_r - 1)) / Real(4 * m_n * m_n);
  72. Real rqrpp = 2*trm1*m_qrm1p + trm1*m_x*m_qrm1pp - s*m_qrm2pp;
  73. Real tmp1 = ((2 * m_r - 1) * m_x * m_qrm1 - s * m_qrm2) / m_r;
  74. Real tmp2 = ((2 * m_r - 1) * (m_qrm1 + m_x * m_qrm1p) - s * m_qrm2p) / m_r;
  75. m_qrm2 = m_qrm1;
  76. m_qrm1 = tmp1;
  77. m_qrm2p = m_qrm1p;
  78. m_qrm1p = tmp2;
  79. m_qrm2pp = m_qrm1pp;
  80. m_qrm1pp = rqrpp/m_r;
  81. ++m_r;
  82. return m_qrm1pp;
  83. }
  84. Real operator()(Real x, std::size_t k)
  85. {
  86. BOOST_MATH_ASSERT_MSG(k <= 2 * m_n, "r <= 2n is required.");
  87. if (k == 0)
  88. {
  89. return 1;
  90. }
  91. if (k == 1)
  92. {
  93. return x;
  94. }
  95. Real qrm2 = 1;
  96. Real qrm1 = x;
  97. Real N = 2 * m_n + 1;
  98. for (std::size_t r = 2; r <= k; ++r) {
  99. Real num = (r - 1) * (N * N - (r - 1) * (r - 1)) * qrm2;
  100. Real tmp = (2 * r - 1) * x * qrm1 - num / Real(4 * m_n * m_n);
  101. qrm2 = qrm1;
  102. qrm1 = tmp / r;
  103. }
  104. return qrm1;
  105. }
  106. Real prime(Real x, std::size_t k) {
  107. BOOST_MATH_ASSERT_MSG(k <= 2 * m_n, "r <= 2n is required.");
  108. if (k == 0) {
  109. return 0;
  110. }
  111. if (k == 1) {
  112. return 1;
  113. }
  114. Real qrm2 = 1;
  115. Real qrm1 = x;
  116. Real qrm2p = 0;
  117. Real qrm1p = 1;
  118. Real N = 2 * m_n + 1;
  119. for (std::size_t r = 2; r <= k; ++r) {
  120. Real s =
  121. (r - 1) * (N * N - (r - 1) * (r - 1)) / Real(4 * m_n * m_n);
  122. Real tmp1 = ((2 * r - 1) * x * qrm1 - s * qrm2) / r;
  123. Real tmp2 = ((2 * r - 1) * (qrm1 + x * qrm1p) - s * qrm2p) / r;
  124. qrm2 = qrm1;
  125. qrm1 = tmp1;
  126. qrm2p = qrm1p;
  127. qrm1p = tmp2;
  128. }
  129. return qrm1p;
  130. }
  131. private:
  132. std::size_t m_n;
  133. std::size_t m_r;
  134. Real m_x;
  135. Real m_qrm2;
  136. Real m_qrm1;
  137. Real m_qrm2p;
  138. Real m_qrm1p;
  139. Real m_qrm2pp;
  140. Real m_qrm1pp;
  141. };
  142. template <class Real>
  143. std::vector<Real> interior_velocity_filter(std::size_t n, std::size_t p) {
  144. auto dlp = discrete_legendre<Real>(n, 0);
  145. std::vector<Real> coeffs(p+1);
  146. coeffs[1] = 1/dlp.norm_sq(1);
  147. for (std::size_t l = 3; l < p + 1; l += 2)
  148. {
  149. dlp.next_prime();
  150. coeffs[l] = dlp.next_prime()/ dlp.norm_sq(l);
  151. }
  152. // We could make the filter length n, as f[0] = 0,
  153. // but that'd make the indexing awkward when applying the filter.
  154. std::vector<Real> f(n + 1);
  155. // This value should never be read, but this is the correct value *if it is read*.
  156. // Hmm, should it be a nan then? I'm not gonna agonize.
  157. f[0] = 0;
  158. for (std::size_t j = 1; j < f.size(); ++j)
  159. {
  160. Real arg = Real(j) / Real(n);
  161. dlp = discrete_legendre<Real>(n, arg);
  162. f[j] = coeffs[1]*arg;
  163. for (std::size_t l = 3; l <= p; l += 2)
  164. {
  165. dlp.next();
  166. f[j] += coeffs[l]*dlp.next();
  167. }
  168. f[j] /= (n * n);
  169. }
  170. return f;
  171. }
  172. template <class Real>
  173. std::vector<Real> boundary_velocity_filter(std::size_t n, std::size_t p, int64_t s)
  174. {
  175. std::vector<Real> coeffs(p+1, std::numeric_limits<Real>::quiet_NaN());
  176. Real sn = Real(s) / Real(n);
  177. auto dlp = discrete_legendre<Real>(n, sn);
  178. coeffs[0] = 0;
  179. coeffs[1] = 1/dlp.norm_sq(1);
  180. for (std::size_t l = 2; l < p + 1; ++l)
  181. {
  182. // Calculation of the norms is common to all filters,
  183. // so it seems like an obvious optimization target.
  184. // I tried this: The spent in computing the norms time is not negligible,
  185. // but still a small fraction of the total compute time.
  186. // Hence I'm not refactoring out these norm calculations.
  187. coeffs[l] = dlp.next_prime()/ dlp.norm_sq(l);
  188. }
  189. std::vector<Real> f(2*n + 1);
  190. for (std::size_t k = 0; k < f.size(); ++k)
  191. {
  192. Real j = Real(k) - Real(n);
  193. Real arg = j/Real(n);
  194. dlp = discrete_legendre<Real>(n, arg);
  195. f[k] = coeffs[1]*arg;
  196. for (std::size_t l = 2; l <= p; ++l)
  197. {
  198. f[k] += coeffs[l]*dlp.next();
  199. }
  200. f[k] /= (n * n);
  201. }
  202. return f;
  203. }
  204. template <class Real>
  205. std::vector<Real> acceleration_filter(std::size_t n, std::size_t p, int64_t s)
  206. {
  207. BOOST_MATH_ASSERT_MSG(p <= 2*n, "Approximation order must be <= 2*n");
  208. BOOST_MATH_ASSERT_MSG(p > 2, "Approximation order must be > 2");
  209. std::vector<Real> coeffs(p+1, std::numeric_limits<Real>::quiet_NaN());
  210. Real sn = Real(s) / Real(n);
  211. auto dlp = discrete_legendre<Real>(n, sn);
  212. coeffs[0] = 0;
  213. coeffs[1] = 0;
  214. for (std::size_t l = 2; l < p + 1; ++l)
  215. {
  216. coeffs[l] = dlp.next_dbl_prime()/ dlp.norm_sq(l);
  217. }
  218. std::vector<Real> f(2*n + 1, 0);
  219. for (std::size_t k = 0; k < f.size(); ++k)
  220. {
  221. Real j = Real(k) - Real(n);
  222. Real arg = j/Real(n);
  223. dlp = discrete_legendre<Real>(n, arg);
  224. for (std::size_t l = 2; l <= p; ++l)
  225. {
  226. f[k] += coeffs[l]*dlp.next();
  227. }
  228. f[k] /= (n * n * n);
  229. }
  230. return f;
  231. }
  232. } // namespace detail
  233. template <typename Real, std::size_t order = 1>
  234. class discrete_lanczos_derivative {
  235. public:
  236. discrete_lanczos_derivative(Real const & spacing,
  237. std::size_t n = 18,
  238. std::size_t approximation_order = 3)
  239. : m_dt{spacing}
  240. {
  241. static_assert(!std::is_integral_v<Real>,
  242. "Spacing must be a floating point type.");
  243. BOOST_MATH_ASSERT_MSG(spacing > 0,
  244. "Spacing between samples must be > 0.");
  245. if constexpr (order == 1)
  246. {
  247. BOOST_MATH_ASSERT_MSG(approximation_order <= 2 * n,
  248. "The approximation order must be <= 2n");
  249. BOOST_MATH_ASSERT_MSG(approximation_order >= 2,
  250. "The approximation order must be >= 2");
  251. if constexpr (std::is_same_v<Real, float> || std::is_same_v<Real, double>)
  252. {
  253. auto interior = detail::interior_velocity_filter<long double>(n, approximation_order);
  254. m_f.resize(interior.size());
  255. for (std::size_t j = 0; j < interior.size(); ++j)
  256. {
  257. m_f[j] = static_cast<Real>(interior[j])/m_dt;
  258. }
  259. }
  260. else
  261. {
  262. m_f = detail::interior_velocity_filter<Real>(n, approximation_order);
  263. for (auto & x : m_f)
  264. {
  265. x /= m_dt;
  266. }
  267. }
  268. m_boundary_filters.resize(n);
  269. // This for loop is a natural candidate for parallelization.
  270. // But does it matter? Probably not.
  271. for (std::size_t i = 0; i < n; ++i)
  272. {
  273. if constexpr (std::is_same_v<Real, float> || std::is_same_v<Real, double>)
  274. {
  275. int64_t s = static_cast<int64_t>(i) - static_cast<int64_t>(n);
  276. auto bf = detail::boundary_velocity_filter<long double>(n, approximation_order, s);
  277. m_boundary_filters[i].resize(bf.size());
  278. for (std::size_t j = 0; j < bf.size(); ++j)
  279. {
  280. m_boundary_filters[i][j] = static_cast<Real>(bf[j])/m_dt;
  281. }
  282. }
  283. else
  284. {
  285. int64_t s = static_cast<int64_t>(i) - static_cast<int64_t>(n);
  286. m_boundary_filters[i] = detail::boundary_velocity_filter<Real>(n, approximation_order, s);
  287. for (auto & bf : m_boundary_filters[i])
  288. {
  289. bf /= m_dt;
  290. }
  291. }
  292. }
  293. }
  294. else if constexpr (order == 2)
  295. {
  296. // High precision isn't warranted for small p; only for large p.
  297. // (The computation appears stable for large n.)
  298. // But given that the filters are reusable for many vectors,
  299. // it's better to do a high precision computation and then cast back,
  300. // since the resulting cost is a factor of 2, and the cost of the filters not working is hours of debugging.
  301. if constexpr (std::is_same_v<Real, double> || std::is_same_v<Real, float>)
  302. {
  303. auto f = detail::acceleration_filter<long double>(n, approximation_order, 0);
  304. m_f.resize(n+1);
  305. for (std::size_t i = 0; i < m_f.size(); ++i)
  306. {
  307. m_f[i] = static_cast<Real>(f[i+n])/(m_dt*m_dt);
  308. }
  309. m_boundary_filters.resize(n);
  310. for (std::size_t i = 0; i < n; ++i)
  311. {
  312. int64_t s = static_cast<int64_t>(i) - static_cast<int64_t>(n);
  313. auto bf = detail::acceleration_filter<long double>(n, approximation_order, s);
  314. m_boundary_filters[i].resize(bf.size());
  315. for (std::size_t j = 0; j < bf.size(); ++j)
  316. {
  317. m_boundary_filters[i][j] = static_cast<Real>(bf[j])/(m_dt*m_dt);
  318. }
  319. }
  320. }
  321. else
  322. {
  323. // Given that the purpose is denoising, for higher precision calculations,
  324. // the default precision should be fine.
  325. auto f = detail::acceleration_filter<Real>(n, approximation_order, 0);
  326. m_f.resize(n+1);
  327. for (std::size_t i = 0; i < m_f.size(); ++i)
  328. {
  329. m_f[i] = f[i+n]/(m_dt*m_dt);
  330. }
  331. m_boundary_filters.resize(n);
  332. for (std::size_t i = 0; i < n; ++i)
  333. {
  334. int64_t s = static_cast<int64_t>(i) - static_cast<int64_t>(n);
  335. m_boundary_filters[i] = detail::acceleration_filter<Real>(n, approximation_order, s);
  336. for (auto & bf : m_boundary_filters[i])
  337. {
  338. bf /= (m_dt*m_dt);
  339. }
  340. }
  341. }
  342. }
  343. else
  344. {
  345. BOOST_MATH_ASSERT_MSG(false, "Derivatives of order 3 and higher are not implemented.");
  346. }
  347. }
  348. Real get_spacing() const
  349. {
  350. return m_dt;
  351. }
  352. template<class RandomAccessContainer>
  353. Real operator()(RandomAccessContainer const & v, std::size_t i) const
  354. {
  355. static_assert(std::is_same_v<typename RandomAccessContainer::value_type, Real>,
  356. "The type of the values in the vector provided does not match the type in the filters.");
  357. BOOST_MATH_ASSERT_MSG(std::size(v) >= m_boundary_filters[0].size(),
  358. "Vector must be at least as long as the filter length");
  359. if constexpr (order==1)
  360. {
  361. if (i >= m_f.size() - 1 && i <= std::size(v) - m_f.size())
  362. {
  363. // The filter has length >= 1:
  364. Real dvdt = m_f[1] * (v[i + 1] - v[i - 1]);
  365. for (std::size_t j = 2; j < m_f.size(); ++j)
  366. {
  367. dvdt += m_f[j] * (v[i + j] - v[i - j]);
  368. }
  369. return dvdt;
  370. }
  371. // m_f.size() = N+1
  372. if (i < m_f.size() - 1)
  373. {
  374. auto &bf = m_boundary_filters[i];
  375. Real dvdt = bf[0]*v[0];
  376. for (std::size_t j = 1; j < bf.size(); ++j)
  377. {
  378. dvdt += bf[j] * v[j];
  379. }
  380. return dvdt;
  381. }
  382. if (i > std::size(v) - m_f.size() && i < std::size(v))
  383. {
  384. int k = std::size(v) - 1 - i;
  385. auto &bf = m_boundary_filters[k];
  386. Real dvdt = bf[0]*v[std::size(v)-1];
  387. for (std::size_t j = 1; j < bf.size(); ++j)
  388. {
  389. dvdt += bf[j] * v[std::size(v) - 1 - j];
  390. }
  391. return -dvdt;
  392. }
  393. }
  394. else if constexpr (order==2)
  395. {
  396. if (i >= m_f.size() - 1 && i <= std::size(v) - m_f.size())
  397. {
  398. Real d2vdt2 = m_f[0]*v[i];
  399. for (std::size_t j = 1; j < m_f.size(); ++j)
  400. {
  401. d2vdt2 += m_f[j] * (v[i + j] + v[i - j]);
  402. }
  403. return d2vdt2;
  404. }
  405. // m_f.size() = N+1
  406. if (i < m_f.size() - 1)
  407. {
  408. auto &bf = m_boundary_filters[i];
  409. Real d2vdt2 = bf[0]*v[0];
  410. for (std::size_t j = 1; j < bf.size(); ++j)
  411. {
  412. d2vdt2 += bf[j] * v[j];
  413. }
  414. return d2vdt2;
  415. }
  416. if (i > std::size(v) - m_f.size() && i < std::size(v))
  417. {
  418. int k = std::size(v) - 1 - i;
  419. auto &bf = m_boundary_filters[k];
  420. Real d2vdt2 = bf[0] * v[std::size(v) - 1];
  421. for (std::size_t j = 1; j < bf.size(); ++j)
  422. {
  423. d2vdt2 += bf[j] * v[std::size(v) - 1 - j];
  424. }
  425. return d2vdt2;
  426. }
  427. }
  428. // OOB access:
  429. std::string msg = "Out of bounds access in Lanczos derivative.";
  430. msg += "Input vector has length " + std::to_string(std::size(v)) + ", but user requested access at index " + std::to_string(i) + ".";
  431. throw std::out_of_range(msg);
  432. return std::numeric_limits<Real>::quiet_NaN();
  433. }
  434. template<class RandomAccessContainer>
  435. void operator()(RandomAccessContainer const & v, RandomAccessContainer & w) const
  436. {
  437. static_assert(std::is_same_v<typename RandomAccessContainer::value_type, Real>,
  438. "The type of the values in the vector provided does not match the type in the filters.");
  439. if (&w[0] == &v[0])
  440. {
  441. throw std::logic_error("This transform cannot be performed in-place.");
  442. }
  443. if (std::size(v) < m_boundary_filters[0].size())
  444. {
  445. std::string msg = "The input vector must be at least as long as the filter length. ";
  446. msg += "The input vector has length = " + std::to_string(std::size(v)) + ", the filter has length " + std::to_string(m_boundary_filters[0].size());
  447. throw std::length_error(msg);
  448. }
  449. if (std::size(w) < std::size(v))
  450. {
  451. std::string msg = "The output vector (containing the derivative) must be at least as long as the input vector.";
  452. msg += "The output vector has length = " + std::to_string(std::size(w)) + ", the input vector has length " + std::to_string(std::size(v));
  453. throw std::length_error(msg);
  454. }
  455. if constexpr (order==1)
  456. {
  457. for (std::size_t i = 0; i < m_f.size() - 1; ++i)
  458. {
  459. auto &bf = m_boundary_filters[i];
  460. Real dvdt = bf[0] * v[0];
  461. for (std::size_t j = 1; j < bf.size(); ++j)
  462. {
  463. dvdt += bf[j] * v[j];
  464. }
  465. w[i] = dvdt;
  466. }
  467. for(std::size_t i = m_f.size() - 1; i <= std::size(v) - m_f.size(); ++i)
  468. {
  469. Real dvdt = m_f[1] * (v[i + 1] - v[i - 1]);
  470. for (std::size_t j = 2; j < m_f.size(); ++j)
  471. {
  472. dvdt += m_f[j] *(v[i + j] - v[i - j]);
  473. }
  474. w[i] = dvdt;
  475. }
  476. for(std::size_t i = std::size(v) - m_f.size() + 1; i < std::size(v); ++i)
  477. {
  478. int k = std::size(v) - 1 - i;
  479. auto &f = m_boundary_filters[k];
  480. Real dvdt = f[0] * v[std::size(v) - 1];;
  481. for (std::size_t j = 1; j < f.size(); ++j)
  482. {
  483. dvdt += f[j] * v[std::size(v) - 1 - j];
  484. }
  485. w[i] = -dvdt;
  486. }
  487. }
  488. else if constexpr (order==2)
  489. {
  490. // m_f.size() = N+1
  491. for (std::size_t i = 0; i < m_f.size() - 1; ++i)
  492. {
  493. auto &bf = m_boundary_filters[i];
  494. Real d2vdt2 = 0;
  495. for (std::size_t j = 0; j < bf.size(); ++j)
  496. {
  497. d2vdt2 += bf[j] * v[j];
  498. }
  499. w[i] = d2vdt2;
  500. }
  501. for (std::size_t i = m_f.size() - 1; i <= std::size(v) - m_f.size(); ++i)
  502. {
  503. Real d2vdt2 = m_f[0]*v[i];
  504. for (std::size_t j = 1; j < m_f.size(); ++j)
  505. {
  506. d2vdt2 += m_f[j] * (v[i + j] + v[i - j]);
  507. }
  508. w[i] = d2vdt2;
  509. }
  510. for (std::size_t i = std::size(v) - m_f.size() + 1; i < std::size(v); ++i)
  511. {
  512. int k = std::size(v) - 1 - i;
  513. auto &bf = m_boundary_filters[k];
  514. Real d2vdt2 = bf[0] * v[std::size(v) - 1];
  515. for (std::size_t j = 1; j < bf.size(); ++j)
  516. {
  517. d2vdt2 += bf[j] * v[std::size(v) - 1 - j];
  518. }
  519. w[i] = d2vdt2;
  520. }
  521. }
  522. }
  523. template<class RandomAccessContainer>
  524. RandomAccessContainer operator()(RandomAccessContainer const & v) const
  525. {
  526. RandomAccessContainer w(std::size(v));
  527. this->operator()(v, w);
  528. return w;
  529. }
  530. // Don't copy; too big.
  531. discrete_lanczos_derivative( const discrete_lanczos_derivative & ) = delete;
  532. discrete_lanczos_derivative& operator=(const discrete_lanczos_derivative&) = delete;
  533. // Allow moves:
  534. discrete_lanczos_derivative(discrete_lanczos_derivative&&) noexcept = default;
  535. discrete_lanczos_derivative& operator=(discrete_lanczos_derivative&&) noexcept = default;
  536. private:
  537. std::vector<Real> m_f;
  538. std::vector<std::vector<Real>> m_boundary_filters;
  539. Real m_dt;
  540. };
  541. } // namespaces
  542. #endif