fft_wrapper.hpp 25 KB

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