fft_wrapper.hpp 25 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707
  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 {
  21. Complex<ValueType> z;
  22. z.real = -real;
  23. z.imag = -imag;
  24. return z;
  25. }
  26. Complex<ValueType> conj() const {
  27. Complex<ValueType> z;
  28. z.real = real;
  29. z.imag = -imag;
  30. return z;
  31. }
  32. bool operator==(const Complex<ValueType>& x) const {
  33. return real == x.real && imag == x.imag;
  34. }
  35. bool operator!=(const Complex<ValueType>& x) const {
  36. return !((*this) == x);;
  37. }
  38. template <class ScalarType> void operator+=(const Complex<ScalarType>& x) {
  39. (*this) = (*this) + x;
  40. }
  41. template <class ScalarType> void operator-=(const Complex<ScalarType>& x) {
  42. (*this) = (*this) - x;
  43. }
  44. template <class ScalarType> void operator*=(const Complex<ScalarType>& x) {
  45. (*this) = (*this) * x;
  46. }
  47. template <class ScalarType> void operator/=(const Complex<ScalarType>& x) {
  48. (*this) = (*this) / x;
  49. }
  50. template <class ScalarType> Complex<ValueType> operator+(const ScalarType& x) const {
  51. Complex<ValueType> z;
  52. z.real = real + x;
  53. z.imag = imag;
  54. return z;
  55. }
  56. template <class ScalarType> Complex<ValueType> operator-(const ScalarType& x) const {
  57. Complex<ValueType> z;
  58. z.real = real - x;
  59. z.imag = imag;
  60. return z;
  61. }
  62. template <class ScalarType> Complex<ValueType> operator*(const ScalarType& x) const {
  63. Complex<ValueType> z;
  64. z.real = real * x;
  65. z.imag = imag * x;
  66. return z;
  67. }
  68. template <class ScalarType> Complex<ValueType> operator/(const ScalarType& y) const {
  69. Complex<ValueType> z;
  70. z.real = real / y;
  71. z.imag = imag / y;
  72. return z;
  73. }
  74. Complex<ValueType> operator+(const Complex<ValueType>& x) const {
  75. Complex<ValueType> z;
  76. z.real = real + x.real;
  77. z.imag = imag + x.imag;
  78. return z;
  79. }
  80. Complex<ValueType> operator-(const Complex<ValueType>& x) const {
  81. Complex<ValueType> z;
  82. z.real = real - x.real;
  83. z.imag = imag - x.imag;
  84. return z;
  85. }
  86. Complex<ValueType> operator*(const Complex<ValueType>& x) const {
  87. Complex<ValueType> z;
  88. z.real = real * x.real - imag * x.imag;
  89. z.imag = imag * x.real + real * x.imag;
  90. return z;
  91. }
  92. Complex<ValueType> operator/(const Complex<ValueType>& y) const {
  93. Complex<ValueType> z;
  94. ValueType y_inv = 1 / (y.real * y.real + y.imag * y.imag);
  95. z.real = (y.real * real + y.imag * imag) * y_inv;
  96. z.imag = (y.real * imag - y.imag * real) * y_inv;
  97. return z;
  98. }
  99. ValueType real;
  100. ValueType imag;
  101. };
  102. template <class ScalarType, class ValueType> Complex<ValueType> operator*(const ScalarType& x, const Complex<ValueType>& y) {
  103. Complex<ValueType> z;
  104. z.real = y.real * x;
  105. z.imag = y.imag * x;
  106. return z;
  107. }
  108. template <class ScalarType, class ValueType> Complex<ValueType> operator+(const ScalarType& x, const Complex<ValueType>& y) {
  109. Complex<ValueType> z;
  110. z.real = y.real + x;
  111. z.imag = y.imag;
  112. return z;
  113. }
  114. template <class ScalarType, class ValueType> Complex<ValueType> operator-(const ScalarType& x, const Complex<ValueType>& y) {
  115. Complex<ValueType> z;
  116. z.real = y.real - x;
  117. z.imag = y.imag;
  118. return z;
  119. }
  120. template <class ScalarType, class ValueType> Complex<ValueType> operator/(const ScalarType& x, const Complex<ValueType>& y) {
  121. Complex<ValueType> z;
  122. ValueType y_inv = 1 / (y.real * y.real + y.imag * y.imag);
  123. z.real = (y.real * x) * y_inv;
  124. z.imag = -(y.imag * x) * y_inv;
  125. return z;
  126. }
  127. enum class FFT_Type {R2C, C2C, C2C_INV, C2R};
  128. template <class ValueType, class FFT_Derived> class FFT_Generic {
  129. typedef Complex<ValueType> ComplexType;
  130. struct FFTPlan {
  131. std::vector<Matrix<ValueType>> M;
  132. };
  133. public:
  134. FFT_Generic() {
  135. dim[0] = 0;
  136. dim[1] = 0;
  137. }
  138. FFT_Generic(const FFT_Generic&) {
  139. dim[0]=0;
  140. dim[1]=0;
  141. }
  142. FFT_Generic& operator=(const FFT_Generic&) {
  143. dim[0]=0;
  144. dim[1]=0;
  145. return *this;
  146. };
  147. Long Dim(Integer i) const {
  148. return dim[i];
  149. }
  150. void Setup(FFT_Type fft_type_, Long howmany_, const Vector<Long>& dim_vec) {
  151. Long rank = dim_vec.Dim();
  152. fft_type = fft_type_;
  153. howmany = howmany_;
  154. plan.M.resize(0);
  155. if (fft_type == FFT_Type::R2C) {
  156. plan.M.push_back(fft_r2c(dim_vec[rank - 1]));
  157. for (Long i = rank - 2; i >= 0; i--) plan.M.push_back(fft_c2c(dim_vec[i]));
  158. } else if (fft_type == FFT_Type::C2C) {
  159. for (Long i = rank - 1; i >= 0; i--) plan.M.push_back(fft_c2c(dim_vec[i]));
  160. } else if (fft_type == FFT_Type::C2C_INV) {
  161. for (Long i = rank - 1; i >= 0; i--) plan.M.push_back(fft_c2c(dim_vec[i]).Transpose());
  162. } else if (fft_type == FFT_Type::C2R) {
  163. for (Long i = rank - 2; i >= 0; i--) plan.M.push_back(fft_c2c(dim_vec[i]).Transpose());
  164. plan.M.push_back(fft_c2r(dim_vec[rank - 1]));
  165. }
  166. Long N0 = howmany * 2;
  167. Long N1 = howmany * 2;
  168. for (const auto M : plan.M) {
  169. N0 = N0 * M.Dim(0) / 2;
  170. N1 = N1 * M.Dim(1) / 2;
  171. }
  172. dim[0] = N0;
  173. dim[1] = N1;
  174. }
  175. void Execute(const Vector<ValueType>& in, Vector<ValueType>& out) const {
  176. Long N0 = Dim(0);
  177. Long N1 = Dim(1);
  178. SCTL_ASSERT_MSG(in.Dim() == N0, "FFT: Wrong input size.");
  179. if (out.Dim() != N1) out.ReInit(N1);
  180. check_align(in, out);
  181. Vector<ValueType> buff0(N0 + N1);
  182. Vector<ValueType> buff1(N0 + N1);
  183. Long rank = plan.M.size();
  184. if (rank <= 0) return;
  185. Long N = N0;
  186. if (fft_type == FFT_Type::C2R) {
  187. const Matrix<ValueType>& M = plan.M[rank - 1];
  188. transpose<ComplexType>(buff0.begin(), in.begin(), N / M.Dim(0), M.Dim(0) / 2);
  189. for (Long i = 0; i < rank - 1; i++) {
  190. const Matrix<ValueType>& M = plan.M[i];
  191. Matrix<ValueType> vi(N / M.Dim(0), M.Dim(0), buff0.begin(), false);
  192. Matrix<ValueType> vo(N / M.Dim(0), M.Dim(1), buff1.begin(), false);
  193. Matrix<ValueType>::GEMM(vo, vi, M);
  194. N = N * M.Dim(1) / M.Dim(0);
  195. transpose<ComplexType>(buff0.begin(), buff1.begin(), N / M.Dim(1), M.Dim(1) / 2);
  196. }
  197. transpose<ComplexType>(buff1.begin(), buff0.begin(), N / howmany / 2, howmany);
  198. Matrix<ValueType> vi(N / M.Dim(0), M.Dim(0), buff1.begin(), false);
  199. Matrix<ValueType> vo(N / M.Dim(0), M.Dim(1), out.begin(), false);
  200. Matrix<ValueType>::GEMM(vo, vi, M);
  201. } else {
  202. memcopy(buff0.begin(), in.begin(), in.Dim());
  203. for (Long i = 0; i < rank; i++) {
  204. const Matrix<ValueType>& M = plan.M[i];
  205. Matrix<ValueType> vi(N / M.Dim(0), M.Dim(0), buff0.begin(), false);
  206. Matrix<ValueType> vo(N / M.Dim(0), M.Dim(1), buff1.begin(), false);
  207. Matrix<ValueType>::GEMM(vo, vi, M);
  208. N = N * M.Dim(1) / M.Dim(0);
  209. transpose<ComplexType>(buff0.begin(), buff1.begin(), N / M.Dim(1), M.Dim(1) / 2);
  210. }
  211. transpose<ComplexType>(out.begin(), buff0.begin(), N / howmany / 2, howmany);
  212. }
  213. }
  214. static void test() {
  215. Vector<Long> fft_dim;
  216. fft_dim.PushBack(2);
  217. fft_dim.PushBack(5);
  218. fft_dim.PushBack(3);
  219. Long howmany = 3;
  220. if (1){ // R2C, C2R
  221. FFT_Derived myfft0, myfft1;
  222. myfft0.Setup(FFT_Type::R2C, howmany, fft_dim);
  223. myfft1.Setup(FFT_Type::C2R, howmany, fft_dim);
  224. Vector<ValueType> v0(myfft0.Dim(0)), v1, v2;
  225. for (int i = 0; i < v0.Dim(); i++) v0[i] = 1 + i;
  226. myfft0.Execute(v0, v1);
  227. myfft1.Execute(v1, v2);
  228. { // Print error
  229. ValueType err = 0;
  230. SCTL_ASSERT(v0.Dim() == v2.Dim());
  231. for (Long i = 0; i < v0.Dim(); i++) err = std::max(err, fabs(v0[i] - v2[i]));
  232. std::cout<<"Error : "<<err<<'\n';
  233. }
  234. }
  235. std::cout<<'\n';
  236. { // C2C, C2C_INV
  237. FFT_Derived myfft0, myfft1;
  238. myfft0.Setup(FFT_Type::C2C, howmany, fft_dim);
  239. myfft1.Setup(FFT_Type::C2C_INV, howmany, fft_dim);
  240. Vector<ValueType> v0(myfft0.Dim(0)), v1, v2;
  241. for (int i = 0; i < v0.Dim(); i++) v0[i] = 1 + i;
  242. myfft0.Execute(v0, v1);
  243. myfft1.Execute(v1, v2);
  244. { // Print error
  245. ValueType err = 0;
  246. SCTL_ASSERT(v0.Dim() == v2.Dim());
  247. for (Long i = 0; i < v0.Dim(); i++) err = std::max(err, fabs(v0[i] - v2[i]));
  248. std::cout<<"Error : "<<err<<'\n';
  249. }
  250. }
  251. std::cout<<'\n';
  252. }
  253. protected:
  254. static Matrix<ValueType> fft_r2c(Long N0) {
  255. ValueType s = 1 / sqrt<ValueType>(N0);
  256. Long N1 = (N0 / 2 + 1);
  257. Matrix<ValueType> M(N0, 2 * N1);
  258. for (Long j = 0; j < N0; j++)
  259. for (Long i = 0; i < N1; i++) {
  260. M[j][2 * i + 0] = cos<ValueType>(2 * const_pi<ValueType>() * j * i / N0)*s;
  261. M[j][2 * i + 1] = -sin<ValueType>(2 * const_pi<ValueType>() * j * i / N0)*s;
  262. }
  263. return M;
  264. }
  265. static Matrix<ValueType> fft_c2c(Long N0) {
  266. ValueType s = 1 / sqrt<ValueType>(N0);
  267. Matrix<ValueType> M(2 * N0, 2 * N0);
  268. for (Long i = 0; i < N0; i++)
  269. for (Long j = 0; j < N0; j++) {
  270. M[2 * i + 0][2 * j + 0] = cos<ValueType>(2 * const_pi<ValueType>() * j * i / N0)*s;
  271. M[2 * i + 1][2 * j + 0] = sin<ValueType>(2 * const_pi<ValueType>() * j * i / N0)*s;
  272. M[2 * i + 0][2 * j + 1] = -sin<ValueType>(2 * const_pi<ValueType>() * j * i / N0)*s;
  273. M[2 * i + 1][2 * j + 1] = cos<ValueType>(2 * const_pi<ValueType>() * j * i / N0)*s;
  274. }
  275. return M;
  276. }
  277. static Matrix<ValueType> fft_c2r(Long N0) {
  278. ValueType s = 1 / sqrt<ValueType>(N0);
  279. Long N1 = (N0 / 2 + 1);
  280. Matrix<ValueType> M(2 * N1, N0);
  281. for (Long i = 0; i < N1; i++) {
  282. for (Long j = 0; j < N0; j++) {
  283. M[2 * i + 0][j] = 2 * cos<ValueType>(2 * const_pi<ValueType>() * j * i / N0)*s;
  284. M[2 * i + 1][j] = -2 * sin<ValueType>(2 * const_pi<ValueType>() * j * i / N0)*s;
  285. }
  286. }
  287. if (N1 > 0) {
  288. for (Long j = 0; j < N0; j++) {
  289. M[0][j] = M[0][j] * (ValueType)0.5;
  290. M[1][j] = M[1][j] * (ValueType)0.5;
  291. }
  292. }
  293. if (N0 % 2 == 0) {
  294. for (Long j = 0; j < N0; j++) {
  295. M[2 * N1 - 2][j] = M[2 * N1 - 2][j] * (ValueType)0.5;
  296. M[2 * N1 - 1][j] = M[2 * N1 - 1][j] * (ValueType)0.5;
  297. }
  298. }
  299. return M;
  300. }
  301. template <class T> static void transpose(Iterator<ValueType> out, ConstIterator<ValueType> in, Long N0, Long N1) {
  302. Matrix<T> M0(N0, N1, (Iterator<T>)in, false);
  303. Matrix<T> M1(N1, N0, (Iterator<T>)out, false);
  304. M1 = M0.Transpose();
  305. }
  306. static void check_align(const Vector<ValueType>& in, const Vector<ValueType>& out) {
  307. //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!");
  308. //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!");
  309. // TODO: copy to auxiliary array if unaligned
  310. }
  311. StaticArray<Long,2> dim;
  312. FFT_Type fft_type;
  313. Long howmany;
  314. FFTPlan plan;
  315. };
  316. template <class ValueType> class FFT : public FFT_Generic<ValueType, FFT<ValueType>> {};
  317. #ifdef SCTL_HAVE_FFTW
  318. template <> class FFT<double> : public FFT_Generic<double, FFT<double>> {
  319. typedef double ValueType;
  320. public:
  321. ~FFT() { if (this->Dim(0) && this->Dim(1)) fftw_destroy_plan(plan); }
  322. void Setup(FFT_Type fft_type_, Long howmany_, const Vector<Long>& dim_vec) {
  323. if (Dim(0) && Dim(1)) fftw_destroy_plan(plan);
  324. this->fft_type = fft_type_;
  325. this->howmany = howmany_;
  326. plan = NULL;
  327. Long rank = dim_vec.Dim();
  328. Vector<int> dim_vec_(rank);
  329. for (Integer i = 0; i < rank; i++) {
  330. dim_vec_[i] = dim_vec[i];
  331. }
  332. Long N0 = 0, N1 = 0;
  333. { // Set N0, N1
  334. Long N = howmany;
  335. for (auto ni : dim_vec) N *= ni;
  336. if (fft_type == FFT_Type::R2C) {
  337. N0 = N;
  338. N1 = (N / dim_vec[rank - 1]) * (dim_vec[rank - 1] / 2 + 1) * 2;
  339. } else if (fft_type == FFT_Type::C2C) {
  340. N0 = N * 2;
  341. N1 = N * 2;
  342. } else if (fft_type == FFT_Type::C2C_INV) {
  343. N0 = N * 2;
  344. N1 = N * 2;
  345. } else if (fft_type == FFT_Type::C2R) {
  346. N0 = (N / dim_vec[rank - 1]) * (dim_vec[rank - 1] / 2 + 1) * 2;
  347. N1 = N;
  348. } else {
  349. N0 = 0;
  350. N1 = 0;
  351. }
  352. this->dim[0] = N0;
  353. this->dim[1] = N1;
  354. }
  355. if (!N0 || !N1) return;
  356. Vector<ValueType> in(N0), out(N1);
  357. if (fft_type == FFT_Type::R2C) {
  358. 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);
  359. } else if (fft_type == FFT_Type::C2C) {
  360. 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);
  361. } else if (fft_type == FFT_Type::C2C_INV) {
  362. 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);
  363. } else if (fft_type == FFT_Type::C2R) {
  364. 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);
  365. }
  366. if (!plan) { // Build plan without FFTW_PRESERVE_INPUT
  367. if (fft_type == FFT_Type::R2C) {
  368. 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);
  369. } else if (fft_type == FFT_Type::C2C) {
  370. 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);
  371. } else if (fft_type == FFT_Type::C2C_INV) {
  372. 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);
  373. } else if (fft_type == FFT_Type::C2R) {
  374. 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);
  375. }
  376. tmp.Swap(in);
  377. } else {
  378. tmp.ReInit(0);
  379. }
  380. SCTL_ASSERT(plan);
  381. }
  382. void Execute(const Vector<ValueType>& in, Vector<ValueType>& out) const {
  383. Long N0 = this->Dim(0);
  384. Long N1 = this->Dim(1);
  385. if (!N0 || !N1) return;
  386. SCTL_ASSERT_MSG(in.Dim() == N0, "FFT: Wrong input size.");
  387. if (out.Dim() != N1) out.ReInit(N1);
  388. check_align(in, out);
  389. ValueType s = 0;
  390. auto in_ptr = in.begin();
  391. if (tmp.Dim()) { // Save input
  392. assert(tmp.Dim() == N0);
  393. in_ptr = tmp.begin();
  394. tmp = in;
  395. }
  396. if (fft_type == FFT_Type::R2C) {
  397. s = 1 / sqrt<ValueType>(N0 / howmany);
  398. fftw_execute_dft_r2c(plan, (double*)&in_ptr[0], (fftw_complex*)&out[0]);
  399. } else if (fft_type == FFT_Type::C2C) {
  400. s = 1 / sqrt<ValueType>(N0 / howmany * (ValueType)0.5);
  401. fftw_execute_dft(plan, (fftw_complex*)&in_ptr[0], (fftw_complex*)&out[0]);
  402. } else if (fft_type == FFT_Type::C2C_INV) {
  403. s = 1 / sqrt<ValueType>(N1 / howmany * (ValueType)0.5);
  404. fftw_execute_dft(plan, (fftw_complex*)&in_ptr[0], (fftw_complex*)&out[0]);
  405. } else if (fft_type == FFT_Type::C2R) {
  406. s = 1 / sqrt<ValueType>(N1 / howmany);
  407. fftw_execute_dft_c2r(plan, (fftw_complex*)&in_ptr[0], (double*)&out[0]);
  408. }
  409. for (auto& x : out) x *= s;
  410. }
  411. private:
  412. mutable Vector<ValueType> tmp;
  413. fftw_plan plan;
  414. };
  415. #endif
  416. #ifdef SCTL_HAVE_FFTWF
  417. template <> class FFT<float> : public FFT_Generic<float, FFT<float>> {
  418. typedef float ValueType;
  419. public:
  420. ~FFT() { if (this->Dim(0) && this->Dim(1)) fftwf_destroy_plan(plan); }
  421. void Setup(FFT_Type fft_type_, Long howmany_, const Vector<Long>& dim_vec) {
  422. if (Dim(0) && Dim(1)) fftwf_destroy_plan(plan);
  423. this->fft_type = fft_type_;
  424. this->howmany = howmany_;
  425. plan = NULL;
  426. Long rank = dim_vec.Dim();
  427. Vector<int> dim_vec_(rank);
  428. for (Integer i = 0; i < rank; i++) {
  429. dim_vec_[i] = dim_vec[i];
  430. }
  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. } else {
  448. N0 = 0;
  449. N1 = 0;
  450. }
  451. this->dim[0] = N0;
  452. this->dim[1] = N1;
  453. }
  454. if (!N0 || !N1) return;
  455. Vector<ValueType> in (N0), out(N1);
  456. if (fft_type == FFT_Type::R2C) {
  457. 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);
  458. } else if (fft_type == FFT_Type::C2C) {
  459. 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);
  460. } else if (fft_type == FFT_Type::C2C_INV) {
  461. 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);
  462. } else if (fft_type == FFT_Type::C2R) {
  463. 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);
  464. }
  465. if (!plan) { // Build plan without FFTW_PRESERVE_INPUT
  466. if (fft_type == FFT_Type::R2C) {
  467. 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);
  468. } else if (fft_type == FFT_Type::C2C) {
  469. 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);
  470. } else if (fft_type == FFT_Type::C2C_INV) {
  471. 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);
  472. } else if (fft_type == FFT_Type::C2R) {
  473. 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);
  474. }
  475. tmp.Swap(in);
  476. } else {
  477. tmp.ReInit(0);
  478. }
  479. SCTL_ASSERT(plan);
  480. }
  481. void Execute(const Vector<ValueType>& in, Vector<ValueType>& out) const {
  482. Long N0 = this->Dim(0);
  483. Long N1 = this->Dim(1);
  484. if (!N0 || !N1) return;
  485. SCTL_ASSERT_MSG(in.Dim() == N0, "FFT: Wrong input size.");
  486. if (out.Dim() != N1) out.ReInit(N1);
  487. check_align(in, out);
  488. ValueType s = 0;
  489. auto in_ptr = in.begin();
  490. if (tmp.Dim()) { // Save input
  491. assert(tmp.Dim() == N0);
  492. in_ptr = tmp.begin();
  493. tmp = in;
  494. }
  495. if (fft_type == FFT_Type::R2C) {
  496. s = 1 / sqrt<ValueType>(N0 / howmany);
  497. fftwf_execute_dft_r2c(plan, (float*)&in_ptr[0], (fftwf_complex*)&out[0]);
  498. } else if (fft_type == FFT_Type::C2C) {
  499. s = 1 / sqrt<ValueType>(N0 / howmany * (ValueType)0.5);
  500. fftwf_execute_dft(plan, (fftwf_complex*)&in_ptr[0], (fftwf_complex*)&out[0]);
  501. } else if (fft_type == FFT_Type::C2C_INV) {
  502. s = 1 / sqrt<ValueType>(N1 / howmany * (ValueType)0.5);
  503. fftwf_execute_dft(plan, (fftwf_complex*)&in_ptr[0], (fftwf_complex*)&out[0]);
  504. } else if (fft_type == FFT_Type::C2R) {
  505. s = 1 / sqrt<ValueType>(N1 / howmany);
  506. fftwf_execute_dft_c2r(plan, (fftwf_complex*)&in_ptr[0], (float*)&out[0]);
  507. }
  508. for (auto& x : out) x *= s;
  509. }
  510. private:
  511. mutable Vector<ValueType> tmp;
  512. fftwf_plan plan;
  513. };
  514. #endif
  515. #ifdef SCTL_HAVE_FFTWL
  516. template <> class FFT<long double> : public FFT_Generic<long double, FFT<long double>> {
  517. typedef long double ValueType;
  518. public:
  519. ~FFT() { if (this->Dim(0) && this->Dim(1)) fftwl_destroy_plan(plan); }
  520. void Setup(FFT_Type fft_type_, Long howmany_, const Vector<Long>& dim_vec) {
  521. if (Dim(0) && Dim(1)) fftwl_destroy_plan(plan);
  522. this->fft_type = fft_type_;
  523. this->howmany = howmany_;
  524. plan = NULL;
  525. Long rank = dim_vec.Dim();
  526. Vector<int> dim_vec_(rank);
  527. for (Integer i = 0; i < rank; i++) dim_vec_[i] = dim_vec[i];
  528. Long N0, N1;
  529. { // Set N0, N1
  530. Long N = howmany;
  531. for (auto ni : dim_vec) N *= ni;
  532. if (fft_type == FFT_Type::R2C) {
  533. N0 = N;
  534. N1 = (N / dim_vec[rank - 1]) * (dim_vec[rank - 1] / 2 + 1) * 2;
  535. } else if (fft_type == FFT_Type::C2C) {
  536. N0 = N * 2;
  537. N1 = N * 2;
  538. } else if (fft_type == FFT_Type::C2C_INV) {
  539. N0 = N * 2;
  540. N1 = N * 2;
  541. } else if (fft_type == FFT_Type::C2R) {
  542. N0 = (N / dim_vec[rank - 1]) * (dim_vec[rank - 1] / 2 + 1) * 2;
  543. N1 = N;
  544. } else {
  545. N0 = 0;
  546. N1 = 0;
  547. }
  548. this->dim[0] = N0;
  549. this->dim[1] = N1;
  550. }
  551. if (!N0 || !N1) return;
  552. Vector<ValueType> in (N0), out(N1);
  553. if (fft_type == FFT_Type::R2C) {
  554. 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);
  555. } else if (fft_type == FFT_Type::C2C) {
  556. 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);
  557. } else if (fft_type == FFT_Type::C2C_INV) {
  558. 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);
  559. } else if (fft_type == FFT_Type::C2R) {
  560. 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);
  561. }
  562. if (!plan) { // Build plan without FFTW_PRESERVE_INPUT
  563. if (fft_type == FFT_Type::R2C) {
  564. 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);
  565. } else if (fft_type == FFT_Type::C2C) {
  566. 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);
  567. } else if (fft_type == FFT_Type::C2C_INV) {
  568. 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);
  569. } else if (fft_type == FFT_Type::C2R) {
  570. 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);
  571. }
  572. tmp.Swap(in);
  573. } else {
  574. tmp.ReInit(0);
  575. }
  576. SCTL_ASSERT(plan);
  577. }
  578. void Execute(const Vector<ValueType>& in, Vector<ValueType>& out) const {
  579. Long N0 = this->Dim(0);
  580. Long N1 = this->Dim(1);
  581. if (!N0 || !N1) return;
  582. SCTL_ASSERT_MSG(in.Dim() == N0, "FFT: Wrong input size.");
  583. if (out.Dim() != N1) out.ReInit(N1);
  584. check_align(in, out);
  585. ValueType s = 0;
  586. auto in_ptr = in.begin();
  587. if (tmp.Dim()) { // Save input
  588. assert(tmp.Dim() == N0);
  589. in_ptr = tmp.begin();
  590. tmp = in;
  591. }
  592. if (fft_type == FFT_Type::R2C) {
  593. s = 1 / sqrt<ValueType>(N0 / howmany);
  594. fftwl_execute_dft_r2c(plan, (long double*)&in_ptr[0], (fftwl_complex*)&out[0]);
  595. } else if (fft_type == FFT_Type::C2C) {
  596. s = 1 / sqrt<ValueType>(N0 / howmany * (ValueType)0.5);
  597. fftwl_execute_dft(plan, (fftwl_complex*)&in_ptr[0], (fftwl_complex*)&out[0]);
  598. } else if (fft_type == FFT_Type::C2C_INV) {
  599. s = 1 / sqrt<ValueType>(N1 / howmany * (ValueType)0.5);
  600. fftwl_execute_dft(plan, (fftwl_complex*)&in_ptr[0], (fftwl_complex*)&out[0]);
  601. } else if (fft_type == FFT_Type::C2R) {
  602. s = 1 / sqrt<ValueType>(N1 / howmany);
  603. fftwl_execute_dft_c2r(plan, (fftwl_complex*)&in_ptr[0], (long double*)&out[0]);
  604. }
  605. for (auto& x : out) x *= s;
  606. }
  607. private:
  608. mutable Vector<ValueType> tmp;
  609. fftwl_plan plan;
  610. };
  611. #endif
  612. } // end namespace
  613. #endif //_SCTL_FFT_WRAPPER_