isce3  0.1.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Pages
FFTPlan.h
1 #pragma once
2 
3 #include <thrust/complex.h>
4 
5 #include "detail/FFTPlanBase.h"
6 
7 namespace isce3 { namespace cuda { namespace fft {
8 
10 template<typename T>
11 class FwdFFTPlan final : public detail::FFTPlanBase<CUFFT_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 
35  thrust::complex<T> * in,
36  int n,
37  int batch = 1);
38 
48  template<int Rank>
50  thrust::complex<T> * in,
51  const int (&n)[Rank],
52  int batch = 1);
53 
66  thrust::complex<T> * in,
67  int n,
68  int nembed,
69  int stride,
70  int dist,
71  int batch = 1);
72 
85  template<int Rank>
87  thrust::complex<T> * in,
88  const int (&n)[Rank],
89  const int (&nembed)[Rank],
90  int stride,
91  int dist,
92  int batch = 1);
93 
109  thrust::complex<T> * in,
110  int n,
111  int inembed,
112  int istride,
113  int idist,
114  int onembed,
115  int ostride,
116  int odist,
117  int batch = 1);
118 
134  template<int Rank>
136  thrust::complex<T> * in,
137  const int (&n)[Rank],
138  const int (&inembed)[Rank],
139  int istride,
140  int idist,
141  const int (&onembed)[Rank],
142  int ostride,
143  int odist,
144  int batch = 1);
145 #endif
146 
156  T * in,
157  int n,
158  int batch = 1);
159 
169  template<int Rank>
171  T * in,
172  const int (&n)[Rank],
173  int batch = 1);
174 
187  T * in,
188  int n,
189  int nembed,
190  int stride,
191  int dist,
192  int batch = 1);
193 
206  template<int Rank>
208  T * in,
209  const int (&n)[Rank],
210  const int (&nembed)[Rank],
211  int stride,
212  int dist,
213  int batch = 1);
214 
230  T * in,
231  int n,
232  int inembed,
233  int istride,
234  int idist,
235  int onembed,
236  int ostride,
237  int odist,
238  int batch = 1);
239 
255  template<int Rank>
257  T * in,
258  const int (&n)[Rank],
259  const int (&inembed)[Rank],
260  int istride,
261  int idist,
262  const int (&onembed)[Rank],
263  int ostride,
264  int odist,
265  int batch = 1);
266 
267  // the following are declared here just for the purpose of doxygenating them
268  // they are actually defined in the base class
269 #if 0
270 
271  explicit operator bool() const;
272 
274  void execute() const;
275 #endif
276 };
277 
279 template<typename T>
280 class InvFFTPlan final : public detail::FFTPlanBase<CUFFT_INVERSE, T> {
281 public:
283  using super_t::super_t;
284 
291 
292  // the following are declared here just for the purpose of doxygenating them
293  // they are actually defined in the base class
294 #if 0
295 
304  thrust::complex<T> * in,
305  int n,
306  int batch = 1);
307 
317  template<int Rank>
319  thrust::complex<T> * in,
320  const int (&n)[Rank],
321  int batch = 1);
322 
335  thrust::complex<T> * in,
336  int n,
337  int nembed,
338  int stride,
339  int dist,
340  int batch = 1);
341 
354  template<int Rank>
356  thrust::complex<T> * in,
357  const int (&n)[Rank],
358  const int (&nembed)[Rank],
359  int stride,
360  int dist,
361  int batch = 1);
362 
378  thrust::complex<T> * in,
379  int n,
380  int inembed,
381  int istride,
382  int idist,
383  int onembed,
384  int ostride,
385  int odist,
386  int batch = 1);
387 
403  template<int Rank>
405  thrust::complex<T> * in,
406  const int (&n)[Rank],
407  const int (&inembed)[Rank],
408  int istride,
409  int idist,
410  const int (&onembed)[Rank],
411  int ostride,
412  int odist,
413  int batch = 1);
414 #endif
415 
424  InvFFTPlan(T * out,
425  thrust::complex<T> * in,
426  int n,
427  int batch = 1);
428 
438  template<int Rank>
439  InvFFTPlan(T * out,
440  thrust::complex<T> * in,
441  const int (&n)[Rank],
442  int batch = 1);
443 
455  InvFFTPlan(T * out,
456  thrust::complex<T> * in,
457  int n,
458  int nembed,
459  int stride,
460  int dist,
461  int batch = 1);
462 
475  template<int Rank>
476  InvFFTPlan(T * out,
477  thrust::complex<T> * in,
478  const int (&n)[Rank],
479  const int (&nembed)[Rank],
480  int stride,
481  int dist,
482  int batch = 1);
483 
498  InvFFTPlan(T * out,
499  thrust::complex<T> * in,
500  int n,
501  int inembed,
502  int istride,
503  int idist,
504  int onembed,
505  int ostride,
506  int odist,
507  int batch = 1);
508 
524  template<int Rank>
525  InvFFTPlan(T * out,
526  thrust::complex<T> * in,
527  const int (&n)[Rank],
528  const int (&inembed)[Rank],
529  int istride,
530  int idist,
531  const int (&onembed)[Rank],
532  int ostride,
533  int odist,
534  int batch = 1);
535 
536  // the following are declared here just for the purpose of doxygenating them
537  // they are actually defined in the base class
538 #if 0
539 
540  explicit operator bool() const;
541 
543  void execute() const;
544 #endif
545 };
546 
547 }}}
548 
549 #define ISCE_CUDA_FFT_FFTPLAN_ICC
550 #include "FFTPlan.icc"
551 #undef ISCE_CUDA_FFT_FFTPLAN_ICC
InvFFTPlan()
Construct an invalid plan.
Definition: FFTPlan.h:290
FwdFFTPlan()
Construct an invalid plan.
Definition: FFTPlan.h:21
RAII wrapper encapsulating cuFFT plan for inverse FFT execution.
Definition: FFTPlan.h:280
RAII wrapper encapsulating cuFFT plan for forward FFT execution.
Definition: FFTPlan.h:11
Definition: FFTPlanBase.h:12

Generated for ISCE3.0 by doxygen 1.8.5.