isce3  0.1.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Pages
FFTPlan.h
1 #pragma once
2 
3 #include <complex>
4 
5 #include "detail/FFTPlanBase.h"
6 
7 namespace isce3 { namespace fft {
8 
10 template<typename T>
11 class FwdFFTPlan final : public detail::FFTPlanBase<FFTW_FORWARD, T> {
12 public:
14  using super_t::super_t;
15 
22 
23  // the following are declared here just for the purpose of doxygenating them
24  // they are actually defined in the base class
25 #if 0
26 
36  FwdFFTPlan(std::complex<T> * out,
37  std::complex<T> * in,
38  int n,
39  int batch = 1,
40  unsigned flags = FFTW_MEASURE,
41  int threads = detail::getMaxThreads());
42 
54  template<int Rank>
55  FwdFFTPlan(std::complex<T> * out,
56  std::complex<T> * in,
57  const int (&n)[Rank],
58  int batch = 1,
59  unsigned flags = FFTW_MEASURE,
60  int threads = detail::getMaxThreads());
61 
75  FwdFFTPlan(std::complex<T> * out,
76  std::complex<T> * in,
77  int n,
78  int nembed,
79  int stride,
80  int dist,
81  int batch = 1,
82  unsigned flags = FFTW_MEASURE,
83  int threads = detail::getMaxThreads());
84 
99  template<int Rank>
100  FwdFFTPlan(std::complex<T> * out,
101  std::complex<T> * in,
102  const int (&n)[Rank],
103  const int (&nembed)[Rank],
104  int stride,
105  int dist,
106  int batch = 1,
107  unsigned flags = FFTW_MEASURE,
108  int threads = detail::getMaxThreads());
109 
126  FwdFFTPlan(std::complex<T> * out,
127  std::complex<T> * in,
128  int n,
129  int inembed,
130  int istride,
131  int idist,
132  int onembed,
133  int ostride,
134  int odist,
135  int batch = 1,
136  unsigned flags = FFTW_MEASURE,
137  int threads = detail::getMaxThreads());
138 
156  template<int Rank>
157  FwdFFTPlan(std::complex<T> * out,
158  std::complex<T> * in,
159  const int (&n)[Rank],
160  const int (&inembed)[Rank],
161  int istride,
162  int idist,
163  const int (&onembed)[Rank],
164  int ostride,
165  int odist,
166  int batch = 1,
167  unsigned flags = FFTW_MEASURE,
168  int threads = detail::getMaxThreads());
169 #endif
170 
181  FwdFFTPlan(std::complex<T> * out,
182  T * in,
183  int n,
184  int batch = 1,
185  unsigned flags = FFTW_MEASURE,
186  int threads = detail::getMaxThreads());
187 
199  template<int Rank>
200  FwdFFTPlan(std::complex<T> * out,
201  T * in,
202  const int (&n)[Rank],
203  int batch = 1,
204  unsigned flags = FFTW_MEASURE,
205  int threads = detail::getMaxThreads());
206 
220  FwdFFTPlan(std::complex<T> * out,
221  T * in,
222  int n,
223  int nembed,
224  int stride,
225  int dist,
226  int batch = 1,
227  unsigned flags = FFTW_MEASURE,
228  int threads = detail::getMaxThreads());
229 
244  template<int Rank>
245  FwdFFTPlan(std::complex<T> * out,
246  T * in,
247  const int (&n)[Rank],
248  const int (&nembed)[Rank],
249  int stride,
250  int dist,
251  int batch = 1,
252  unsigned flags = FFTW_MEASURE,
253  int threads = detail::getMaxThreads());
254 
271  FwdFFTPlan(std::complex<T> * out,
272  T * in,
273  int n,
274  int inembed,
275  int istride,
276  int idist,
277  int onembed,
278  int ostride,
279  int odist,
280  int batch = 1,
281  unsigned flags = FFTW_MEASURE,
282  int threads = detail::getMaxThreads());
283 
301  template<int Rank>
302  FwdFFTPlan(std::complex<T> * out,
303  T * in,
304  const int (&n)[Rank],
305  const int (&inembed)[Rank],
306  int istride,
307  int idist,
308  const int (&onembed)[Rank],
309  int ostride,
310  int odist,
311  int batch = 1,
312  unsigned flags = FFTW_MEASURE,
313  int threads = detail::getMaxThreads());
314 
315  // the following are declared here just for the purpose of doxygenating them
316  // they are actually defined in the base class
317 #if 0
318 
319  explicit operator bool() const;
320 
322  void execute() const;
323 #endif
324 };
325 
327 template<typename T>
328 class InvFFTPlan final : public detail::FFTPlanBase<FFTW_BACKWARD, T> {
329 public:
331  using super_t::super_t;
332 
339 
340  // the following are declared here just for the purpose of doxygenating them
341  // they are actually defined in the base class
342 #if 0
343 
353  InvFFTPlan(std::complex<T> * out,
354  std::complex<T> * in,
355  int n,
356  int batch = 1,
357  unsigned flags = FFTW_MEASURE,
358  int threads = detail::getMaxThreads());
359 
371  template<int Rank>
372  InvFFTPlan(std::complex<T> * out,
373  std::complex<T> * in,
374  const int (&n)[Rank],
375  int batch = 1,
376  unsigned flags = FFTW_MEASURE,
377  int threads = detail::getMaxThreads());
378 
392  InvFFTPlan(std::complex<T> * out,
393  std::complex<T> * in,
394  int n,
395  int nembed,
396  int stride,
397  int dist,
398  int batch = 1,
399  unsigned flags = FFTW_MEASURE,
400  int threads = detail::getMaxThreads());
401 
416  template<int Rank>
417  InvFFTPlan(std::complex<T> * out,
418  std::complex<T> * in,
419  const int (&n)[Rank],
420  const int (&nembed)[Rank],
421  int stride,
422  int dist,
423  int batch = 1,
424  unsigned flags = FFTW_MEASURE,
425  int threads = detail::getMaxThreads());
426 
443  InvFFTPlan(std::complex<T> * out,
444  std::complex<T> * in,
445  int n,
446  int inembed,
447  int istride,
448  int idist,
449  int onembed,
450  int ostride,
451  int odist,
452  int batch = 1,
453  unsigned flags = FFTW_MEASURE,
454  int threads = detail::getMaxThreads());
455 
473  template<int Rank>
474  InvFFTPlan(std::complex<T> * out,
475  std::complex<T> * in,
476  const int (&n)[Rank],
477  const int (&inembed)[Rank],
478  int istride,
479  int idist,
480  const int (&onembed)[Rank],
481  int ostride,
482  int odist,
483  int batch = 1,
484  unsigned flags = FFTW_MEASURE,
485  int threads = detail::getMaxThreads());
486 #endif
487 
498  InvFFTPlan(T * out,
499  std::complex<T> * in,
500  int n,
501  int batch = 1,
502  unsigned flags = FFTW_MEASURE,
503  int threads = detail::getMaxThreads());
504 
516  template<int Rank>
517  InvFFTPlan(T * out,
518  std::complex<T> * in,
519  const int (&n)[Rank],
520  int batch = 1,
521  unsigned flags = FFTW_MEASURE,
522  int threads = detail::getMaxThreads());
523 
537  InvFFTPlan(T * out,
538  std::complex<T> * in,
539  int n,
540  int nembed,
541  int stride,
542  int dist,
543  int batch = 1,
544  unsigned flags = FFTW_MEASURE,
545  int threads = detail::getMaxThreads());
546 
561  template<int Rank>
562  InvFFTPlan(T * out,
563  std::complex<T> * in,
564  const int (&n)[Rank],
565  const int (&nembed)[Rank],
566  int stride,
567  int dist,
568  int batch = 1,
569  unsigned flags = FFTW_MEASURE,
570  int threads = detail::getMaxThreads());
571 
588  InvFFTPlan(T * out,
589  std::complex<T> * in,
590  int n,
591  int inembed,
592  int istride,
593  int idist,
594  int onembed,
595  int ostride,
596  int odist,
597  int batch = 1,
598  unsigned flags = FFTW_MEASURE,
599  int threads = detail::getMaxThreads());
600 
618  template<int Rank>
619  InvFFTPlan(T * out,
620  std::complex<T> * in,
621  const int (&n)[Rank],
622  const int (&inembed)[Rank],
623  int istride,
624  int idist,
625  const int (&onembed)[Rank],
626  int ostride,
627  int odist,
628  int batch = 1,
629  unsigned flags = FFTW_MEASURE,
630  int threads = detail::getMaxThreads());
631 
632  // the following are declared here just for the purpose of doxygenating them
633  // they are actually defined in the base class
634 #if 0
635 
636  explicit operator bool() const;
637 
639  void execute() const;
640 #endif
641 };
642 
643 }}
644 
645 #define ISCE_FFT_FFTPLAN_ICC
646 #include "FFTPlan.icc"
647 #undef ISCE_FFT_FFTPLAN_ICC
RAII wrapper encapsulating FFTW plan for forward FFT execution.
Definition: FFTPlan.h:11
FwdFFTPlan()
Construct an invalid plan.
Definition: FFTPlan.h:21
Definition: FFTPlanBase.h:13
RAII wrapper encapsulating FFTW plan for inverse FFT execution.
Definition: FFTPlan.h:328
InvFFTPlan()
Construct an invalid plan.
Definition: FFTPlan.h:338

Generated for ISCE3.0 by doxygen 1.8.5.