isce3 0.25.0
Loading...
Searching...
No Matches
FFT.h
1#pragma once
2
3#include "FFTPlan.h"
4
5namespace isce3 { namespace cuda { namespace fft {
6
15template<typename T>
16FwdFFTPlan<T> planfft1d(thrust::complex<T> * out, thrust::complex<T> * in, int n);
17
19template<typename T>
20FwdFFTPlan<T> planfft1d(thrust::complex<T> * out, T * in, int n);
21
34template<typename T>
35FwdFFTPlan<T> planfft1d(thrust::complex<T> * out, thrust::complex<T> * in, const int (&dims)[2], int axis);
36
38template<typename T>
39FwdFFTPlan<T> planfft1d(thrust::complex<T> * out, T * in, const int (&dims)[2], int axis);
40
51template<typename T>
52FwdFFTPlan<T> planfft2d(thrust::complex<T> * out, thrust::complex<T> * in, const int (&dims)[2]);
53
55template<typename T>
56FwdFFTPlan<T> planfft2d(thrust::complex<T> * out, T * in, const int (&dims)[2]);
57
66template<typename T>
67InvFFTPlan<T> planifft1d(thrust::complex<T> * out, thrust::complex<T> * in, int n);
68
70template<typename T>
71InvFFTPlan<T> planifft1d(T * out, thrust::complex<T> * in, int n);
72
85template<typename T>
86InvFFTPlan<T> planifft1d(thrust::complex<T> * out, thrust::complex<T> * in, const int (&dims)[2], int axis);
87
89template<typename T>
90InvFFTPlan<T> planifft1d(T * out, thrust::complex<T> * in, const int (&dims)[2], int axis);
91
102template<typename T>
103InvFFTPlan<T> planifft2d(thrust::complex<T> * out, thrust::complex<T> * in, const int (&dims)[2]);
104
106template<typename T>
107InvFFTPlan<T> planifft2d(T * out, thrust::complex<T> * in, const int (&dims)[2]);
108
116template<typename T>
117void fft1d(thrust::complex<T> * out, thrust::complex<T> * in, int n);
118
120template<typename T>
121void fft1d(thrust::complex<T> * out, T * in, int n);
122
134template<typename T>
135void fft1d(thrust::complex<T> * out, thrust::complex<T> * in, const int (&dims)[2], int axis);
136
138template<typename T>
139void fft1d(thrust::complex<T> * out, T * in, const int (&dims)[2], int axis);
140
150template<typename T>
151void fft2d(thrust::complex<T> * out, thrust::complex<T> * in, const int (&dims)[2]);
152
154template<typename T>
155void fft2d(thrust::complex<T> * out, T * in, const int (&dims)[2]);
156
164template<typename T>
165void ifft1d(thrust::complex<T> * out, thrust::complex<T> * in, int n);
166
168template<typename T>
169void ifft1d(T * out, thrust::complex<T> * in, int n);
170
182template<typename T>
183void ifft1d(thrust::complex<T> * out, thrust::complex<T> * in, const int (&dims)[2], int axis);
184
186template<typename T>
187void ifft1d(T * out, thrust::complex<T> * in, const int (&dims)[2], int axis);
188
198template<typename T>
199void ifft2d(thrust::complex<T> * out, thrust::complex<T> * in, const int (&dims)[2]);
200
202template<typename T>
203void ifft2d(T * out, thrust::complex<T> * in, int (&dims)[2]);
204
205}}}
206
207#define ISCE_CUDA_FFT_FFT_ICC
208#include "FFT.icc"
209#undef ISCE_CUDA_FFT_FFT_ICC
RAII wrapper encapsulating cuFFT plan for forward FFT execution.
Definition FFTPlan.h:11
RAII wrapper encapsulating cuFFT plan for inverse FFT execution.
Definition FFTPlan.h:280
base interpolator is an abstract base class
Definition BinarySearchFunc.cpp:5

Generated for ISCE3.0 by doxygen 1.13.2.