fft_wrapper.hpp 18 KB

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