LTL  2.0.x
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
ltl::FourierTransform< T, N > Class Template Reference

Public Member Functions

 FourierTransform ()
 constructor More...
 
virtual ~FourierTransform ()
 destructor More...
 
void FFT_Real2Complex (MArray< T, N > &A, MArray< complex< T >, N > &FFT_A)
 Real to complex (forward) fourier transform. Plan and execute the transform. Output array is not normalized. More...
 
void FFT_Complex2Real (MArray< complex< T >, N > &FFT_A, MArray< T, N > &A)
 Complex to real (inverse) fourier transform. Plan and execute the transform. Output array is not normalized. More...
 
void reset ()
 dispose plan and call fftw_cleanup(). More...
 
void FFT (MArray< complex< T >, N > &A, MArray< complex< T >, N > &FFT_A)
 Forward and inverse fourier transform. Plan and execute the transform. Output array is not normalized. More...
 
void iFFT (MArray< complex< T >, N > &FFT_A, MArray< complex< T >, N > &A)
 
void normalize (MArray< complex< T >, N > &A)
 normalize the array by dividing by the number of elements. More...
 
void normalize (MArray< T, N > &A)
 
MArray< T, N > shiftDC (MArray< T, N > &A)
 shift DC component to center of MArray. More...
 
MArray< complex< T >, N > shiftDC (MArray< complex< T >, N > &A)
 

Protected Member Functions

void execute ()
 

Protected Attributes

unsigned mode
 
fftw_plan plan
 

Detailed Description

template<typename T, int N>
class ltl::FourierTransform< T, N >

FourierTransform class

Interface to fftw3 FFT library. Currently, each call to one of the FFT methods causes a new plan to be generated. This is highly inefficient, but safe. It guarantees that pointers to the data strorage of MArray do not escape by being stored in a plan that might survive the MArray object.

All transforms leave the output array un-normalized. To normalize (divide by number of elements) call normalize().

* #include <ltl/fftw.h>
*
* // example of 1-dimensional FFT using fftw3 library:
* MArray<double, 1> in(size);
* MArray<std::complex<double>, 1> out(size);
* FourierTransform<double,1> FFT;
*
* in = 10.0 +
* 20.0 * sin(indexPosDbl(in,1) / double(size) * 2.0 * M_PI * 3.0) +
* 30.0 * cos(indexPosDbl(in,1) / double(size) * 2.0 * M_PI * 4.0);
* FFT.FFT_Real2Complex(in,out);
* FFT.normalize(out);
* out = merge(real(out * conj(out)) > 1e-9, out, 0);
* std::cout << out << std::endl;
*

Constructor & Destructor Documentation

template<typename T , int N>
ltl::FourierTransform< T, N >::FourierTransform ( )

constructor

template<typename T , int N>
virtual ltl::FourierTransform< T, N >::~FourierTransform ( )
virtual

destructor

Member Function Documentation

template<typename T , int N>
void ltl::FourierTransform< T, N >::FFT ( MArray< complex< T >, N > &  A,
MArray< complex< T >, N > &  FFT_A 
)

Forward and inverse fourier transform. Plan and execute the transform. Output array is not normalized.

template<typename T , int N>
void ltl::FourierTransform< T, N >::iFFT ( MArray< complex< T >, N > &  FFT_A,
MArray< complex< T >, N > &  A 
)
template<typename T , int N>
void ltl::FourierTransform< T, N >::FFT_Real2Complex ( MArray< T, N > &  A,
MArray< complex< T >, N > &  FFT_A 
)

Real to complex (forward) fourier transform. Plan and execute the transform. Output array is not normalized.

template<typename T , int N>
void ltl::FourierTransform< T, N >::FFT_Complex2Real ( MArray< complex< T >, N > &  FFT_A,
MArray< T, N > &  A 
)

Complex to real (inverse) fourier transform. Plan and execute the transform. Output array is not normalized.

template<typename T , int N>
void ltl::FourierTransform< T, N >::normalize ( MArray< complex< T >, N > &  A)

normalize the array by dividing by the number of elements.

template<typename T , int N>
void ltl::FourierTransform< T, N >::normalize ( MArray< T, N > &  A)
template<typename T , int N>
MArray<T,N> ltl::FourierTransform< T, N >::shiftDC ( MArray< T, N > &  A)

shift DC component to center of MArray.

template<typename T , int N>
MArray<complex<T>,N> ltl::FourierTransform< T, N >::shiftDC ( MArray< complex< T >, N > &  A)
template<typename T , int N>
void ltl::FourierTransform< T, N >::reset ( )
inline

dispose plan and call fftw_cleanup().

References ltl::FourierTransform< T, N >::mode.

template<typename T , int N>
void ltl::FourierTransform< T, N >::execute ( )
protected

Member Data Documentation

template<typename T , int N>
unsigned ltl::FourierTransform< T, N >::mode
protected
template<typename T , int N>
fftw_plan ltl::FourierTransform< T, N >::plan
protected