cardinal_trigonometric_detail.hpp 17 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676
  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_INTERPOLATORS_DETAIL_CARDINAL_TRIGONOMETRIC_HPP
  6. #define BOOST_MATH_INTERPOLATORS_DETAIL_CARDINAL_TRIGONOMETRIC_HPP
  7. #include <cmath>
  8. #include <stdexcept>
  9. #include <fftw3.h>
  10. #include <boost/math/constants/constants.hpp>
  11. #ifdef BOOST_HAS_FLOAT128
  12. #include <quadmath.h>
  13. #endif
  14. namespace boost { namespace math { namespace interpolators { namespace detail {
  15. template<typename Real>
  16. class cardinal_trigonometric_detail {
  17. public:
  18. cardinal_trigonometric_detail(const Real* data, size_t length, Real t0, Real h)
  19. {
  20. m_data = data;
  21. m_length = length;
  22. m_t0 = t0;
  23. m_h = h;
  24. throw std::domain_error("Not implemented.");
  25. }
  26. private:
  27. size_t m_length;
  28. Real m_t0;
  29. Real m_h;
  30. Real* m_data;
  31. };
  32. template<>
  33. class cardinal_trigonometric_detail<float> {
  34. public:
  35. cardinal_trigonometric_detail(const float* data, size_t length, float t0, float h) : m_t0{t0}, m_h{h}
  36. {
  37. if (length == 0)
  38. {
  39. throw std::logic_error("At least one sample is required.");
  40. }
  41. if (h <= 0)
  42. {
  43. throw std::logic_error("The step size must be > 0");
  44. }
  45. // The period sadly must be stored, since the complex vector has length that cannot be used to recover the period:
  46. m_T = m_h*length;
  47. m_complex_vector_size = length/2 + 1;
  48. m_gamma = fftwf_alloc_complex(m_complex_vector_size);
  49. // The const_cast is legitimate: FFTW does not change the data as long as FFTW_ESTIMATE is provided.
  50. fftwf_plan plan = fftwf_plan_dft_r2c_1d(length, const_cast<float*>(data), m_gamma, FFTW_ESTIMATE);
  51. // FFTW says a null plan is impossible with the basic interface we are using, and I have no reason to doubt them.
  52. // But it just feels weird not to check this:
  53. if (!plan)
  54. {
  55. throw std::logic_error("A null fftw plan was created.");
  56. }
  57. fftwf_execute(plan);
  58. fftwf_destroy_plan(plan);
  59. float denom = length;
  60. for (size_t k = 0; k < m_complex_vector_size; ++k)
  61. {
  62. m_gamma[k][0] /= denom;
  63. m_gamma[k][1] /= denom;
  64. }
  65. if (length % 2 == 0)
  66. {
  67. m_gamma[m_complex_vector_size -1][0] /= 2;
  68. // numerically, m_gamma[m_complex_vector_size -1][1] should be zero . . .
  69. // I believe, but need to check, that FFTW guarantees that it is identically zero.
  70. }
  71. }
  72. cardinal_trigonometric_detail(const cardinal_trigonometric_detail& old) = delete;
  73. cardinal_trigonometric_detail& operator=(const cardinal_trigonometric_detail&) = delete;
  74. cardinal_trigonometric_detail(cardinal_trigonometric_detail &&) = delete;
  75. float operator()(float t) const
  76. {
  77. using std::sin;
  78. using std::cos;
  79. using boost::math::constants::two_pi;
  80. using std::exp;
  81. float s = m_gamma[0][0];
  82. float x = two_pi<float>()*(t - m_t0)/m_T;
  83. fftwf_complex z;
  84. // boost::math::cos_pi with a redefinition of x? Not now . . .
  85. z[0] = cos(x);
  86. z[1] = sin(x);
  87. fftwf_complex b{0, 0};
  88. // u = b*z
  89. fftwf_complex u;
  90. for (size_t k = m_complex_vector_size - 1; k >= 1; --k) {
  91. u[0] = b[0]*z[0] - b[1]*z[1];
  92. u[1] = b[0]*z[1] + b[1]*z[0];
  93. b[0] = m_gamma[k][0] + u[0];
  94. b[1] = m_gamma[k][1] + u[1];
  95. }
  96. s += 2*(b[0]*z[0] - b[1]*z[1]);
  97. return s;
  98. }
  99. float prime(float t) const
  100. {
  101. using std::sin;
  102. using std::cos;
  103. using boost::math::constants::two_pi;
  104. using std::exp;
  105. float x = two_pi<float>()*(t - m_t0)/m_T;
  106. fftwf_complex z;
  107. z[0] = cos(x);
  108. z[1] = sin(x);
  109. fftwf_complex b{0, 0};
  110. // u = b*z
  111. fftwf_complex u;
  112. for (size_t k = m_complex_vector_size - 1; k >= 1; --k)
  113. {
  114. u[0] = b[0]*z[0] - b[1]*z[1];
  115. u[1] = b[0]*z[1] + b[1]*z[0];
  116. b[0] = k*m_gamma[k][0] + u[0];
  117. b[1] = k*m_gamma[k][1] + u[1];
  118. }
  119. // b*z = (b[0]*z[0] - b[1]*z[1]) + i(b[1]*z[0] + b[0]*z[1])
  120. return -2*two_pi<float>()*(b[1]*z[0] + b[0]*z[1])/m_T;
  121. }
  122. float double_prime(float t) const
  123. {
  124. using std::sin;
  125. using std::cos;
  126. using boost::math::constants::two_pi;
  127. using std::exp;
  128. float x = two_pi<float>()*(t - m_t0)/m_T;
  129. fftwf_complex z;
  130. z[0] = cos(x);
  131. z[1] = sin(x);
  132. fftwf_complex b{0, 0};
  133. // u = b*z
  134. fftwf_complex u;
  135. for (size_t k = m_complex_vector_size - 1; k >= 1; --k)
  136. {
  137. u[0] = b[0]*z[0] - b[1]*z[1];
  138. u[1] = b[0]*z[1] + b[1]*z[0];
  139. b[0] = k*k*m_gamma[k][0] + u[0];
  140. b[1] = k*k*m_gamma[k][1] + u[1];
  141. }
  142. // b*z = (b[0]*z[0] - b[1]*z[1]) + i(b[1]*z[0] + b[0]*z[1])
  143. return -2*two_pi<float>()*two_pi<float>()*(b[0]*z[0] - b[1]*z[1])/(m_T*m_T);
  144. }
  145. float period() const
  146. {
  147. return m_T;
  148. }
  149. float integrate() const
  150. {
  151. return m_T*m_gamma[0][0];
  152. }
  153. float squared_l2() const
  154. {
  155. float s = 0;
  156. // Always add smallest to largest for accuracy.
  157. for (size_t i = m_complex_vector_size - 1; i >= 1; --i)
  158. {
  159. s += (m_gamma[i][0]*m_gamma[i][0] + m_gamma[i][1]*m_gamma[i][1]);
  160. }
  161. s *= 2;
  162. s += m_gamma[0][0]*m_gamma[0][0];
  163. return s*m_T;
  164. }
  165. ~cardinal_trigonometric_detail()
  166. {
  167. if (m_gamma)
  168. {
  169. fftwf_free(m_gamma);
  170. m_gamma = nullptr;
  171. }
  172. }
  173. private:
  174. float m_t0;
  175. float m_h;
  176. float m_T;
  177. fftwf_complex* m_gamma;
  178. size_t m_complex_vector_size;
  179. };
  180. template<>
  181. class cardinal_trigonometric_detail<double> {
  182. public:
  183. cardinal_trigonometric_detail(const double* data, size_t length, double t0, double h) : m_t0{t0}, m_h{h}
  184. {
  185. if (length == 0)
  186. {
  187. throw std::logic_error("At least one sample is required.");
  188. }
  189. if (h <= 0)
  190. {
  191. throw std::logic_error("The step size must be > 0");
  192. }
  193. m_T = m_h*length;
  194. m_complex_vector_size = length/2 + 1;
  195. m_gamma = fftw_alloc_complex(m_complex_vector_size);
  196. fftw_plan plan = fftw_plan_dft_r2c_1d(length, const_cast<double*>(data), m_gamma, FFTW_ESTIMATE);
  197. if (!plan)
  198. {
  199. throw std::logic_error("A null fftw plan was created.");
  200. }
  201. fftw_execute(plan);
  202. fftw_destroy_plan(plan);
  203. double denom = length;
  204. for (size_t k = 0; k < m_complex_vector_size; ++k)
  205. {
  206. m_gamma[k][0] /= denom;
  207. m_gamma[k][1] /= denom;
  208. }
  209. if (length % 2 == 0)
  210. {
  211. m_gamma[m_complex_vector_size -1][0] /= 2;
  212. }
  213. }
  214. cardinal_trigonometric_detail(const cardinal_trigonometric_detail& old) = delete;
  215. cardinal_trigonometric_detail& operator=(const cardinal_trigonometric_detail&) = delete;
  216. cardinal_trigonometric_detail(cardinal_trigonometric_detail &&) = delete;
  217. double operator()(double t) const
  218. {
  219. using std::sin;
  220. using std::cos;
  221. using boost::math::constants::two_pi;
  222. using std::exp;
  223. double s = m_gamma[0][0];
  224. double x = two_pi<double>()*(t - m_t0)/m_T;
  225. fftw_complex z;
  226. z[0] = cos(x);
  227. z[1] = sin(x);
  228. fftw_complex b{0, 0};
  229. // u = b*z
  230. fftw_complex u;
  231. for (size_t k = m_complex_vector_size - 1; k >= 1; --k)
  232. {
  233. u[0] = b[0]*z[0] - b[1]*z[1];
  234. u[1] = b[0]*z[1] + b[1]*z[0];
  235. b[0] = m_gamma[k][0] + u[0];
  236. b[1] = m_gamma[k][1] + u[1];
  237. }
  238. s += 2*(b[0]*z[0] - b[1]*z[1]);
  239. return s;
  240. }
  241. double prime(double t) const
  242. {
  243. using std::sin;
  244. using std::cos;
  245. using boost::math::constants::two_pi;
  246. using std::exp;
  247. double x = two_pi<double>()*(t - m_t0)/m_T;
  248. fftw_complex z;
  249. z[0] = cos(x);
  250. z[1] = sin(x);
  251. fftw_complex b{0, 0};
  252. // u = b*z
  253. fftw_complex u;
  254. for (size_t k = m_complex_vector_size - 1; k >= 1; --k)
  255. {
  256. u[0] = b[0]*z[0] - b[1]*z[1];
  257. u[1] = b[0]*z[1] + b[1]*z[0];
  258. b[0] = k*m_gamma[k][0] + u[0];
  259. b[1] = k*m_gamma[k][1] + u[1];
  260. }
  261. // b*z = (b[0]*z[0] - b[1]*z[1]) + i(b[1]*z[0] + b[0]*z[1])
  262. return -2*two_pi<double>()*(b[1]*z[0] + b[0]*z[1])/m_T;
  263. }
  264. double double_prime(double t) const
  265. {
  266. using std::sin;
  267. using std::cos;
  268. using boost::math::constants::two_pi;
  269. using std::exp;
  270. double x = two_pi<double>()*(t - m_t0)/m_T;
  271. fftw_complex z;
  272. z[0] = cos(x);
  273. z[1] = sin(x);
  274. fftw_complex b{0, 0};
  275. // u = b*z
  276. fftw_complex u;
  277. for (size_t k = m_complex_vector_size - 1; k >= 1; --k)
  278. {
  279. u[0] = b[0]*z[0] - b[1]*z[1];
  280. u[1] = b[0]*z[1] + b[1]*z[0];
  281. b[0] = k*k*m_gamma[k][0] + u[0];
  282. b[1] = k*k*m_gamma[k][1] + u[1];
  283. }
  284. // b*z = (b[0]*z[0] - b[1]*z[1]) + i(b[1]*z[0] + b[0]*z[1])
  285. return -2*two_pi<double>()*two_pi<double>()*(b[0]*z[0] - b[1]*z[1])/(m_T*m_T);
  286. }
  287. double period() const
  288. {
  289. return m_T;
  290. }
  291. double integrate() const
  292. {
  293. return m_T*m_gamma[0][0];
  294. }
  295. double squared_l2() const
  296. {
  297. double s = 0;
  298. for (size_t i = m_complex_vector_size - 1; i >= 1; --i)
  299. {
  300. s += (m_gamma[i][0]*m_gamma[i][0] + m_gamma[i][1]*m_gamma[i][1]);
  301. }
  302. s *= 2;
  303. s += m_gamma[0][0]*m_gamma[0][0];
  304. return s*m_T;
  305. }
  306. ~cardinal_trigonometric_detail()
  307. {
  308. if (m_gamma)
  309. {
  310. fftw_free(m_gamma);
  311. m_gamma = nullptr;
  312. }
  313. }
  314. private:
  315. double m_t0;
  316. double m_h;
  317. double m_T;
  318. fftw_complex* m_gamma;
  319. size_t m_complex_vector_size;
  320. };
  321. template<>
  322. class cardinal_trigonometric_detail<long double> {
  323. public:
  324. cardinal_trigonometric_detail(const long double* data, size_t length, long double t0, long double h) : m_t0{t0}, m_h{h}
  325. {
  326. if (length == 0)
  327. {
  328. throw std::logic_error("At least one sample is required.");
  329. }
  330. if (h <= 0)
  331. {
  332. throw std::logic_error("The step size must be > 0");
  333. }
  334. m_T = m_h*length;
  335. m_complex_vector_size = length/2 + 1;
  336. m_gamma = fftwl_alloc_complex(m_complex_vector_size);
  337. fftwl_plan plan = fftwl_plan_dft_r2c_1d(length, const_cast<long double*>(data), m_gamma, FFTW_ESTIMATE);
  338. if (!plan)
  339. {
  340. throw std::logic_error("A null fftw plan was created.");
  341. }
  342. fftwl_execute(plan);
  343. fftwl_destroy_plan(plan);
  344. long double denom = length;
  345. for (size_t k = 0; k < m_complex_vector_size; ++k)
  346. {
  347. m_gamma[k][0] /= denom;
  348. m_gamma[k][1] /= denom;
  349. }
  350. if (length % 2 == 0) {
  351. m_gamma[m_complex_vector_size -1][0] /= 2;
  352. }
  353. }
  354. cardinal_trigonometric_detail(const cardinal_trigonometric_detail& old) = delete;
  355. cardinal_trigonometric_detail& operator=(const cardinal_trigonometric_detail&) = delete;
  356. cardinal_trigonometric_detail(cardinal_trigonometric_detail &&) = delete;
  357. long double operator()(long double t) const
  358. {
  359. using std::sin;
  360. using std::cos;
  361. using boost::math::constants::two_pi;
  362. using std::exp;
  363. long double s = m_gamma[0][0];
  364. long double x = two_pi<long double>()*(t - m_t0)/m_T;
  365. fftwl_complex z;
  366. z[0] = cos(x);
  367. z[1] = sin(x);
  368. fftwl_complex b{0, 0};
  369. fftwl_complex u;
  370. for (size_t k = m_complex_vector_size - 1; k >= 1; --k)
  371. {
  372. u[0] = b[0]*z[0] - b[1]*z[1];
  373. u[1] = b[0]*z[1] + b[1]*z[0];
  374. b[0] = m_gamma[k][0] + u[0];
  375. b[1] = m_gamma[k][1] + u[1];
  376. }
  377. s += 2*(b[0]*z[0] - b[1]*z[1]);
  378. return s;
  379. }
  380. long double prime(long double t) const
  381. {
  382. using std::sin;
  383. using std::cos;
  384. using boost::math::constants::two_pi;
  385. using std::exp;
  386. long double x = two_pi<long double>()*(t - m_t0)/m_T;
  387. fftwl_complex z;
  388. z[0] = cos(x);
  389. z[1] = sin(x);
  390. fftwl_complex b{0, 0};
  391. // u = b*z
  392. fftwl_complex u;
  393. for (size_t k = m_complex_vector_size - 1; k >= 1; --k)
  394. {
  395. u[0] = b[0]*z[0] - b[1]*z[1];
  396. u[1] = b[0]*z[1] + b[1]*z[0];
  397. b[0] = k*m_gamma[k][0] + u[0];
  398. b[1] = k*m_gamma[k][1] + u[1];
  399. }
  400. // b*z = (b[0]*z[0] - b[1]*z[1]) + i(b[1]*z[0] + b[0]*z[1])
  401. return -2*two_pi<long double>()*(b[1]*z[0] + b[0]*z[1])/m_T;
  402. }
  403. long double double_prime(long double t) const
  404. {
  405. using std::sin;
  406. using std::cos;
  407. using boost::math::constants::two_pi;
  408. using std::exp;
  409. long double x = two_pi<long double>()*(t - m_t0)/m_T;
  410. fftwl_complex z;
  411. z[0] = cos(x);
  412. z[1] = sin(x);
  413. fftwl_complex b{0, 0};
  414. // u = b*z
  415. fftwl_complex u;
  416. for (size_t k = m_complex_vector_size - 1; k >= 1; --k)
  417. {
  418. u[0] = b[0]*z[0] - b[1]*z[1];
  419. u[1] = b[0]*z[1] + b[1]*z[0];
  420. b[0] = k*k*m_gamma[k][0] + u[0];
  421. b[1] = k*k*m_gamma[k][1] + u[1];
  422. }
  423. // b*z = (b[0]*z[0] - b[1]*z[1]) + i(b[1]*z[0] + b[0]*z[1])
  424. return -2*two_pi<long double>()*two_pi<long double>()*(b[0]*z[0] - b[1]*z[1])/(m_T*m_T);
  425. }
  426. long double period() const
  427. {
  428. return m_T;
  429. }
  430. long double integrate() const
  431. {
  432. return m_T*m_gamma[0][0];
  433. }
  434. long double squared_l2() const
  435. {
  436. long double s = 0;
  437. for (size_t i = m_complex_vector_size - 1; i >= 1; --i)
  438. {
  439. s += (m_gamma[i][0]*m_gamma[i][0] + m_gamma[i][1]*m_gamma[i][1]);
  440. }
  441. s *= 2;
  442. s += m_gamma[0][0]*m_gamma[0][0];
  443. return s*m_T;
  444. }
  445. ~cardinal_trigonometric_detail()
  446. {
  447. if (m_gamma)
  448. {
  449. fftwl_free(m_gamma);
  450. m_gamma = nullptr;
  451. }
  452. }
  453. private:
  454. long double m_t0;
  455. long double m_h;
  456. long double m_T;
  457. fftwl_complex* m_gamma;
  458. size_t m_complex_vector_size;
  459. };
  460. #ifdef BOOST_HAS_FLOAT128
  461. template<>
  462. class cardinal_trigonometric_detail<__float128> {
  463. public:
  464. cardinal_trigonometric_detail(const __float128* data, size_t length, __float128 t0, __float128 h) : m_t0{t0}, m_h{h}
  465. {
  466. if (length == 0)
  467. {
  468. throw std::logic_error("At least one sample is required.");
  469. }
  470. if (h <= 0)
  471. {
  472. throw std::logic_error("The step size must be > 0");
  473. }
  474. m_T = m_h*length;
  475. m_complex_vector_size = length/2 + 1;
  476. m_gamma = fftwq_alloc_complex(m_complex_vector_size);
  477. fftwq_plan plan = fftwq_plan_dft_r2c_1d(length, reinterpret_cast<__float128*>(const_cast<__float128*>(data)), m_gamma, FFTW_ESTIMATE);
  478. if (!plan)
  479. {
  480. throw std::logic_error("A null fftw plan was created.");
  481. }
  482. fftwq_execute(plan);
  483. fftwq_destroy_plan(plan);
  484. __float128 denom = length;
  485. for (size_t k = 0; k < m_complex_vector_size; ++k)
  486. {
  487. m_gamma[k][0] /= denom;
  488. m_gamma[k][1] /= denom;
  489. }
  490. if (length % 2 == 0)
  491. {
  492. m_gamma[m_complex_vector_size -1][0] /= 2;
  493. }
  494. }
  495. cardinal_trigonometric_detail(const cardinal_trigonometric_detail& old) = delete;
  496. cardinal_trigonometric_detail& operator=(const cardinal_trigonometric_detail&) = delete;
  497. cardinal_trigonometric_detail(cardinal_trigonometric_detail &&) = delete;
  498. __float128 operator()(__float128 t) const
  499. {
  500. using std::sin;
  501. using std::cos;
  502. using boost::math::constants::two_pi;
  503. using std::exp;
  504. __float128 s = m_gamma[0][0];
  505. __float128 x = two_pi<__float128>()*(t - m_t0)/m_T;
  506. fftwq_complex z;
  507. z[0] = cosq(x);
  508. z[1] = sinq(x);
  509. fftwq_complex b{0, 0};
  510. fftwq_complex u;
  511. for (size_t k = m_complex_vector_size - 1; k >= 1; --k)
  512. {
  513. u[0] = b[0]*z[0] - b[1]*z[1];
  514. u[1] = b[0]*z[1] + b[1]*z[0];
  515. b[0] = m_gamma[k][0] + u[0];
  516. b[1] = m_gamma[k][1] + u[1];
  517. }
  518. s += 2*(b[0]*z[0] - b[1]*z[1]);
  519. return s;
  520. }
  521. __float128 prime(__float128 t) const
  522. {
  523. using std::sin;
  524. using std::cos;
  525. using boost::math::constants::two_pi;
  526. using std::exp;
  527. __float128 x = two_pi<__float128>()*(t - m_t0)/m_T;
  528. fftwq_complex z;
  529. z[0] = cosq(x);
  530. z[1] = sinq(x);
  531. fftwq_complex b{0, 0};
  532. // u = b*z
  533. fftwq_complex u;
  534. for (size_t k = m_complex_vector_size - 1; k >= 1; --k)
  535. {
  536. u[0] = b[0]*z[0] - b[1]*z[1];
  537. u[1] = b[0]*z[1] + b[1]*z[0];
  538. b[0] = k*m_gamma[k][0] + u[0];
  539. b[1] = k*m_gamma[k][1] + u[1];
  540. }
  541. // b*z = (b[0]*z[0] - b[1]*z[1]) + i(b[1]*z[0] + b[0]*z[1])
  542. return -2*two_pi<__float128>()*(b[1]*z[0] + b[0]*z[1])/m_T;
  543. }
  544. __float128 double_prime(__float128 t) const
  545. {
  546. using std::sin;
  547. using std::cos;
  548. using boost::math::constants::two_pi;
  549. using std::exp;
  550. __float128 x = two_pi<__float128>()*(t - m_t0)/m_T;
  551. fftwq_complex z;
  552. z[0] = cosq(x);
  553. z[1] = sinq(x);
  554. fftwq_complex b{0, 0};
  555. // u = b*z
  556. fftwq_complex u;
  557. for (size_t k = m_complex_vector_size - 1; k >= 1; --k)
  558. {
  559. u[0] = b[0]*z[0] - b[1]*z[1];
  560. u[1] = b[0]*z[1] + b[1]*z[0];
  561. b[0] = k*k*m_gamma[k][0] + u[0];
  562. b[1] = k*k*m_gamma[k][1] + u[1];
  563. }
  564. // b*z = (b[0]*z[0] - b[1]*z[1]) + i(b[1]*z[0] + b[0]*z[1])
  565. return -2*two_pi<__float128>()*two_pi<__float128>()*(b[0]*z[0] - b[1]*z[1])/(m_T*m_T);
  566. }
  567. __float128 period() const
  568. {
  569. return m_T;
  570. }
  571. __float128 integrate() const
  572. {
  573. return m_T*m_gamma[0][0];
  574. }
  575. __float128 squared_l2() const
  576. {
  577. __float128 s = 0;
  578. for (size_t i = m_complex_vector_size - 1; i >= 1; --i)
  579. {
  580. s += (m_gamma[i][0]*m_gamma[i][0] + m_gamma[i][1]*m_gamma[i][1]);
  581. }
  582. s *= 2;
  583. s += m_gamma[0][0]*m_gamma[0][0];
  584. return s*m_T;
  585. }
  586. ~cardinal_trigonometric_detail()
  587. {
  588. if (m_gamma)
  589. {
  590. fftwq_free(m_gamma);
  591. m_gamma = nullptr;
  592. }
  593. }
  594. private:
  595. __float128 m_t0;
  596. __float128 m_h;
  597. __float128 m_T;
  598. fftwq_complex* m_gamma;
  599. size_t m_complex_vector_size;
  600. };
  601. #endif
  602. }}}}
  603. #endif