fft_wrapper.hpp 19 KB

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