fft_wrapper.hpp 23 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644
  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. check_align(in, out);
  143. Vector<ValueType> buff0(N0 + N1);
  144. Vector<ValueType> buff1(N0 + N1);
  145. Long rank = plan.M.size();
  146. if (rank <= 0) return;
  147. Long N = N0;
  148. if (fft_type == FFT_Type::C2R) {
  149. const Matrix<ValueType>& M = plan.M[rank - 1];
  150. transpose<ComplexType>(buff0.begin(), in.begin(), N / M.Dim(0), M.Dim(0) / 2);
  151. for (Long i = 0; i < rank - 1; 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>(buff1.begin(), buff0.begin(), N / howmany / 2, howmany);
  160. Matrix<ValueType> vi(N / M.Dim(0), M.Dim(0), buff1.begin(), false);
  161. Matrix<ValueType> vo(N / M.Dim(0), M.Dim(1), out.begin(), false);
  162. Matrix<ValueType>::GEMM(vo, vi, M);
  163. } else {
  164. memcopy(buff0.begin(), in.begin(), in.Dim());
  165. for (Long i = 0; i < rank; i++) {
  166. const Matrix<ValueType>& M = plan.M[i];
  167. Matrix<ValueType> vi(N / M.Dim(0), M.Dim(0), buff0.begin(), false);
  168. Matrix<ValueType> vo(N / M.Dim(0), M.Dim(1), buff1.begin(), false);
  169. Matrix<ValueType>::GEMM(vo, vi, M);
  170. N = N * M.Dim(1) / M.Dim(0);
  171. transpose<ComplexType>(buff0.begin(), buff1.begin(), N / M.Dim(1), M.Dim(1) / 2);
  172. }
  173. transpose<ComplexType>(out.begin(), buff0.begin(), N / howmany / 2, howmany);
  174. }
  175. }
  176. static void test() {
  177. Vector<Long> fft_dim;
  178. fft_dim.PushBack(2);
  179. fft_dim.PushBack(5);
  180. fft_dim.PushBack(3);
  181. Long howmany = 3;
  182. if (1){ // R2C, C2R
  183. FFT_Derived myfft0, myfft1;
  184. myfft0.Setup(FFT_Type::R2C, howmany, fft_dim);
  185. myfft1.Setup(FFT_Type::C2R, howmany, fft_dim);
  186. Vector<ValueType> v0(myfft0.Dim(0)), v1, v2;
  187. for (int i = 0; i < v0.Dim(); i++) v0[i] = 1 + i;
  188. myfft0.Execute(v0, v1);
  189. myfft1.Execute(v1, v2);
  190. { // Print error
  191. ValueType err = 0;
  192. SCTL_ASSERT(v0.Dim() == v2.Dim());
  193. for (Long i = 0; i < v0.Dim(); i++) err = std::max(err, fabs(v0[i] - v2[i]));
  194. std::cout<<"Error : "<<err<<'\n';
  195. }
  196. }
  197. std::cout<<'\n';
  198. { // C2C, C2C_INV
  199. FFT_Derived myfft0, myfft1;
  200. myfft0.Setup(FFT_Type::C2C, howmany, fft_dim);
  201. myfft1.Setup(FFT_Type::C2C_INV, howmany, fft_dim);
  202. Vector<ValueType> v0(myfft0.Dim(0)), v1, v2;
  203. for (int i = 0; i < v0.Dim(); i++) v0[i] = 1 + i;
  204. myfft0.Execute(v0, v1);
  205. myfft1.Execute(v1, v2);
  206. { // Print error
  207. ValueType err = 0;
  208. SCTL_ASSERT(v0.Dim() == v2.Dim());
  209. for (Long i = 0; i < v0.Dim(); i++) err = std::max(err, fabs(v0[i] - v2[i]));
  210. std::cout<<"Error : "<<err<<'\n';
  211. }
  212. }
  213. std::cout<<'\n';
  214. }
  215. protected:
  216. static Matrix<ValueType> fft_r2c(Long N0) {
  217. ValueType s = 1 / sqrt<ValueType>(N0);
  218. Long N1 = (N0 / 2 + 1);
  219. Matrix<ValueType> M(N0, 2 * N1);
  220. for (Long j = 0; j < N0; j++)
  221. for (Long i = 0; i < N1; i++) {
  222. M[j][2 * i + 0] = cos<ValueType>(2 * const_pi<ValueType>() * j * i / N0)*s;
  223. M[j][2 * i + 1] = -sin<ValueType>(2 * const_pi<ValueType>() * j * i / N0)*s;
  224. }
  225. return M;
  226. }
  227. static Matrix<ValueType> fft_c2c(Long N0) {
  228. ValueType s = 1 / sqrt<ValueType>(N0);
  229. Matrix<ValueType> M(2 * N0, 2 * N0);
  230. for (Long i = 0; i < N0; i++)
  231. for (Long j = 0; j < N0; j++) {
  232. M[2 * i + 0][2 * j + 0] = cos<ValueType>(2 * const_pi<ValueType>() * j * i / N0)*s;
  233. M[2 * i + 1][2 * j + 0] = sin<ValueType>(2 * const_pi<ValueType>() * j * i / N0)*s;
  234. M[2 * i + 0][2 * j + 1] = -sin<ValueType>(2 * const_pi<ValueType>() * j * i / N0)*s;
  235. M[2 * i + 1][2 * j + 1] = cos<ValueType>(2 * const_pi<ValueType>() * j * i / N0)*s;
  236. }
  237. return M;
  238. }
  239. static Matrix<ValueType> fft_c2r(Long N0) {
  240. ValueType s = 1 / sqrt<ValueType>(N0);
  241. Long N1 = (N0 / 2 + 1);
  242. Matrix<ValueType> M(2 * N1, N0);
  243. for (Long i = 0; i < N1; i++) {
  244. for (Long j = 0; j < N0; j++) {
  245. M[2 * i + 0][j] = 2 * cos<ValueType>(2 * const_pi<ValueType>() * j * i / N0)*s;
  246. M[2 * i + 1][j] = -2 * sin<ValueType>(2 * const_pi<ValueType>() * j * i / N0)*s;
  247. }
  248. }
  249. if (N1 > 0) {
  250. for (Long j = 0; j < N0; j++) {
  251. M[0][j] = M[0][j] * (ValueType)0.5;
  252. M[1][j] = M[1][j] * (ValueType)0.5;
  253. }
  254. }
  255. if (N0 % 2 == 0) {
  256. for (Long j = 0; j < N0; j++) {
  257. M[2 * N1 - 2][j] = M[2 * N1 - 2][j] * (ValueType)0.5;
  258. M[2 * N1 - 1][j] = M[2 * N1 - 1][j] * (ValueType)0.5;
  259. }
  260. }
  261. return M;
  262. }
  263. template <class T> static void transpose(Iterator<ValueType> out, ConstIterator<ValueType> in, Long N0, Long N1) {
  264. Matrix<T> M0(N0, N1, (Iterator<T>)in, false);
  265. Matrix<T> M1(N1, N0, (Iterator<T>)out, false);
  266. M1 = M0.Transpose();
  267. }
  268. static void check_align(const Vector<ValueType>& in, const Vector<ValueType>& out) {
  269. SCTL_ASSERT_MSG((((uintptr_t)& in[0]) & ((uintptr_t)(SCTL_MEM_ALIGN - 1))) == 0, "sctl::FFT: Input vector not aligned to " <<SCTL_MEM_ALIGN<<" bits!");
  270. SCTL_ASSERT_MSG((((uintptr_t)&out[0]) & ((uintptr_t)(SCTL_MEM_ALIGN - 1))) == 0, "sctl::FFT: Output vector not aligned to "<<SCTL_MEM_ALIGN<<" bits!");
  271. }
  272. StaticArray<Long,2> dim;
  273. FFT_Type fft_type;
  274. Long howmany;
  275. FFTPlan plan;
  276. };
  277. template <class ValueType> class FFT : public FFT_Generic<ValueType, FFT<ValueType>> {};
  278. #ifdef SCTL_HAVE_FFTW
  279. template <> class FFT<double> : public FFT_Generic<double, FFT<double>> {
  280. typedef double ValueType;
  281. public:
  282. ~FFT() { if (this->Dim(0) && this->Dim(1)) fftw_destroy_plan(plan); }
  283. void Setup(FFT_Type fft_type_, Long howmany_, const Vector<Long>& dim_vec) {
  284. if (Dim(0) && Dim(1)) fftw_destroy_plan(plan); plan = NULL;
  285. this->fft_type = fft_type_;
  286. this->howmany = howmany_;
  287. Long rank = dim_vec.Dim();
  288. Vector<int> dim_vec_(rank);
  289. for (Integer i = 0; i < rank; i++) {
  290. dim_vec_[i] = dim_vec[i];
  291. }
  292. Long N0 = 0, N1 = 0;
  293. { // Set N0, N1
  294. Long N = howmany;
  295. for (auto ni : dim_vec) N *= ni;
  296. if (fft_type == FFT_Type::R2C) {
  297. N0 = N;
  298. N1 = (N / dim_vec[rank - 1]) * (dim_vec[rank - 1] / 2 + 1) * 2;
  299. } else if (fft_type == FFT_Type::C2C) {
  300. N0 = N * 2;
  301. N1 = N * 2;
  302. } else if (fft_type == FFT_Type::C2C_INV) {
  303. N0 = N * 2;
  304. N1 = N * 2;
  305. } else if (fft_type == FFT_Type::C2R) {
  306. N0 = (N / dim_vec[rank - 1]) * (dim_vec[rank - 1] / 2 + 1) * 2;
  307. N1 = N;
  308. }
  309. this->dim[0] = N0;
  310. this->dim[1] = N1;
  311. }
  312. if (!N0 || !N1) return;
  313. Vector<ValueType> in(N0), out(N1);
  314. if (fft_type == FFT_Type::R2C) {
  315. 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);
  316. } else if (fft_type == FFT_Type::C2C) {
  317. 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);
  318. } else if (fft_type == FFT_Type::C2C_INV) {
  319. 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);
  320. } else if (fft_type == FFT_Type::C2R) {
  321. 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);
  322. }
  323. if (!plan) { // Build plan without FFTW_PRESERVE_INPUT
  324. if (fft_type == FFT_Type::R2C) {
  325. 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);
  326. } else if (fft_type == FFT_Type::C2C) {
  327. 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);
  328. } else if (fft_type == FFT_Type::C2C_INV) {
  329. 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);
  330. } else if (fft_type == FFT_Type::C2R) {
  331. 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);
  332. }
  333. tmp.Swap(in);
  334. } else {
  335. tmp.ReInit(0);
  336. }
  337. SCTL_ASSERT(plan);
  338. }
  339. void Execute(const Vector<ValueType>& in, Vector<ValueType>& out) const {
  340. // TODO: check alignment of in and out
  341. Long N0 = this->Dim(0);
  342. Long N1 = this->Dim(1);
  343. if (!N0 || !N1) return;
  344. SCTL_ASSERT_MSG(in.Dim() == N0, "FFT: Wrong input size.");
  345. if (out.Dim() != N1) out.ReInit(N1);
  346. check_align(in, out);
  347. ValueType s = 0;
  348. auto in_ptr = in.begin();
  349. if (tmp.Dim()) { // Save input
  350. assert(tmp.Dim() == N0);
  351. in_ptr = tmp.begin();
  352. tmp = in;
  353. }
  354. if (fft_type == FFT_Type::R2C) {
  355. s = 1 / sqrt<ValueType>(N0 / howmany);
  356. fftw_execute_dft_r2c(plan, (double*)&in_ptr[0], (fftw_complex*)&out[0]);
  357. } else if (fft_type == FFT_Type::C2C) {
  358. s = 1 / sqrt<ValueType>(N0 / howmany * (ValueType)0.5);
  359. fftw_execute_dft(plan, (fftw_complex*)&in_ptr[0], (fftw_complex*)&out[0]);
  360. } else if (fft_type == FFT_Type::C2C_INV) {
  361. s = 1 / sqrt<ValueType>(N1 / howmany * (ValueType)0.5);
  362. fftw_execute_dft(plan, (fftw_complex*)&in_ptr[0], (fftw_complex*)&out[0]);
  363. } else if (fft_type == FFT_Type::C2R) {
  364. s = 1 / sqrt<ValueType>(N1 / howmany);
  365. fftw_execute_dft_c2r(plan, (fftw_complex*)&in_ptr[0], (double*)&out[0]);
  366. }
  367. for (auto& x : out) x *= s;
  368. }
  369. private:
  370. mutable Vector<ValueType> tmp;
  371. fftw_plan plan;
  372. };
  373. #endif
  374. #ifdef SCTL_HAVE_FFTWF
  375. template <> class FFT<float> : public FFT_Generic<float, FFT<float>> {
  376. typedef float ValueType;
  377. public:
  378. ~FFT() { if (this->Dim(0) && this->Dim(1)) fftwf_destroy_plan(plan); }
  379. void Setup(FFT_Type fft_type_, Long howmany_, const Vector<Long>& dim_vec) {
  380. if (Dim(0) && Dim(1)) fftwf_destroy_plan(plan); plan = NULL;
  381. this->fft_type = fft_type_;
  382. this->howmany = howmany_;
  383. Long rank = dim_vec.Dim();
  384. Vector<int> dim_vec_(rank);
  385. for (Integer i = 0; i < rank; i++) {
  386. dim_vec_[i] = dim_vec[i];
  387. }
  388. Long N0, N1;
  389. { // Set N0, N1
  390. Long N = howmany;
  391. for (auto ni : dim_vec) N *= ni;
  392. if (fft_type == FFT_Type::R2C) {
  393. N0 = N;
  394. N1 = (N / dim_vec[rank - 1]) * (dim_vec[rank - 1] / 2 + 1) * 2;
  395. } else if (fft_type == FFT_Type::C2C) {
  396. N0 = N * 2;
  397. N1 = N * 2;
  398. } else if (fft_type == FFT_Type::C2C_INV) {
  399. N0 = N * 2;
  400. N1 = N * 2;
  401. } else if (fft_type == FFT_Type::C2R) {
  402. N0 = (N / dim_vec[rank - 1]) * (dim_vec[rank - 1] / 2 + 1) * 2;
  403. N1 = N;
  404. }
  405. this->dim[0] = N0;
  406. this->dim[1] = N1;
  407. }
  408. if (!N0 || !N1) return;
  409. Vector<ValueType> in (N0), out(N1);
  410. if (fft_type == FFT_Type::R2C) {
  411. 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);
  412. } else if (fft_type == FFT_Type::C2C) {
  413. 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);
  414. } else if (fft_type == FFT_Type::C2C_INV) {
  415. 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);
  416. } else if (fft_type == FFT_Type::C2R) {
  417. 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);
  418. }
  419. if (!plan) { // Build plan without FFTW_PRESERVE_INPUT
  420. if (fft_type == FFT_Type::R2C) {
  421. 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);
  422. } else if (fft_type == FFT_Type::C2C) {
  423. 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);
  424. } else if (fft_type == FFT_Type::C2C_INV) {
  425. 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);
  426. } else if (fft_type == FFT_Type::C2R) {
  427. 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);
  428. }
  429. tmp.Swap(in);
  430. } else {
  431. tmp.ReInit(0);
  432. }
  433. SCTL_ASSERT(plan);
  434. }
  435. void Execute(const Vector<ValueType>& in, Vector<ValueType>& out) const {
  436. Long N0 = this->Dim(0);
  437. Long N1 = this->Dim(1);
  438. if (!N0 || !N1) return;
  439. SCTL_ASSERT_MSG(in.Dim() == N0, "FFT: Wrong input size.");
  440. if (out.Dim() != N1) out.ReInit(N1);
  441. check_align(in, out);
  442. ValueType s = 0;
  443. auto in_ptr = in.begin();
  444. if (tmp.Dim()) { // Save input
  445. assert(tmp.Dim() == N0);
  446. in_ptr = tmp.begin();
  447. tmp = in;
  448. }
  449. if (fft_type == FFT_Type::R2C) {
  450. s = 1 / sqrt<ValueType>(N0 / howmany);
  451. fftwf_execute_dft_r2c(plan, (float*)&in_ptr[0], (fftwf_complex*)&out[0]);
  452. } else if (fft_type == FFT_Type::C2C) {
  453. s = 1 / sqrt<ValueType>(N0 / howmany * (ValueType)0.5);
  454. fftwf_execute_dft(plan, (fftwf_complex*)&in_ptr[0], (fftwf_complex*)&out[0]);
  455. } else if (fft_type == FFT_Type::C2C_INV) {
  456. s = 1 / sqrt<ValueType>(N1 / howmany * (ValueType)0.5);
  457. fftwf_execute_dft(plan, (fftwf_complex*)&in_ptr[0], (fftwf_complex*)&out[0]);
  458. } else if (fft_type == FFT_Type::C2R) {
  459. s = 1 / sqrt<ValueType>(N1 / howmany);
  460. fftwf_execute_dft_c2r(plan, (fftwf_complex*)&in_ptr[0], (float*)&out[0]);
  461. }
  462. for (auto& x : out) x *= s;
  463. }
  464. private:
  465. mutable Vector<ValueType> tmp;
  466. fftwf_plan plan;
  467. };
  468. #endif
  469. #ifdef SCTL_HAVE_FFTWL
  470. template <> class FFT<long double> : public FFT_Generic<long double, FFT<long double>> {
  471. typedef long double ValueType;
  472. public:
  473. ~FFT() { if (this->Dim(0) && this->Dim(1)) fftwl_destroy_plan(plan); }
  474. void Setup(FFT_Type fft_type_, Long howmany_, const Vector<Long>& dim_vec) {
  475. if (Dim(0) && Dim(1)) fftwl_destroy_plan(plan); plan = NULL;
  476. this->fft_type = fft_type_;
  477. this->howmany = howmany_;
  478. Long rank = dim_vec.Dim();
  479. Vector<int> dim_vec_(rank);
  480. for (Integer i = 0; i < rank; i++) dim_vec_[i] = dim_vec[i];
  481. Long N0, N1;
  482. { // Set N0, N1
  483. Long N = howmany;
  484. for (auto ni : dim_vec) N *= ni;
  485. if (fft_type == FFT_Type::R2C) {
  486. N0 = N;
  487. N1 = (N / dim_vec[rank - 1]) * (dim_vec[rank - 1] / 2 + 1) * 2;
  488. } else if (fft_type == FFT_Type::C2C) {
  489. N0 = N * 2;
  490. N1 = N * 2;
  491. } else if (fft_type == FFT_Type::C2C_INV) {
  492. N0 = N * 2;
  493. N1 = N * 2;
  494. } else if (fft_type == FFT_Type::C2R) {
  495. N0 = (N / dim_vec[rank - 1]) * (dim_vec[rank - 1] / 2 + 1) * 2;
  496. N1 = N;
  497. }
  498. this->dim[0] = N0;
  499. this->dim[1] = N1;
  500. }
  501. if (!N0 || !N1) return;
  502. Vector<ValueType> in (N0), out(N1);
  503. if (fft_type == FFT_Type::R2C) {
  504. 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);
  505. } else if (fft_type == FFT_Type::C2C) {
  506. 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);
  507. } else if (fft_type == FFT_Type::C2C_INV) {
  508. 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);
  509. } else if (fft_type == FFT_Type::C2R) {
  510. 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);
  511. }
  512. if (!plan) { // Build plan without FFTW_PRESERVE_INPUT
  513. if (fft_type == FFT_Type::R2C) {
  514. 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);
  515. } else if (fft_type == FFT_Type::C2C) {
  516. 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);
  517. } else if (fft_type == FFT_Type::C2C_INV) {
  518. 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);
  519. } else if (fft_type == FFT_Type::C2R) {
  520. 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);
  521. }
  522. tmp.Swap(in);
  523. } else {
  524. tmp.ReInit(0);
  525. }
  526. SCTL_ASSERT(plan);
  527. }
  528. void Execute(const Vector<ValueType>& in, Vector<ValueType>& out) const {
  529. Long N0 = this->Dim(0);
  530. Long N1 = this->Dim(1);
  531. if (!N0 || !N1) return;
  532. SCTL_ASSERT_MSG(in.Dim() == N0, "FFT: Wrong input size.");
  533. if (out.Dim() != N1) out.ReInit(N1);
  534. check_align(in, out);
  535. ValueType s = 0;
  536. auto in_ptr = in.begin();
  537. if (tmp.Dim()) { // Save input
  538. assert(tmp.Dim() == N0);
  539. in_ptr = tmp.begin();
  540. tmp = in;
  541. }
  542. if (fft_type == FFT_Type::R2C) {
  543. s = 1 / sqrt<ValueType>(N0 / howmany);
  544. fftwl_execute_dft_r2c(plan, (long double*)&in_ptr[0], (fftwl_complex*)&out[0]);
  545. } else if (fft_type == FFT_Type::C2C) {
  546. s = 1 / sqrt<ValueType>(N0 / howmany * (ValueType)0.5);
  547. fftwl_execute_dft(plan, (fftwl_complex*)&in_ptr[0], (fftwl_complex*)&out[0]);
  548. } else if (fft_type == FFT_Type::C2C_INV) {
  549. s = 1 / sqrt<ValueType>(N1 / howmany * (ValueType)0.5);
  550. fftwl_execute_dft(plan, (fftwl_complex*)&in_ptr[0], (fftwl_complex*)&out[0]);
  551. } else if (fft_type == FFT_Type::C2R) {
  552. s = 1 / sqrt<ValueType>(N1 / howmany);
  553. fftwl_execute_dft_c2r(plan, (fftwl_complex*)&in_ptr[0], (long double*)&out[0]);
  554. }
  555. for (auto& x : out) x *= s;
  556. }
  557. private:
  558. mutable Vector<ValueType> tmp;
  559. fftwl_plan plan;
  560. };
  561. #endif
  562. } // end namespace
  563. #endif //_SCTL_FFT_WRAPPER_