Fourier library: Operations in the Fourier domain
|
#include <fftwaff.h>
Public Types | |
typedef double | Tsample |
typedef std::complex< Tsample > | Tcoeff |
typedef aff::Series< Tsample > | Tseries |
typedef aff::Series< Tcoeff > | Tspectrum |
Public Member Functions | |
DRFFTWAFF (const unsigned int &n, const bool &deletearrays=false) | |
DRFFTWAFF (const bool &deletearrays=false) | |
~DRFFTWAFF () | |
delete plan. More... | |
Tspectrum | operator() (const Tseries::Tcoc &s, const bool &debug=false) const |
Transform time series to Fourier coefficients. More... | |
Tseries | operator() (const Tspectrum::Tcoc &s, const bool &debug=false) const |
Transform Fourier coefficients to time series. More... | |
Tspectrum | operator() (const Tseries::Tcoc &s, const Tsample &dt, const bool &debug=false) const |
Transform time series to Fourier coefficients and scale. More... | |
Tseries | operator() (const Tspectrum::Tcoc &s, const Tsample &dt, const bool &debug=false) const |
Transform Fourier coefficients to time series and scale. More... | |
Tsample | scale_series (const Tsample &dt) const |
Return appropriate scaling factor for sampling interval dt. More... | |
Tsample | scale_spectrum (const Tsample &dt) const |
Return appropriate scaling factor for sampling interval dt. More... | |
unsigned int | size () const |
void | size (const unsigned int &s) const |
Static Public Member Functions | |
static unsigned int | spectrumsize (const unsigned int &n) |
return number of coefficients for given number of samples More... | |
static unsigned int | seriessize (const unsigned int &n) |
return number of samples for given number of coefficients More... | |
Private Member Functions | |
void | create_plan_forward () const |
create plan. More... | |
void | create_plan_backward () const |
create plan. More... | |
void | delete_plans () const |
delete plans. More... | |
void | create_arrays () const |
create plans. More... | |
void | delete_arrays () const |
delete arrays. More... | |
unsigned int | ssize () const |
void | set_size (const unsigned int &n) const |
prepare FFT settings for size n. More... | |
Private Attributes | |
unsigned int | Msize |
fftw_plan | Mplan_forward |
fftw_plan | Mplan_backward |
double * | Mseriesarray |
fftw_complex * | Mspectrumarray |
Tspectrum | Mspectrum |
Tseries | Mseries |
bool | Mdeletearrays |
A rigid class to do simple transforms using libdrfftw.a.
uses real double arrays
How to use this class:
You may create one instance of this class and use it to transform several signals in both directions in turn. The class itself takes care of the transform size and creates a new plan if necessary. FFTs are invoked by the bracket operators. You should use the class object like a function. The scaling operators (taking the sampling interval as one of their arguments) return a series or spectrum scaled appropriately to match the values of samples from the corresponding Fourier integral transform (usual convention with and
integrals - not
).