isce3 0.25.0
|
A class to handle 2D FFT or 1D FFT in range or azimuth directions. More...
#include <Signal.h>
Public Member Functions | |
Signal () | |
Default constructor. | |
Signal (int nthreads) | |
Constructor with number of threads. | |
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 fftw_plan_many_dft | |
void | fftPlanForward (std::complex< T > *input, 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 fftw_plan_many_dft | |
void | fftPlanForward (std::valarray< 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) |
initiate forward FFTW3 plan for a block of real data input parameters follow FFTW3 interface for fftw_plan_many_dft_r2c | |
void | fftPlanForward (T *input, std::complex< T > *output, int rank, int *n, int howmany, int *inembed, int istride, int idist, int *onembed, int ostride, int odist) |
initiate forward FFTW3 plan for a block of real data input parameters follow FFTW3 interface for fftw_plan_many_dft_r2c | |
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_plan_many_dft | |
void | fftPlanBackward (std::complex< T > *input, 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_plan_many_dft | |
void | fftPlanBackward (std::valarray< std::complex< T > > &input, std::valarray< T > &output, int rank, int *n, int howmany, int *inembed, int istride, int idist, int *onembed, int ostride, int odist) |
initiate iverse FFTW3 plan for a block of spectrum to be transformed to real data input parameters follow FFTW3 interface for fftw_plan_many_dft_c2r | |
void | fftPlanBackward (std::complex< T > *input, T *output, int rank, int *n, int howmany, int *inembed, int istride, int idist, int *onembed, int ostride, int odist) |
initiate iverse FFTW3 plan for a block of spectrum to be transformed to real data input parameters follow FFTW3 interface for fftw_plan_many_dft_c2r | |
void | forward (std::valarray< std::complex< T > > &input, std::valarray< std::complex< T > > &output) |
forward transform | |
void | forward (std::complex< T > *input, std::complex< T > *output) |
forward transform | |
void | forward (std::valarray< T > &input, std::valarray< std::complex< T > > &output) |
forward transform | |
void | forward (T *input, std::complex< T > *output) |
forward transform | |
void | inverse (std::valarray< std::complex< T > > &input, std::valarray< std::complex< T > > &output) |
inverse transform | |
void | inverse (std::complex< T > *input, std::complex< T > *output) |
inverse transform | |
void | inverse (std::valarray< std::complex< T > > &input, std::valarray< T > &output) |
inverse transform | |
void | inverse (std::complex< T > *input, T *output) |
inverse transform | |
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. | |
void | forwardRangeFFT (std::complex< T > *signal, std::complex< T > *spectrum, int ncolumns, int nrows) |
initiate plan for forward FFT in range direction for a block of complex data. | |
void | forwardRangeFFT (std::valarray< T > &signal, std::valarray< std::complex< T > > &spectrum, int ncolumns, int nrows) |
initiate plan for forward FFT in range direction for a block of real data. | |
void | forwardRangeFFT (T *signal, std::complex< T > *spectrum, int ncolumns, int nrows) |
initiate plan for forward FFT in range direction for a block of real data. | |
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. | |
void | forwardAzimuthFFT (std::complex< T > *signal, std::complex< T > *spectrum, int ncolumns, int nrows) |
initiate plan for forward FFT in azimuth direction for a block of complex data. | |
void | forwardAzimuthFFT (std::valarray< T > &signal, std::valarray< std::complex< T > > &spectrum, int ncolumns, int nrows) |
initiate plan for forward FFT in azimuth direction for a block of real data. | |
void | forwardAzimuthFFT (T *signal, std::complex< T > *spectrum, int ncolumns, int nrows) |
initiate plan for forward FFT in azimuth direction for a block of real data. | |
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 | |
void | forward2DFFT (std::valarray< std::complex< T > > &signal, std::valarray< std::complex< T > > &spectrum, int incolumns, int inrows, int oncolumns, int onrows) |
initiate plan for forward two imensional FFT for a block of complex data | |
void | forward2DFFT (std::complex< T > *signal, std::complex< T > *spectrum, int ncolumns, int nrows) |
initiate plan for forward two imensional FFT for a block of complex data | |
void | forward2DFFT (std::complex< T > *signal, std::complex< T > *spectrum, int incolumns, int inrows, int oncolumns, int onrows) |
initiate plan for forward two imensional FFT for a block of complex data | |
void | forward2DFFT (std::valarray< T > &signal, std::valarray< std::complex< T > > &spectrum, int ncolumns, int nrows) |
initiate plan for forward two imensional FFT for a block of real data | |
void | forward2DFFT (std::valarray< T > &signal, std::valarray< std::complex< T > > &spectrum, int incolumns, int inrows, int oncolumns, int onrows) |
initiate plan for forward two imensional FFT for a block of real data | |
void | forward2DFFT (T *signal, std::complex< T > *spectrum, int ncolumns, int nrows) |
initiate plan for forward two imensional FFT for a block of real data | |
void | forward2DFFT (T *signal, std::complex< T > *spectrum, int incolumns, int inrows, int oncolumns, int onrows) |
initiate plan for forward two imensional FFT for a block of real data | |
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 | |
void | inverseRangeFFT (std::complex< T > *spectrum, std::complex< T > *signal, int ncolumns, int nrows) |
initiate plan for backward FFT in range direction for a block of data | |
void | inverseRangeFFT (std::valarray< std::complex< T > > &spectrum, std::valarray< T > &signal, int ncolumns, int nrows) |
initiate plan for backward FFT in range direction for a block of data | |
void | inverseRangeFFT (std::complex< T > *spectrum, T *signal, int ncolumns, int nrows) |
initiate plan for backward FFT in range direction for a block of data | |
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 | |
void | inverseAzimuthFFT (std::complex< T > *spectrum, std::complex< T > *signal, int ncolumns, int nrows) |
initiate plan for inverse FFT in azimuth direction for a block of data | |
void | inverseAzimuthFFT (std::valarray< std::complex< T > > &spectrum, std::valarray< T > &signal, int ncolumns, int nrows) |
initiate plan for inverse FFT in azimuth direction for a block of data | |
void | inverseAzimuthFFT (std::complex< T > *spectrum, T *signal, int ncolumns, int nrows) |
initiate plan for inverse FFT in azimuth direction for a block of data | |
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 | |
void | inverse2DFFT (std::complex< T > *spectrum, std::complex< T > *signal, int ncolumns, int nrows) |
initiate plan for inverse two dimensional FFT for a block of data | |
void | inverse2DFFT (std::valarray< std::complex< T > > &spectrum, std::valarray< T > &signal, int ncolumns, int nrows) |
initiate plan for inverse two dimensional FFT for a block of data | |
void | inverse2DFFT (std::complex< T > *spectrum, T *signal, int ncolumns, int nrows) |
initiate plan for inverse two dimensional FFT for a block of data | |
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 | |
void | upsample (std::valarray< std::complex< T > > &signal, std::valarray< std::complex< T > > &signalOversampled, int rows, int fft_size, int oversampleFactor, std::valarray< std::complex< T > > shiftImpact) |
upsampling a block of data in range direction and shifting the upsampled signal by a constant. | |
void | upsample (isce3::core::EArray2D< std::complex< T > > &signal, isce3::core::EArray2D< std::complex< T > > &signalUpsampled, const isce3::core::EArray2D< std::complex< T > > &shiftImpact={}) |
upsampling a basebanded block of data in range (columns) direction and shifting the upsampled signal by a constant. | |
void | upsample2D (std::valarray< std::complex< T > > &signal, std::valarray< std::complex< T > > &signalOversampled, int oversampleFactor) |
2D upsampling a block of 2D data | |
void | upsample2D (std::valarray< std::complex< T > > &signal, std::valarray< std::complex< T > > &signalOversampled, int oversampleFactor, std::valarray< std::complex< T > > &shiftImpact) |
2D upsampling a block of 2D data and shifting the upsampled signal by a constant. | |
void | upsample2D (std::complex< T > *signal, std::complex< T > *signalOversampled, int oversampleFactor) |
2D upsampling a block of 2D data | |
void | upsample2D (std::complex< T > *signal, std::complex< T > *signalOversampled, int oversampleFactor, std::complex< T > *shiftImpact) |
2D upsampling a block of 2D data and shifting the upsampled signal by a constant. | |
void | nextPowerOfTwo (size_t N, size_t &fftLength) |
next power of two | |
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 | |
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 | |
void | _fwd_configureRangeFFT (int ncolumns, int nrows) |
determine the required parameters for setting range FFT plans | |
void | _fwd_configureAzimuthFFT (int ncolumns, int nrows) |
determine the required parameters for setting azimuth FFT plans | |
void | _fwd_configure2DFFT (int incolumns, int inrows, int oncolumns, int onrows) |
determine the required parameters for setting 2D FFT plans | |
void | _rev_configureRangeFFT (int ncolumns, int nrows) |
determine the required parameters for setting range FFT plans | |
void | _rev_configureAzimuthFFT (int ncolumns, int nrows) |
determine the required parameters for setting azimuth FFT plans | |
void | _rev_configure2DFFT (int ncolumns, int nrows) |
determine the required parameters for setting 2D FFT plans | |
A class to handle 2D FFT or 1D FFT in range or azimuth directions.
isce3::signal::Signal< T >::Signal | ( | int | nthreads | ) |
Constructor with number of threads.
This uses the Multi-threaded FFTW
|
inline |
save FFT plan parameters
[in] | rank | dataset number of dimensions |
[in] | n | size of each dimension (array) |
[in] | howmany | number of forward FFT |
[in] | inembed | layout of the input container |
[in] | istride | packing of the input data |
[in] | idist | number of elements between two input dataset |
[in] | onembed | layout of the output container |
[in] | ostride | packing of the output data |
[in] | odist | number of elements between two output dataset |
|
inline |
determine the required parameters for setting 2D FFT plans
[in] | incolumns | number of columns |
[in] | inrows | number of rows |
[in] | oncolumns | number of columns in output container |
[in] | onrows | number of rows in output container |
|
inline |
determine the required parameters for setting azimuth FFT plans
[in] | ncolumns | number of columns |
[in] | nrows | number of rows |
|
inline |
determine the required parameters for setting range FFT plans
[in] | ncolumns | number of columns |
[in] | nrows | number of rows |
|
inline |
save FFT plan parameters
[in] | rank | dataset number of dimensions |
[in] | n | size of each dimension (array) |
[in] | howmany | number of forward FFT |
[in] | inembed | layout of the input container |
[in] | istride | packing of the input data |
[in] | idist | number of elements between two input dataset |
[in] | onembed | layout of the output container |
[in] | ostride | packing of the output data |
[in] | odist | number of elements between two output dataset |
|
inline |
determine the required parameters for setting 2D FFT plans
[in] | ncolumns | number of columns |
[in] | nrows | number of rows |
|
inline |
determine the required parameters for setting azimuth FFT plans
[in] | ncolumns | number of columns |
[in] | nrows | number of rows |
|
inline |
determine the required parameters for setting range FFT plans
[in] | ncolumns | number of columns |
[in] | nrows | number of rows |
void isce3::signal::Signal< T >::fftPlanBackward | ( | std::complex< T > * | input, |
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_plan_many_dft
[in] | input | block of data |
[out] | output | block of data |
[in] | rank | |
[in] | size | |
[in] | howmany | |
[in] | inembed | |
[in] | istride | |
[in] | idist | |
[in] | onembed | |
[in] | ostride | |
[in] | odist | |
[in] | sign |
void isce3::signal::Signal< T >::fftPlanBackward | ( | std::complex< T > * | input, |
T * | output, | ||
int | rank, | ||
int * | n, | ||
int | howmany, | ||
int * | inembed, | ||
int | istride, | ||
int | idist, | ||
int * | onembed, | ||
int | ostride, | ||
int | odist ) |
initiate iverse FFTW3 plan for a block of spectrum to be transformed to real data input parameters follow FFTW3 interface for fftw_plan_many_dft_c2r
[in] | input | block of data |
[out] | output | block of data |
[in] | rank | |
[in] | size | |
[in] | howmany | |
[in] | inembed | |
[in] | istride | |
[in] | idist | |
[in] | onembed | |
[in] | ostride | |
[in] | odist |
void isce3::signal::Signal< T >::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_plan_many_dft
[in] | input | block of data |
[out] | output | block of data |
[in] | rank | |
[in] | size | |
[in] | howmany | |
[in] | inembed | |
[in] | istride | |
[in] | idist | |
[in] | onembed | |
[in] | ostride | |
[in] | odist | |
[in] | sign |
void isce3::signal::Signal< T >::fftPlanBackward | ( | std::valarray< std::complex< T > > & | input, |
std::valarray< T > & | output, | ||
int | rank, | ||
int * | n, | ||
int | howmany, | ||
int * | inembed, | ||
int | istride, | ||
int | idist, | ||
int * | onembed, | ||
int | ostride, | ||
int | odist ) |
initiate iverse FFTW3 plan for a block of spectrum to be transformed to real data input parameters follow FFTW3 interface for fftw_plan_many_dft_c2r
[in] | input | block of data |
[out] | output | block of data |
[in] | rank | |
[in] | size | |
[in] | howmany | |
[in] | inembed | |
[in] | istride | |
[in] | idist | |
[in] | onembed | |
[in] | ostride | |
[in] | odist |
void isce3::signal::Signal< T >::fftPlanForward | ( | std::complex< T > * | input, |
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 fftw_plan_many_dft
[in] | input | block of data |
[out] | output | block of data |
[in] | rank | rank of the transform (1: for one dimensional and 2: for two dimensional transform) |
[in] | size | size of each transform (ncols: for range FFT, nrows: for azimuth FFT) |
[in] | howmany | number of FFT transforms for a block of data (nrows: for range FFT, ncols: for azimuth FFT) |
[in] | inembed | |
[in] | istride | |
[in] | idist | |
[in] | onembed | |
[in] | ostride | |
[in] | odist | |
[in] | sign |
void isce3::signal::Signal< T >::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 fftw_plan_many_dft
[in] | input | block of data |
[out] | output | block of data |
[in] | rank | rank of the transform (1: for one dimensional and 2: for two dimensional transform) |
[in] | size | size of each transform (ncols: for range FFT, nrows: for azimuth FFT) |
[in] | howmany | number of FFT transforms for a block of data (nrows: for range FFT, ncols: for azimuth FFT) |
[in] | inembed | |
[in] | istride | |
[in] | idist | |
[in] | onembed | |
[in] | ostride | |
[in] | odist | |
[in] | sign |
void isce3::signal::Signal< T >::fftPlanForward | ( | std::valarray< 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 ) |
initiate forward FFTW3 plan for a block of real data input parameters follow FFTW3 interface for fftw_plan_many_dft_r2c
[in] | input | block of data |
[out] | output | block of data |
[in] | rank | rank of the transform (1: for one dimensional and 2: for two dimensional transform) |
[in] | size | size of each transform (ncols: for range FFT, nrows: for azimuth FFT) |
[in] | howmany | number of FFT transforms for a block of data (nrows: for range FFT, ncols: for azimuth FFT) |
[in] | inembed | |
[in] | istride | |
[in] | idist | |
[in] | onembed | |
[in] | ostride | |
[in] | odist |
void isce3::signal::Signal< T >::fftPlanForward | ( | T * | input, |
std::complex< T > * | output, | ||
int | rank, | ||
int * | n, | ||
int | howmany, | ||
int * | inembed, | ||
int | istride, | ||
int | idist, | ||
int * | onembed, | ||
int | ostride, | ||
int | odist ) |
initiate forward FFTW3 plan for a block of real data input parameters follow FFTW3 interface for fftw_plan_many_dft_r2c
[in] | input | block of data |
[out] | output | block of data |
[in] | rank | rank of the transform (1: for one dimensional and 2: for two dimensional transform) |
[in] | size | size of each transform (ncols: for range FFT, nrows: for azimuth FFT) |
[in] | howmany | number of FFT transforms for a block of data (nrows: for range FFT, ncols: for azimuth FFT) |
[in] | inembed | |
[in] | istride | |
[in] | idist | |
[in] | onembed | |
[in] | ostride | |
[in] | odist |
void isce3::signal::Signal< T >::forward | ( | std::complex< T > * | input, |
std::complex< T > * | output ) |
forward transform
unnormalized forward transform
[in] | input | block of data |
[out] | output | block of spectrum |
void isce3::signal::Signal< T >::forward | ( | std::valarray< std::complex< T > > & | input, |
std::valarray< std::complex< T > > & | output ) |
forward transform
unnormalized forward transform
[in] | input | block of data |
[out] | output | block of spectrum |
void isce3::signal::Signal< T >::forward | ( | std::valarray< T > & | input, |
std::valarray< std::complex< T > > & | output ) |
forward transform
unnormalized forward transform
[in] | input | block of data |
[out] | output | block of spectrum |
void isce3::signal::Signal< T >::forward | ( | T * | input, |
std::complex< T > * | output ) |
forward transform
unnormalized forward transform
[in] | input | block of data |
[out] | output | block of spectrum |
void isce3::signal::Signal< T >::forward2DFFT | ( | std::complex< T > * | signal, |
std::complex< T > * | spectrum, | ||
int | incolumns, | ||
int | inrows, | ||
int | oncolumns, | ||
int | onrows ) |
initiate plan for forward two imensional FFT for a block of complex data
[in] | signal | input block of data |
[out] | spectrum | output block of spectrum |
[in] | incolumns | number of columns of the block of data in input container |
[in] | inrows | number of rows of the block of data in input container |
[in] | oncolumns | number of columns of the block of data in output container |
[in] | oupnrows | number of rows of the block of data in output container |
If doing an out-of-place FFT, the output can be stored in a container whose number of columns and rows could be different than the input. This would be necessary for instance in the case of 2D upsampling
void isce3::signal::Signal< T >::forward2DFFT | ( | std::complex< T > * | signal, |
std::complex< T > * | spectrum, | ||
int | ncolumns, | ||
int | nrows ) |
initiate plan for forward two imensional FFT for a block of complex data
[in] | signal | input block of data |
[out] | spectrum | output block of spectrum |
[in] | ncolumns | number of columns of the block of data |
[in] | nrows | number of rows of the block of data |
void isce3::signal::Signal< T >::forward2DFFT | ( | std::valarray< std::complex< T > > & | signal, |
std::valarray< std::complex< T > > & | spectrum, | ||
int | incolumns, | ||
int | inrows, | ||
int | oncolumns, | ||
int | onrows ) |
initiate plan for forward two imensional FFT for a block of complex data
[in] | signal | input block of data |
[out] | spectrum | output block of spectrum |
[in] | incolumns | number of columns of the block of data |
[in] | inrows | number of rows of the block of data |
[in] | oncolumns | number of columns of the block of data in output container |
[in] | onrows | number of rows of the block of data in output container |
If doing an out-of-place FFT, the output can be stored in a container whose number of columns and rows could be different than the input. This would be necessary for instance in the case of 2D upsampling
void isce3::signal::Signal< T >::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
[in] | signal | input block of data |
[out] | spectrum | output block of spectrum |
[in] | ncolumns | number of columns of the block of data |
[in] | nrows | number of rows of the block of data |
void isce3::signal::Signal< T >::forward2DFFT | ( | std::valarray< T > & | signal, |
std::valarray< std::complex< T > > & | spectrum, | ||
int | incolumns, | ||
int | inrows, | ||
int | oncolumns, | ||
int | onrows ) |
initiate plan for forward two imensional FFT for a block of real data
[in] | signal | input block of data |
[out] | spectrum | output block of spectrum |
[in] | incolumns | number of columns of the block of data in input container |
[in] | inrows | number of rows of the block of data in input container |
[in] | oncolumns | number of columns of the block of data in output container |
[in] | oupnrows | number of rows of the block of data in output container |
If doing an out-of-place FFT, the output can be stored in a container whose number of columns and rows could be different than the input. This would be necessary for instance in the case of 2D upsampling
void isce3::signal::Signal< T >::forward2DFFT | ( | std::valarray< T > & | signal, |
std::valarray< std::complex< T > > & | spectrum, | ||
int | ncolumns, | ||
int | nrows ) |
initiate plan for forward two imensional FFT for a block of real data
[in] | signal | input block of data |
[out] | spectrum | output block of spectrum |
[in] | ncolumns | number of columns of the block of data |
[in] | nrows | number of rows of the block of data |
void isce3::signal::Signal< T >::forward2DFFT | ( | T * | signal, |
std::complex< T > * | spectrum, | ||
int | incolumns, | ||
int | inrows, | ||
int | oncolumns, | ||
int | onrows ) |
initiate plan for forward two imensional FFT for a block of real data
[in] | signal | input block of data |
[out] | spectrum | output block of spectrum |
[in] | incolumns | number of columns of the block of data in input container |
[in] | inrows | number of rows of the block of data in input container |
[in] | oncolumns | number of columns of the block of data in output container |
[in] | oupnrows | number of rows of the block of data in output container |
If doing an out-of-place FFT, the output can be stored in a container whose number of columns and rows could be different than the input. This would be necessary for instance in the case of 2D upsampling.
void isce3::signal::Signal< T >::forward2DFFT | ( | T * | signal, |
std::complex< T > * | spectrum, | ||
int | ncolumns, | ||
int | nrows ) |
initiate plan for forward two imensional FFT for a block of real data
[in] | signal | input block of data |
[out] | spectrum | output block of spectrum |
[in] | ncolumns | number of columns of the block of data |
[in] | nrows | number of rows of the block of data |
void isce3::signal::Signal< T >::forwardAzimuthFFT | ( | std::complex< T > * | signal, |
std::complex< T > * | spectrum, | ||
int | ncolumns, | ||
int | nrows ) |
initiate plan for forward FFT in azimuth direction for a block of complex data.
azimuth direction is assumed to be in the direction of the rows of the array.
[in] | signal | input block of data |
[out] | spectrum | output block of spectrum |
[in] | ncolumns | number of columns of the block of data |
[in] | nrows | number of rows of the block of data |
void isce3::signal::Signal< T >::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.
azimuth direction is assumed to be in the direction of the rows of the array.
[in] | signal | input block of data |
[out] | spectrum | output block of spectrum |
[in] | ncolumns | number of columns of the block of data |
[in] | nrows | number of rows of the block of data |
void isce3::signal::Signal< T >::forwardAzimuthFFT | ( | std::valarray< T > & | signal, |
std::valarray< std::complex< T > > & | spectrum, | ||
int | ncolumns, | ||
int | nrows ) |
initiate plan for forward FFT in azimuth direction for a block of real data.
azimuth direction is assumed to be in the direction of the rows of the array.
[in] | signal | input block of data |
[out] | spectrum | output block of spectrum |
[in] | ncolumns | number of columns of the block of data |
[in] | nrows | number of rows of the block of data |
void isce3::signal::Signal< T >::forwardAzimuthFFT | ( | T * | signal, |
std::complex< T > * | spectrum, | ||
int | ncolumns, | ||
int | nrows ) |
initiate plan for forward FFT in azimuth direction for a block of real data.
azimuth direction is assumed to be in the direction of the rows of the array.
[in] | signal | input block of data |
[out] | spectrum | output block of spectrum |
[in] | ncolumns | number of columns of the block of data |
[in] | nrows | number of rows of the block of data |
void isce3::signal::Signal< T >::forwardRangeFFT | ( | std::complex< T > * | signal, |
std::complex< T > * | spectrum, | ||
int | ncolumns, | ||
int | nrows ) |
initiate plan for forward FFT in range direction for a block of complex data.
range direction is assumed to be in the direction of the columns of the array.
[in] | signal | input block of data |
[out] | spectrum | output block of spectrum |
[in] | ncolumns | number of columns of the block of data |
[in] | nrows | number of rows of the block of data |
void isce3::signal::Signal< T >::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.
range direction is assumed to be in the direction of the columns of the array.
[in] | signal | input block of data |
[out] | spectrum | output block of spectrum |
[in] | ncolumns | number of columns of the block of data |
[in] | nrows | number of rows of the block of data |
void isce3::signal::Signal< T >::forwardRangeFFT | ( | std::valarray< T > & | signal, |
std::valarray< std::complex< T > > & | spectrum, | ||
int | ncolumns, | ||
int | nrows ) |
initiate plan for forward FFT in range direction for a block of real data.
range direction is assumed to be in the direction of the columns of the array.
[in] | signal | input block of data |
[out] | spectrum | output block of spectrum |
[in] | ncolumns | number of columns of the block of data |
[in] | nrows | number of rows of the block of data |
void isce3::signal::Signal< T >::forwardRangeFFT | ( | T * | signal, |
std::complex< T > * | spectrum, | ||
int | ncolumns, | ||
int | nrows ) |
initiate plan for forward FFT in range direction for a block of real data.
range direction is assumed to be in the direction of the columns of the array.
[in] | signal | input block of data |
[out] | spectrum | output block of spectrum |
[in] | ncolumns | number of columns of the block of data |
[in] | nrows | number of rows of the block of data |
void isce3::signal::Signal< T >::inverse | ( | std::complex< T > * | input, |
std::complex< T > * | output ) |
inverse transform
unnormalized inverse transform.
void isce3::signal::Signal< T >::inverse | ( | std::complex< T > * | input, |
T * | output ) |
inverse transform
unnormalized inverse transform.
void isce3::signal::Signal< T >::inverse | ( | std::valarray< std::complex< T > > & | input, |
std::valarray< std::complex< T > > & | output ) |
inverse transform
unnormalized inverse transform.
Note that since the FFTW library does not normalize the DFT computations, computing a forward followed by a backward transform (or vice versa) results in the original array scaled by length of fft.
[in] | input | block of spectrum |
[out] | output | block of data |
void isce3::signal::Signal< T >::inverse | ( | std::valarray< std::complex< T > > & | input, |
std::valarray< T > & | output ) |
inverse transform
unnormalized inverse transform.
void isce3::signal::Signal< T >::inverse2DFFT | ( | std::complex< T > * | spectrum, |
std::complex< T > * | signal, | ||
int | ncolumns, | ||
int | nrows ) |
initiate plan for inverse two dimensional FFT for a block of data
[in] | spectrum | input block of spectrum |
[out] | signal | output block of data |
[in] | ncolumns | number of columns of the block of data |
[in] | nrows | number of rows of the block of data |
void isce3::signal::Signal< T >::inverse2DFFT | ( | std::complex< T > * | spectrum, |
T * | signal, | ||
int | ncolumns, | ||
int | nrows ) |
initiate plan for inverse two dimensional FFT for a block of data
[in] | spectrum | input block of spectrum |
[out] | signal | output block of data |
[in] | ncolumns | number of columns of the block of data |
[in] | nrows | number of rows of the block of data |
void isce3::signal::Signal< T >::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
[in] | spectrum | input block of spectrum |
[out] | signal | output block of data |
[in] | ncolumns | number of columns of the block of data |
[in] | nrows | number of rows of the block of data |
void isce3::signal::Signal< T >::inverse2DFFT | ( | std::valarray< std::complex< T > > & | spectrum, |
std::valarray< T > & | signal, | ||
int | ncolumns, | ||
int | nrows ) |
initiate plan for inverse two dimensional FFT for a block of data
[in] | spectrum | input block of spectrum |
[out] | signal | output block of data |
[in] | ncolumns | number of columns of the block of data |
[in] | nrows | number of rows of the block of data |
void isce3::signal::Signal< T >::inverseAzimuthFFT | ( | std::complex< T > * | spectrum, |
std::complex< T > * | signal, | ||
int | ncolumns, | ||
int | nrows ) |
initiate plan for inverse FFT in azimuth direction for a block of data
[in] | spectrum | input block of spectrum |
[out] | signal | output block of data |
[in] | ncolumns | number of columns of the block of data |
[in] | nrows | number of rows of the block of data |
void isce3::signal::Signal< T >::inverseAzimuthFFT | ( | std::complex< T > * | spectrum, |
T * | signal, | ||
int | ncolumns, | ||
int | nrows ) |
initiate plan for inverse FFT in azimuth direction for a block of data
[in] | spectrum | input block of spectrum |
[out] | signal | output block of data |
[in] | ncolumns | number of columns of the block of data |
[in] | nrows | number of rows of the block of data |
void isce3::signal::Signal< T >::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
[in] | spectrum | input block of spectrum |
[out] | signal | output block of data |
[in] | ncolumns | number of columns of the block of data |
[in] | nrows | number of rows of the block of data |
void isce3::signal::Signal< T >::inverseAzimuthFFT | ( | std::valarray< std::complex< T > > & | spectrum, |
std::valarray< T > & | signal, | ||
int | ncolumns, | ||
int | nrows ) |
initiate plan for inverse FFT in azimuth direction for a block of data
[in] | spectrum | input block of spectrum |
[out] | signal | output block of data |
[in] | ncolumns | number of columns of the block of data |
[in] | nrows | number of rows of the block of data |
void isce3::signal::Signal< T >::inverseRangeFFT | ( | std::complex< T > * | spectrum, |
std::complex< T > * | signal, | ||
int | ncolumns, | ||
int | nrows ) |
initiate plan for backward FFT in range direction for a block of data
[in] | spectrum | input block of spectrum |
[out] | signal | output block of data |
[in] | ncolumns | number of columns of the block of data |
[in] | nrows | number of rows of the block of data |
void isce3::signal::Signal< T >::inverseRangeFFT | ( | std::complex< T > * | spectrum, |
T * | signal, | ||
int | ncolumns, | ||
int | nrows ) |
initiate plan for backward FFT in range direction for a block of data
[in] | spectrum | input block of spectrum |
[out] | signal | output block of data |
[in] | ncolumns | number of columns of the block of data |
[in] | nrows | number of rows of the block of data |
void isce3::signal::Signal< T >::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
[in] | spectrum | input block of spectrum |
[out] | signal | output block of data |
[in] | ncolumns | number of columns of the block of data |
[in] | nrows | number of rows of the block of data |
void isce3::signal::Signal< T >::inverseRangeFFT | ( | std::valarray< std::complex< T > > & | spectrum, |
std::valarray< T > & | signal, | ||
int | ncolumns, | ||
int | nrows ) |
initiate plan for backward FFT in range direction for a block of data
[in] | spectrum | input block of spectrum |
[out] | signal | output block of data |
[in] | ncolumns | number of columns of the block of data |
[in] | nrows | number of rows of the block of data |
|
inline |
next power of two
[in] | N | the actual length of a signal |
[in] | fftLength | next power of two |
void isce3::signal::Signal< T >::upsample | ( | isce3::core::EArray2D< std::complex< T > > & | signal, |
isce3::core::EArray2D< std::complex< T > > & | signalUpsampled, | ||
const isce3::core::EArray2D< std::complex< T > > & | shiftImpact = {} ) |
upsampling a basebanded block of data in range (columns) direction and shifting the upsampled signal by a constant.
The shift is applied by an inout linear phase term in frequency domain.
[in] | signal | input block of data |
[out] | signalUpsampled | output block of oversampled data |
[in] | shiftImpact | a linear phase term equivalent to a constant shift in time domain |
void isce3::signal::Signal< T >::upsample | ( | std::valarray< std::complex< T > > & | signal, |
std::valarray< std::complex< T > > & | signalUpsampled, | ||
int | rows, | ||
int | fft_size, | ||
int | upsampleFactor ) |
upsampling a block of data in range direction
[in] | signal | input block of data |
[in] | signalUpsampled | output block of oversampled data |
[in] | rows | number of rows of the block of input and upsampled data |
[in] | fft_size | number of columns of the block of input data |
[in] | upsampleFactor | upsampling factor |
void isce3::signal::Signal< T >::upsample | ( | std::valarray< std::complex< T > > & | signal, |
std::valarray< std::complex< T > > & | signalUpsampled, | ||
int | rows, | ||
int | fft_size, | ||
int | upsampleFactor, | ||
std::valarray< std::complex< T > > | shiftImpact ) |
upsampling a block of data in range direction and shifting the upsampled signal by a constant.
The shift is applied by an inout linear phase term in frequency domain.
[in] | signal | input block of data |
[out] | signalUpsampled | output block of oversampled data |
[in] | rows | number of rows of the block of input and upsampled data |
[in] | fft_size | number of columns of the block of input data |
[in] | upsampleFactor | upsampling factor |
[out] | shiftImpact | a linear phase term equivalent to a constant shift in time domain |
void isce3::signal::Signal< T >::upsample2D | ( | std::complex< T > * | signal, |
std::complex< T > * | signalUpsampled, | ||
int | upsampleFactor ) |
2D upsampling a block of 2D data
[in] | signal | input block of 2D data |
[out] | signalUpsampled | output block of oversampled data |
[in] | upsampleFactor | upsampling factor |
When doing out-of-place upsampling, i.e., when the container of the input to upsample is different to the one that will contain the upsampled data, the FFT forward plan must be set to out-of-place (with the input and output containers) and the FFT reverse plan to in-place (with the output container). In both case (in-place or out-of-place 2D upsampling), it is the user responsability to provide an output container that is padded (for in-place) or filled (for out-of-place) with zeros.
void isce3::signal::Signal< T >::upsample2D | ( | std::complex< T > * | signal, |
std::complex< T > * | signalUpsampled, | ||
int | upsampleFactor, | ||
std::complex< T > * | shiftImpact ) |
2D upsampling a block of 2D data and shifting the upsampled signal by a constant.
The shift is applied by an inout linear phase term in frequency domain.
[in] | signal | input block of 2D data |
[out] | signalUpsampled | output block of oversampled data |
[in] | upsampleFactor | upsampling factor |
[in] | shiftImpact | a linear phase term equivalent to a constant shift in time domain |
When doing out-of-place upsampling, i.e., when the container of the input to upsample is different to the one that will contain the upsampled data, the FFT forward plan must be set to out-of-place (with the input and output containers) and the FFT reverse plan to in-place (with the output container). In both case (in-place or out-of-place 2D upsampling), it is the user responsability to provide an output container that is padded (for in-place) or filled (for out-of-place) with zeros.
void isce3::signal::Signal< T >::upsample2D | ( | std::valarray< std::complex< T > > & | signal, |
std::valarray< std::complex< T > > & | signalUpsampled, | ||
int | oversampleFactor ) |
2D upsampling a block of 2D data
[in] | signal | input block of 2D data |
[out] | signalUpsampled | output block of oversampled data |
[in] | upsampleFactor | upsampling factor |
When doing out-of-place upsampling, i.e., when the container of the input to upsample is different to the one that will contain the upsampled data, the FFT forward plan must be set to out-of-place (with the input and output containers) and the FFT reverse plan to in-place (with the output container). In both case (in-place or out-of-place 2D upsampling), it is the user responsability to provide an output container that is padded (for in-place) or filled (for out-of-place) with zeros.
void isce3::signal::Signal< T >::upsample2D | ( | std::valarray< std::complex< T > > & | signal, |
std::valarray< std::complex< T > > & | signalUpsampled, | ||
int | upsampleFactor, | ||
std::valarray< std::complex< T > > & | shiftImpact ) |
2D upsampling a block of 2D data and shifting the upsampled signal by a constant.
The shift is applied by an inout linear phase term in frequency domain.
[in] | signal | input block of 2D data |
[out] | signalUpsampled | output block of oversampled data |
[in] | upsampleFactor | upsampling factor |
[in] | shiftImpact | a linear phase term equivalent to a constant shift in time domain |
When doing out-of-place upsampling, i.e., when the container of the input to upsample is different to the one that will contain the upsampled data, the FFT forward plan must be set to out-of-place (with the input and output containers) and the FFT reverse plan to in-place (with the output container). In both case (in-place or out-of-place 2D upsampling), it is the user responsability to provide an output container that is padded (for in-place) or filled (for out-of-place) with zeros.