fft_wrapper.hpp 19 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571
  1. #ifndef _SCTL_FFT_WRAPPER_
  2. #define _SCTL_FFT_WRAPPER_
  3. #include <cmath>
  4. #include <cassert>
  5. #include <cstdlib>
  6. #include <vector>
  7. #if defined(SCTL_HAVE_FFTW) || defined(SCTL_HAVE_FFTWF)
  8. #include <fftw3.h>
  9. #ifdef SCTL_FFTW3_MKL
  10. #include <fftw3_mkl.h>
  11. #endif
  12. #endif
  13. #include SCTL_INCLUDE(common.hpp)
  14. #include SCTL_INCLUDE(mem_mgr.hpp)
  15. #include SCTL_INCLUDE(matrix.hpp)
  16. namespace SCTL_NAMESPACE {
  17. template <class ValueType> class Complex {
  18. public:
  19. Complex<ValueType>(ValueType r=0, ValueType i=0) : real(r), imag(i) {}
  20. Complex<ValueType> operator*(const Complex<ValueType>& x) const {
  21. Complex<ValueType> z;
  22. z.real = real * x.real - imag * x.imag;
  23. z.imag = imag * x.real + real * x.imag;
  24. return z;
  25. }
  26. Complex<ValueType> operator*(const ValueType& x) const {
  27. Complex<ValueType> z;
  28. z.real = real * x;
  29. z.imag = imag * x;
  30. return z;
  31. }
  32. Complex<ValueType> operator+(const Complex<ValueType>& x) const {
  33. Complex<ValueType> z;
  34. z.real = real + x.real;
  35. z.imag = imag + x.imag;
  36. return z;
  37. }
  38. Complex<ValueType> operator+(const ValueType& x) const {
  39. Complex<ValueType> z;
  40. z.real = real + x;
  41. z.imag = imag;
  42. return z;
  43. }
  44. Complex<ValueType> operator-(const Complex<ValueType>& x) const {
  45. Complex<ValueType> z;
  46. z.real = real - x.real;
  47. z.imag = imag - x.imag;
  48. return z;
  49. }
  50. Complex<ValueType> operator-(const ValueType& x) const {
  51. Complex<ValueType> z;
  52. z.real = real - x;
  53. z.imag = imag;
  54. return z;
  55. }
  56. Complex<ValueType> operator-() const {
  57. Complex<ValueType> z;
  58. z.real = -real;
  59. z.imag = -imag;
  60. return z;
  61. }
  62. Complex<ValueType> conj() const {
  63. Complex<ValueType> z;
  64. z.real = real;
  65. z.imag = -imag;
  66. return z;
  67. }
  68. ValueType real;
  69. ValueType imag;
  70. };
  71. template <class ValueType> Complex<ValueType> operator*(const ValueType& x, const Complex<ValueType>& y){
  72. Complex<ValueType> z;
  73. z.real = y.real * x;
  74. z.imag = y.imag * x;
  75. return z;
  76. }
  77. template <class ValueType> Complex<ValueType> operator+(const ValueType& x, const Complex<ValueType>& y){
  78. Complex<ValueType> z;
  79. z.real = y.real + x;
  80. z.imag = y.imag;
  81. return z;
  82. }
  83. template <class ValueType> Complex<ValueType> operator-(const ValueType& x, const Complex<ValueType>& y){
  84. Complex<ValueType> z;
  85. z.real = y.real - x;
  86. z.imag = y.imag;
  87. return z;
  88. }
  89. enum class FFT_Type {R2C, C2C, C2C_INV, C2R};
  90. template <class ValueType, class FFT_Derived> class FFT_Generic {
  91. typedef Complex<ValueType> ComplexType;
  92. struct FFTPlan {
  93. std::vector<Matrix<ValueType>> M;
  94. };
  95. public:
  96. FFT_Generic() {
  97. dim[0] = 0;
  98. dim[1] = 0;
  99. }
  100. FFT_Generic(const FFT_Generic&) {
  101. dim[0]=0;
  102. dim[1]=0;
  103. }
  104. FFT_Generic& operator=(const FFT_Generic&) {
  105. dim[0]=0;
  106. dim[1]=0;
  107. return *this;
  108. };
  109. Long Dim(Integer i) const {
  110. return dim[i];
  111. }
  112. void Setup(FFT_Type fft_type_, Long howmany_, const Vector<Long>& dim_vec) {
  113. Long rank = dim_vec.Dim();
  114. fft_type = fft_type_;
  115. howmany = howmany_;
  116. plan.M.resize(0);
  117. if (fft_type == FFT_Type::R2C) {
  118. plan.M.push_back(fft_r2c(dim_vec[rank - 1]));
  119. for (Long i = rank - 2; i >= 0; i--) plan.M.push_back(fft_c2c(dim_vec[i]));
  120. } else if (fft_type == FFT_Type::C2C) {
  121. for (Long i = rank - 1; i >= 0; i--) plan.M.push_back(fft_c2c(dim_vec[i]));
  122. } else if (fft_type == FFT_Type::C2C_INV) {
  123. for (Long i = rank - 1; i >= 0; i--) plan.M.push_back(fft_c2c(dim_vec[i]).Transpose());
  124. } else if (fft_type == FFT_Type::C2R) {
  125. for (Long i = rank - 2; i >= 0; i--) plan.M.push_back(fft_c2c(dim_vec[i]).Transpose());
  126. plan.M.push_back(fft_c2r(dim_vec[rank - 1]));
  127. }
  128. Long N0 = howmany * 2;
  129. Long N1 = howmany * 2;
  130. for (const auto M : plan.M) {
  131. N0 = N0 * M.Dim(0) / 2;
  132. N1 = N1 * M.Dim(1) / 2;
  133. }
  134. dim[0] = N0;
  135. dim[1] = N1;
  136. }
  137. void Execute(const Vector<ValueType>& in, Vector<ValueType>& out) const {
  138. Long N0 = Dim(0);
  139. Long N1 = Dim(1);
  140. SCTL_ASSERT_MSG(in.Dim() == N0, "FFT: Wrong input size.");
  141. if (out.Dim() != N1) out.ReInit(N1);
  142. Vector<ValueType> buff0(N0 + N1);
  143. Vector<ValueType> buff1(N0 + N1);
  144. Long rank = plan.M.size();
  145. if (rank <= 0) return;
  146. Long N = N0;
  147. if (fft_type == FFT_Type::C2R) {
  148. const Matrix<ValueType>& M = plan.M[rank - 1];
  149. transpose<ComplexType>(buff0.begin(), in.begin(), N / M.Dim(0), M.Dim(0) / 2);
  150. for (Long i = 0; i < rank - 1; i++) {
  151. const Matrix<ValueType>& M = plan.M[i];
  152. Matrix<ValueType> vi(N / M.Dim(0), M.Dim(0), buff0.begin(), false);
  153. Matrix<ValueType> vo(N / M.Dim(0), M.Dim(1), buff1.begin(), false);
  154. Matrix<ValueType>::GEMM(vo, vi, M);
  155. N = N * M.Dim(1) / M.Dim(0);
  156. transpose<ComplexType>(buff0.begin(), buff1.begin(), N / M.Dim(1), M.Dim(1) / 2);
  157. }
  158. transpose<ComplexType>(buff1.begin(), buff0.begin(), N / howmany / 2, howmany);
  159. Matrix<ValueType> vi(N / M.Dim(0), M.Dim(0), buff1.begin(), false);
  160. Matrix<ValueType> vo(N / M.Dim(0), M.Dim(1), out.begin(), false);
  161. Matrix<ValueType>::GEMM(vo, vi, M);
  162. } else {
  163. memcopy(buff0.begin(), in.begin(), in.Dim());
  164. for (Long i = 0; i < rank; i++) {
  165. const Matrix<ValueType>& M = plan.M[i];
  166. Matrix<ValueType> vi(N / M.Dim(0), M.Dim(0), buff0.begin(), false);
  167. Matrix<ValueType> vo(N / M.Dim(0), M.Dim(1), buff1.begin(), false);
  168. Matrix<ValueType>::GEMM(vo, vi, M);
  169. N = N * M.Dim(1) / M.Dim(0);
  170. transpose<ComplexType>(buff0.begin(), buff1.begin(), N / M.Dim(1), M.Dim(1) / 2);
  171. }
  172. transpose<ComplexType>(out.begin(), buff0.begin(), N / howmany / 2, howmany);
  173. }
  174. }
  175. static void test() {
  176. Vector<Long> fft_dim;
  177. fft_dim.PushBack(2);
  178. fft_dim.PushBack(5);
  179. fft_dim.PushBack(3);
  180. Long howmany = 3;
  181. if (1){ // R2C, C2R
  182. FFT_Derived myfft0, myfft1;
  183. myfft0.Setup(FFT_Type::R2C, howmany, fft_dim);
  184. myfft1.Setup(FFT_Type::C2R, howmany, fft_dim);
  185. Vector<ValueType> v0(myfft0.Dim(0)), v1, v2;
  186. for (int i = 0; i < v0.Dim(); i++) v0[i] = 1 + i;
  187. myfft0.Execute(v0, v1);
  188. myfft1.Execute(v1, v2);
  189. { // Print error
  190. ValueType err = 0;
  191. SCTL_ASSERT(v0.Dim() == v2.Dim());
  192. for (Long i = 0; i < v0.Dim(); i++) err = std::max(err, fabs(v0[i] - v2[i]));
  193. std::cout<<"Error : "<<err<<'\n';
  194. }
  195. }
  196. std::cout<<'\n';
  197. { // C2C, C2C_INV
  198. FFT_Derived myfft0, myfft1;
  199. myfft0.Setup(FFT_Type::C2C, howmany, fft_dim);
  200. myfft1.Setup(FFT_Type::C2C_INV, howmany, fft_dim);
  201. Vector<ValueType> v0(myfft0.Dim(0)), v1, v2;
  202. for (int i = 0; i < v0.Dim(); i++) v0[i] = 1 + i;
  203. myfft0.Execute(v0, v1);
  204. myfft1.Execute(v1, v2);
  205. { // Print error
  206. ValueType err = 0;
  207. SCTL_ASSERT(v0.Dim() == v2.Dim());
  208. for (Long i = 0; i < v0.Dim(); i++) err = std::max(err, fabs(v0[i] - v2[i]));
  209. std::cout<<"Error : "<<err<<'\n';
  210. }
  211. }
  212. std::cout<<'\n';
  213. }
  214. protected:
  215. static Matrix<ValueType> fft_r2c(Long N0) {
  216. ValueType s = 1 / sqrt<ValueType>(N0);
  217. Long N1 = (N0 / 2 + 1);
  218. Matrix<ValueType> M(N0, 2 * N1);
  219. for (Long j = 0; j < N0; j++)
  220. for (Long i = 0; i < N1; i++) {
  221. M[j][2 * i + 0] = cos<ValueType>(2 * const_pi<ValueType>() * j * i / N0)*s;
  222. M[j][2 * i + 1] = -sin<ValueType>(2 * const_pi<ValueType>() * j * i / N0)*s;
  223. }
  224. return M;
  225. }
  226. static Matrix<ValueType> fft_c2c(Long N0) {
  227. ValueType s = 1 / sqrt<ValueType>(N0);
  228. Matrix<ValueType> M(2 * N0, 2 * N0);
  229. for (Long i = 0; i < N0; i++)
  230. for (Long j = 0; j < N0; j++) {
  231. M[2 * i + 0][2 * j + 0] = cos<ValueType>(2 * const_pi<ValueType>() * j * i / N0)*s;
  232. M[2 * i + 1][2 * j + 0] = sin<ValueType>(2 * const_pi<ValueType>() * j * i / N0)*s;
  233. M[2 * i + 0][2 * j + 1] = -sin<ValueType>(2 * const_pi<ValueType>() * j * i / N0)*s;
  234. M[2 * i + 1][2 * j + 1] = cos<ValueType>(2 * const_pi<ValueType>() * j * i / N0)*s;
  235. }
  236. return M;
  237. }
  238. static Matrix<ValueType> fft_c2r(Long N0) {
  239. ValueType s = 1 / sqrt<ValueType>(N0);
  240. Long N1 = (N0 / 2 + 1);
  241. Matrix<ValueType> M(2 * N1, N0);
  242. for (Long i = 0; i < N1; i++) {
  243. for (Long j = 0; j < N0; j++) {
  244. M[2 * i + 0][j] = 2 * cos<ValueType>(2 * const_pi<ValueType>() * j * i / N0)*s;
  245. M[2 * i + 1][j] = -2 * sin<ValueType>(2 * const_pi<ValueType>() * j * i / N0)*s;
  246. }
  247. }
  248. if (N1 > 0) {
  249. for (Long j = 0; j < N0; j++) {
  250. M[0][j] = M[0][j] * (ValueType)0.5;
  251. M[1][j] = M[1][j] * (ValueType)0.5;
  252. }
  253. }
  254. if (N0 % 2 == 0) {
  255. for (Long j = 0; j < N0; j++) {
  256. M[2 * N1 - 2][j] = M[2 * N1 - 2][j] * (ValueType)0.5;
  257. M[2 * N1 - 1][j] = M[2 * N1 - 1][j] * (ValueType)0.5;
  258. }
  259. }
  260. return M;
  261. }
  262. template <class T> static void transpose(Iterator<ValueType> out, ConstIterator<ValueType> in, Long N0, Long N1) {
  263. Matrix<T> M0(N0, N1, (Iterator<T>)in, false);
  264. Matrix<T> M1(N1, N0, (Iterator<T>)out, false);
  265. M1 = M0.Transpose();
  266. }
  267. StaticArray<Long,2> dim;
  268. FFT_Type fft_type;
  269. Long howmany;
  270. FFTPlan plan;
  271. };
  272. template <class ValueType> class FFT : public FFT_Generic<ValueType, FFT<ValueType>> {};
  273. #ifdef SCTL_HAVE_FFTW
  274. template <> class FFT<double> : public FFT_Generic<double, FFT<double>> {
  275. typedef double ValueType;
  276. public:
  277. ~FFT() { if (this->Dim(0) && this->Dim(1)) fftw_destroy_plan(plan); }
  278. void Setup(FFT_Type fft_type_, Long howmany_, const Vector<Long>& dim_vec) {
  279. if (Dim(0) && Dim(1)) fftw_destroy_plan(plan);
  280. this->fft_type = fft_type_;
  281. this->howmany = howmany_;
  282. Long rank = dim_vec.Dim();
  283. Vector<int> dim_vec_(rank);
  284. for (Integer i = 0; i < rank; i++) {
  285. dim_vec_[i] = dim_vec[i];
  286. }
  287. Long N0 = 0, N1 = 0;
  288. { // Set N0, N1
  289. Long N = howmany;
  290. for (auto ni : dim_vec) N *= ni;
  291. if (fft_type == FFT_Type::R2C) {
  292. N0 = N;
  293. N1 = (N / dim_vec[rank - 1]) * (dim_vec[rank - 1] / 2 + 1) * 2;
  294. } else if (fft_type == FFT_Type::C2C) {
  295. N0 = N * 2;
  296. N1 = N * 2;
  297. } else if (fft_type == FFT_Type::C2C_INV) {
  298. N0 = N * 2;
  299. N1 = N * 2;
  300. } else if (fft_type == FFT_Type::C2R) {
  301. N0 = (N / dim_vec[rank - 1]) * (dim_vec[rank - 1] / 2 + 1) * 2;
  302. N1 = N;
  303. }
  304. this->dim[0] = N0;
  305. this->dim[1] = N1;
  306. }
  307. if (!N0 || !N1) return;
  308. Vector<ValueType> in(N0), out(N1);
  309. if (fft_type == FFT_Type::R2C) {
  310. plan = fftw_plan_many_dft_r2c(rank, &dim_vec_[0], howmany_, &in[0], NULL, 1, N0 / howmany, (fftw_complex*)&out[0], NULL, 1, N1 / 2 / howmany, FFTW_ESTIMATE | FFTW_PRESERVE_INPUT);
  311. } else if (fft_type == FFT_Type::C2C) {
  312. plan = fftw_plan_many_dft(rank, &dim_vec_[0], howmany_, (fftw_complex*)&in[0], NULL, 1, N0 / 2 / howmany, (fftw_complex*)&out[0], NULL, 1, N1 / 2 / howmany, FFTW_FORWARD, FFTW_ESTIMATE | FFTW_PRESERVE_INPUT);
  313. } else if (fft_type == FFT_Type::C2C_INV) {
  314. plan = fftw_plan_many_dft(rank, &dim_vec_[0], howmany_, (fftw_complex*)&in[0], NULL, 1, N0 / 2 / howmany, (fftw_complex*)&out[0], NULL, 1, N1 / 2 / howmany, FFTW_BACKWARD, FFTW_ESTIMATE | FFTW_PRESERVE_INPUT);
  315. } else if (fft_type == FFT_Type::C2R) {
  316. plan = fftw_plan_many_dft_c2r(rank, &dim_vec_[0], howmany_, (fftw_complex*)&in[0], NULL, 1, N0 / 2 / howmany, &out[0], NULL, 1, N1 / howmany, FFTW_ESTIMATE | FFTW_PRESERVE_INPUT);
  317. }
  318. SCTL_ASSERT(plan);
  319. }
  320. void Execute(const Vector<ValueType>& in, Vector<ValueType>& out) const {
  321. Long N0 = this->Dim(0);
  322. Long N1 = this->Dim(1);
  323. if (!N0 || !N1) return;
  324. SCTL_ASSERT_MSG(in.Dim() == N0, "FFT: Wrong input size.");
  325. if (out.Dim() != N1) out.ReInit(N1);
  326. ValueType s = 0;
  327. if (fft_type == FFT_Type::R2C) {
  328. s = 1 / sqrt<ValueType>(N0 / howmany);
  329. fftw_execute_dft_r2c(plan, (double*)&in[0], (fftw_complex*)&out[0]);
  330. } else if (fft_type == FFT_Type::C2C) {
  331. s = 1 / sqrt<ValueType>(N0 / howmany * (ValueType)0.5);
  332. fftw_execute_dft(plan, (fftw_complex*)&in[0], (fftw_complex*)&out[0]);
  333. } else if (fft_type == FFT_Type::C2C_INV) {
  334. s = 1 / sqrt<ValueType>(N1 / howmany * (ValueType)0.5);
  335. fftw_execute_dft(plan, (fftw_complex*)&in[0], (fftw_complex*)&out[0]);
  336. } else if (fft_type == FFT_Type::C2R) {
  337. s = 1 / sqrt<ValueType>(N1 / howmany);
  338. fftw_execute_dft_c2r(plan, (fftw_complex*)&in[0], (double*)&out[0]);
  339. }
  340. for (auto& x : out) x *= s;
  341. }
  342. private:
  343. fftw_plan plan;
  344. };
  345. #endif
  346. #ifdef SCTL_HAVE_FFTWF
  347. template <> class FFT<float> : public FFT_Generic<float, FFT<float>> {
  348. typedef float ValueType;
  349. public:
  350. ~FFT() { if (this->Dim(0) && this->Dim(1)) fftwf_destroy_plan(plan); }
  351. void Setup(FFT_Type fft_type_, Long howmany_, const Vector<Long>& dim_vec) {
  352. if (Dim(0) && Dim(1)) fftwf_destroy_plan(plan);
  353. this->fft_type = fft_type_;
  354. this->howmany = howmany_;
  355. Long rank = dim_vec.Dim();
  356. Vector<int> dim_vec_(rank);
  357. for (Integer i = 0; i < rank; i++) {
  358. dim_vec_[i] = dim_vec[i];
  359. }
  360. Long N0, N1;
  361. { // Set N0, N1
  362. Long N = howmany;
  363. for (auto ni : dim_vec) N *= ni;
  364. if (fft_type == FFT_Type::R2C) {
  365. N0 = N;
  366. N1 = (N / dim_vec[rank - 1]) * (dim_vec[rank - 1] / 2 + 1) * 2;
  367. } else if (fft_type == FFT_Type::C2C) {
  368. N0 = N * 2;
  369. N1 = N * 2;
  370. } else if (fft_type == FFT_Type::C2C_INV) {
  371. N0 = N * 2;
  372. N1 = N * 2;
  373. } else if (fft_type == FFT_Type::C2R) {
  374. N0 = (N / dim_vec[rank - 1]) * (dim_vec[rank - 1] / 2 + 1) * 2;
  375. N1 = N;
  376. }
  377. this->dim[0] = N0;
  378. this->dim[1] = N1;
  379. }
  380. if (!N0 || !N1) return;
  381. Vector<ValueType> in (N0), out(N1);
  382. if (fft_type == FFT_Type::R2C) {
  383. plan = fftwf_plan_many_dft_r2c(rank, &dim_vec_[0], howmany_, &in[0], NULL, 1, N0 / howmany, (fftwf_complex*)&out[0], NULL, 1, N1 / 2 / howmany, FFTW_ESTIMATE | FFTW_PRESERVE_INPUT);
  384. } else if (fft_type == FFT_Type::C2C) {
  385. plan = fftwf_plan_many_dft(rank, &dim_vec_[0], howmany_, (fftwf_complex*)&in[0], NULL, 1, N0 / 2 / howmany, (fftwf_complex*)&out[0], NULL, 1, N1 / 2 / howmany, FFTW_FORWARD, FFTW_ESTIMATE | FFTW_PRESERVE_INPUT);
  386. } else if (fft_type == FFT_Type::C2C_INV) {
  387. plan = fftwf_plan_many_dft(rank, &dim_vec_[0], howmany_, (fftwf_complex*)&in[0], NULL, 1, N0 / 2 / howmany, (fftwf_complex*)&out[0], NULL, 1, N1 / 2 / howmany, FFTW_BACKWARD, FFTW_ESTIMATE | FFTW_PRESERVE_INPUT);
  388. } else if (fft_type == FFT_Type::C2R) {
  389. plan = fftwf_plan_many_dft_c2r(rank, &dim_vec_[0], howmany_, (fftwf_complex*)&in[0], NULL, 1, N0 / 2 / howmany, &out[0], NULL, 1, N1 / howmany, FFTW_ESTIMATE | FFTW_PRESERVE_INPUT);
  390. }
  391. SCTL_ASSERT(plan);
  392. }
  393. void Execute(const Vector<ValueType>& in, Vector<ValueType>& out) const {
  394. Long N0 = this->Dim(0);
  395. Long N1 = this->Dim(1);
  396. if (!N0 || !N1) return;
  397. SCTL_ASSERT_MSG(in.Dim() == N0, "FFT: Wrong input size.");
  398. if (out.Dim() != N1) out.ReInit(N1);
  399. ValueType s = 0;
  400. if (fft_type == FFT_Type::R2C) {
  401. s = 1 / sqrt<ValueType>(N0 / howmany);
  402. fftwf_execute_dft_r2c(plan, (float*)&in[0], (fftwf_complex*)&out[0]);
  403. } else if (fft_type == FFT_Type::C2C) {
  404. s = 1 / sqrt<ValueType>(N0 / howmany * (ValueType)0.5);
  405. fftwf_execute_dft(plan, (fftwf_complex*)&in[0], (fftwf_complex*)&out[0]);
  406. } else if (fft_type == FFT_Type::C2C_INV) {
  407. s = 1 / sqrt<ValueType>(N1 / howmany * (ValueType)0.5);
  408. fftwf_execute_dft(plan, (fftwf_complex*)&in[0], (fftwf_complex*)&out[0]);
  409. } else if (fft_type == FFT_Type::C2R) {
  410. s = 1 / sqrt<ValueType>(N1 / howmany);
  411. fftwf_execute_dft_c2r(plan, (fftwf_complex*)&in[0], (float*)&out[0]);
  412. }
  413. for (auto& x : out) x *= s;
  414. }
  415. private:
  416. fftwf_plan plan;
  417. };
  418. #endif
  419. #ifdef SCTL_HAVE_FFTWL
  420. template <> class FFT<long double> : public FFT_Generic<long double, FFT<long double>> {
  421. typedef long double ValueType;
  422. public:
  423. ~FFT() { if (this->Dim(0) && this->Dim(1)) fftwl_destroy_plan(plan); }
  424. void Setup(FFT_Type fft_type_, Long howmany_, const Vector<Long>& dim_vec) {
  425. if (Dim(0) && Dim(1)) fftwl_destroy_plan(plan);
  426. this->fft_type = fft_type_;
  427. this->howmany = howmany_;
  428. Long rank = dim_vec.Dim();
  429. Vector<int> dim_vec_(rank);
  430. for (Integer i = 0; i < rank; i++) dim_vec_[i] = dim_vec[i];
  431. Long N0, N1;
  432. { // Set N0, N1
  433. Long N = howmany;
  434. for (auto ni : dim_vec) N *= ni;
  435. if (fft_type == FFT_Type::R2C) {
  436. N0 = N;
  437. N1 = (N / dim_vec[rank - 1]) * (dim_vec[rank - 1] / 2 + 1) * 2;
  438. } else if (fft_type == FFT_Type::C2C) {
  439. N0 = N * 2;
  440. N1 = N * 2;
  441. } else if (fft_type == FFT_Type::C2C_INV) {
  442. N0 = N * 2;
  443. N1 = N * 2;
  444. } else if (fft_type == FFT_Type::C2R) {
  445. N0 = (N / dim_vec[rank - 1]) * (dim_vec[rank - 1] / 2 + 1) * 2;
  446. N1 = N;
  447. }
  448. this->dim[0] = N0;
  449. this->dim[1] = N1;
  450. }
  451. if (!N0 || !N1) return;
  452. Vector<ValueType> in (N0), out(N1);
  453. if (fft_type == FFT_Type::R2C) {
  454. plan = fftwl_plan_many_dft_r2c(rank, &dim_vec_[0], howmany_, &in[0], NULL, 1, N0 / howmany, (fftwl_complex*)&out[0], NULL, 1, N1 / 2 / howmany, FFTW_ESTIMATE | FFTW_PRESERVE_INPUT);
  455. } else if (fft_type == FFT_Type::C2C) {
  456. plan = fftwl_plan_many_dft(rank, &dim_vec_[0], howmany_, (fftwl_complex*)&in[0], NULL, 1, N0 / 2 / howmany, (fftwl_complex*)&out[0], NULL, 1, N1 / 2 / howmany, FFTW_FORWARD, FFTW_ESTIMATE | FFTW_PRESERVE_INPUT);
  457. } else if (fft_type == FFT_Type::C2C_INV) {
  458. plan = fftwl_plan_many_dft(rank, &dim_vec_[0], howmany_, (fftwl_complex*)&in[0], NULL, 1, N0 / 2 / howmany, (fftwl_complex*)&out[0], NULL, 1, N1 / 2 / howmany, FFTW_BACKWARD, FFTW_ESTIMATE | FFTW_PRESERVE_INPUT);
  459. } else if (fft_type == FFT_Type::C2R) {
  460. plan = fftwl_plan_many_dft_c2r(rank, &dim_vec_[0], howmany_, (fftwl_complex*)&in[0], NULL, 1, N0 / 2 / howmany, &out[0], NULL, 1, N1 / howmany, FFTW_ESTIMATE | FFTW_PRESERVE_INPUT);
  461. }
  462. SCTL_ASSERT(plan);
  463. }
  464. void Execute(const Vector<ValueType>& in, Vector<ValueType>& out) const {
  465. Long N0 = this->Dim(0);
  466. Long N1 = this->Dim(1);
  467. if (!N0 || !N1) return;
  468. SCTL_ASSERT_MSG(in.Dim() == N0, "FFT: Wrong input size.");
  469. if (out.Dim() != N1) out.ReInit(N1);
  470. ValueType s = 0;
  471. if (fft_type == FFT_Type::R2C) {
  472. s = 1 / sqrt<ValueType>(N0 / howmany);
  473. fftwl_execute_dft_r2c(plan, (long double*)&in[0], (fftwl_complex*)&out[0]);
  474. } else if (fft_type == FFT_Type::C2C) {
  475. s = 1 / sqrt<ValueType>(N0 / howmany * (ValueType)0.5);
  476. fftwl_execute_dft(plan, (fftwl_complex*)&in[0], (fftwl_complex*)&out[0]);
  477. } else if (fft_type == FFT_Type::C2C_INV) {
  478. s = 1 / sqrt<ValueType>(N1 / howmany * (ValueType)0.5);
  479. fftwl_execute_dft(plan, (fftwl_complex*)&in[0], (fftwl_complex*)&out[0]);
  480. } else if (fft_type == FFT_Type::C2R) {
  481. s = 1 / sqrt<ValueType>(N1 / howmany);
  482. fftwl_execute_dft_c2r(plan, (fftwl_complex*)&in[0], (long double*)&out[0]);
  483. }
  484. for (auto& x : out) x *= s;
  485. }
  486. private:
  487. fftwl_plan plan;
  488. };
  489. #endif
  490. } // end namespace
  491. #endif //_SCTL_FFT_WRAPPER_