16 #include <isce3/core/Constants.h>
35 std::valarray<std::complex<T>> &output,
36 int rank,
int* n,
int howmany,
37 int* inembed,
int istride,
int idist,
38 int* onembed,
int ostride,
int odist,
int sign);
44 std::complex<T>* output,
45 int rank,
int* n,
int howmany,
46 int* inembed,
int istride,
int idist,
47 int* onembed,
int ostride,
int odist,
int sign);
53 std::valarray<std::complex<T>> &output,
54 int rank,
int* n,
int howmany,
55 int* inembed,
int istride,
int idist,
56 int* onembed,
int ostride,
int odist);
62 std::complex<T>* output,
63 int rank,
int* n,
int howmany,
64 int* inembed,
int istride,
int idist,
65 int* onembed,
int ostride,
int odist);
71 std::valarray<std::complex<T>> &output,
72 int rank,
int* n,
int howmany,
73 int* inembed,
int istride,
int idist,
74 int* onembed,
int ostride,
int odist,
int sign);
80 std::complex<T>* output,
81 int rank,
int* n,
int howmany,
82 int* inembed,
int istride,
int idist,
83 int* onembed,
int ostride,
int odist,
int sign);
90 std::valarray<T>& output,
91 int rank,
int* n,
int howmany,
92 int* inembed,
int istride,
int idist,
93 int* onembed,
int ostride,
int odist);
101 int rank,
int* n,
int howmany,
102 int* inembed,
int istride,
int idist,
103 int* onembed,
int ostride,
int odist);
106 void forward(std::valarray<std::complex<T>> &input,
107 std::valarray<std::complex<T>> &output);
110 void forward(std::complex<T>* input,
111 std::complex<T>* output);
114 void forward(std::valarray<T> &input,
115 std::valarray<std::complex<T>> &output);
119 std::complex<T>* output);
124 void inverse(std::valarray<std::complex<T>> &input,
125 std::valarray<std::complex<T>> &output);
128 void inverse(std::complex<T>* input,
129 std::complex<T>* output);
132 void inverse(std::valarray<std::complex<T>>& input,
133 std::valarray<T>& output);
136 void inverse(std::complex<T>* input,
145 std::valarray<std::complex<T>>& spectrum,
146 int ncolumns,
int nrows);
154 std::complex<T>* spectrum,
155 int ncolumns,
int nrows);
163 std::valarray<std::complex<T>>& spectrum,
164 int ncolumns,
int nrows);
172 std::complex<T>* spectrum,
173 int ncolumns,
int nrows);
181 std::valarray<std::complex<T>> &spectrum,
182 int ncolumns,
int nrows);
190 std::complex<T>* spectrum,
191 int ncolumns,
int nrows);
199 std::valarray<std::complex<T>> &spectrum,
200 int ncolumns,
int nrows);
208 std::complex<T>* spectrum,
209 int ncolumns,
int nrows);
213 void forward2DFFT(std::valarray<std::complex<T>> &signal,
214 std::valarray<std::complex<T>> &spectrum,
215 int ncolumns,
int nrows);
219 void forward2DFFT(std::valarray<std::complex<T>> &signal,
220 std::valarray<std::complex<T>> &spectrum,
221 int incolumns,
int inrows,
222 int oncolumns,
int onrows);
227 std::complex<T>* spectrum,
228 int ncolumns,
int nrows);
233 std::complex<T>* spectrum,
234 int incolumns,
int inrows,
235 int oncolumns,
int onrows);
240 std::valarray<std::complex<T>>& spectrum,
241 int ncolumns,
int nrows);
246 std::valarray<std::complex<T>>& spectrum,
247 int incolumns,
int inrows,
248 int oncolumns,
int onrows);
253 std::complex<T>* spectrum,
254 int ncolumns,
int nrows);
259 std::complex<T>* spectrum,
260 int incolumns,
int inrows,
261 int oncolumns,
int onrows);
266 std::valarray<std::complex<T>> &signal,
267 int ncolumns,
int nrows);
272 std::complex<T>* signal,
273 int ncolumns,
int nrows);
278 std::valarray<T> &signal,
279 int ncolumns,
int nrows);
285 int ncolumns,
int nrows);
290 std::valarray<std::complex<T>> &signal,
291 int ncolumns,
int nrows);
296 std::complex<T>* signal,
297 int ncolumns,
int nrows);
301 std::valarray<T> &signal,
302 int ncolumns,
int nrows);
308 int ncolumns,
int nrows);
312 void inverse2DFFT(std::valarray<std::complex<T>> &spectrum,
313 std::valarray<std::complex<T>> &signal,
314 int ncolumns,
int nrows);
319 std::complex<T>* signal,
320 int ncolumns,
int nrows);
324 void inverse2DFFT(std::valarray<std::complex<T>> &spectrum,
325 std::valarray<T> &signal,
326 int ncolumns,
int nrows);
332 int ncolumns,
int nrows);
336 void upsample(std::valarray<std::complex<T>> &signal,
337 std::valarray<std::complex<T>> &signalOversampled,
338 int rows,
int fft_size,
int oversampleFactor);
344 void upsample(std::valarray<std::complex<T>> &signal,
345 std::valarray<std::complex<T>> &signalOversampled,
346 int rows,
int fft_size,
int oversampleFactor,
347 std::valarray<std::complex<T>> shiftImpact);
353 void upsample2D(std::valarray<std::complex<T>> &signal,
354 std::valarray<std::complex<T>> &signalOversampled,
355 int oversampleFactor);
361 void upsample2D(std::valarray<std::complex<T>> &signal,
362 std::valarray<std::complex<T>> &signalOversampled,
363 int oversampleFactor,
364 std::valarray<std::complex<T>> &shiftImpact);
369 std::complex<T> *signalOversampled,
370 int oversampleFactor);
377 std::complex<T> *signalOversampled,
378 int oversampleFactor,
379 std::complex<T> *shiftImpact);
390 int* inembed,
int istride,
int idist,
391 int* onembed,
int ostride,
int odist);
395 int* inembed,
int istride,
int idist,
396 int* onembed,
int ostride,
int odist);
440 std::unique_ptr<impl, void(*)(impl*)> pimpl;
443 #define ISCE_SIGNAL_SIGNAL_ICC
444 #include "Signal.icc"
445 #undef ISCE_SIGNAL_SIGNAL_ICC
void inverse2DFFT(std::valarray< std::complex< T >> &spectrum, std::valarray< std::complex< T >> &signal, int ncolumns, int nrows)
initiate plan for inverse two dimensional FFT for a block of data
Definition: Signal.cpp:871
Signal()
Default constructor.
Definition: Signal.cpp:20
A class to handle 2D FFT or 1D FFT in range or azimuth directions.
Definition: forward.h:10
void forwardAzimuthFFT(std::valarray< std::complex< T >> &signal, std::valarray< std::complex< T >> &spectrum, int ncolumns, int nrows)
initiate plan for forward FFT in azimuth direction for a block of complex data.
Definition: Signal.cpp:459
void _fwd_configure(int rank, int *n, int howmany, int *inembed, int istride, int idist, int *onembed, int ostride, int odist)
save FFT plan parameters
Definition: Signal.icc:221
void inverse(std::valarray< std::complex< T >> &input, std::valarray< std::complex< T >> &output)
inverse transform
Definition: Signal.cpp:341
void inverseAzimuthFFT(std::valarray< std::complex< T >> &spectrum, std::valarray< std::complex< T >> &signal, int ncolumns, int nrows)
initiate plan for inverse FFT in azimuth direction for a block of data
Definition: Signal.cpp:790
void _rev_configure(int rank, int *n, int howmany, int *inembed, int istride, int idist, int *onembed, int ostride, int odist)
save FFT plan parameters
Definition: Signal.icc:271
void inverseRangeFFT(std::valarray< std::complex< T >> &spectrum, std::valarray< std::complex< T >> &signal, int ncolumns, int nrows)
initiate plan for backward FFT in range direction for a block of data
Definition: Signal.cpp:733
void _fwd_configure2DFFT(int incolumns, int inrows, int oncolumns, int onrows)
determine the required parameters for setting 2D FFT plans
Definition: Signal.icc:150
void forward(std::valarray< std::complex< T >> &input, std::valarray< std::complex< T >> &output)
forward transform
Definition: Signal.cpp:288
void _rev_configure2DFFT(int ncolumns, int nrows)
determine the required parameters for setting 2D FFT plans
Definition: Signal.icc:182
void _fwd_configureRangeFFT(int ncolumns, int nrows)
determine the required parameters for setting range FFT plans
Definition: Signal.icc:35
void upsample(std::valarray< std::complex< T >> &signal, std::valarray< std::complex< T >> &signalOversampled, int rows, int fft_size, int oversampleFactor)
upsampling a block of data in range direction
Definition: Signal.cpp:952
void _rev_configureRangeFFT(int ncolumns, int nrows)
determine the required parameters for setting range FFT plans
Definition: Signal.icc:62
void fftPlanBackward(std::valarray< std::complex< T >> &input, std::valarray< std::complex< T >> &output, int rank, int *n, int howmany, int *inembed, int istride, int idist, int *onembed, int ostride, int odist, int sign)
initiate iverse FFTW3 plan for a block of spectrum input parameters follow FFTW3 interface for fftw_p...
Definition: Signal.cpp:171
void forwardRangeFFT(std::valarray< std::complex< T >> &signal, std::valarray< std::complex< T >> &spectrum, int ncolumns, int nrows)
initiate plan for forward FFT in range direction for a block of complex data.
Definition: Signal.cpp:381
void nextPowerOfTwo(size_t N, size_t &fftLength)
next power of two
Definition: Signal.icc:19
void forward2DFFT(std::valarray< std::complex< T >> &signal, std::valarray< std::complex< T >> &spectrum, int ncolumns, int nrows)
initiate plan for forward two imensional FFT for a block of complex data
Definition: Signal.cpp:535
void upsample2D(std::valarray< std::complex< T >> &signal, std::valarray< std::complex< T >> &signalOversampled, int oversampleFactor)
2D upsampling a block of 2D data
Definition: Signal.cpp:1051
void _fwd_configureAzimuthFFT(int ncolumns, int nrows)
determine the required parameters for setting azimuth FFT plans
Definition: Signal.icc:92
void _rev_configureAzimuthFFT(int ncolumns, int nrows)
determine the required parameters for setting azimuth FFT plans
Definition: Signal.icc:120
void fftPlanForward(std::valarray< std::complex< T >> &input, std::valarray< std::complex< T >> &output, int rank, int *n, int howmany, int *inembed, int istride, int idist, int *onembed, int ostride, int odist, int sign)
initiate forward FFTW3 plan for a block of complex data input parameters follow FFTW3 interface for f...
Definition: Signal.cpp:46