GWOLLUM  4.2.0
Tools for gravitational-wave analyses
fft Class Reference

Wrap and optimize FFTW. More...

#include <FFT.h>

Public Member Functions

bool Backward (double *aRe_f)
 Performs the backward FFT of a frequency-domain purely real data vector. More...
 
bool Backward (double *aRe_f, double *aIm_f)
 Performs the backward FFT of a frequency-domain complex data vector. More...
 
void Backward (void)
 Performs the backward FFT. More...
 
bool Forward (double *aRe_t)
 Performs the forward FFT of a time-domain purely real data vector. More...
 
bool Forward (double *aRe_t, double *aIm_t)
 Performs the forward FFT of a time-domain complex data vector. More...
 
void Forward (void)
 Performs the forward FFT. More...
 
double * GetIm_f (const double aNorm=1.0)
 Returns the imaginary part of the frequency-domain data vector. More...
 
double GetIm_f (const unsigned int aIndex)
 Returns the imaginary part of a frequency-domain data vector element. More...
 
double * GetIm_t (const double aNorm=1.0)
 Returns the imaginary part of the time-domain data vector. More...
 
double GetIm_t (const unsigned int aIndex)
 Returns the imaginary part of a time-domain data vector element. More...
 
double * GetNorm2_f (const double aNorm=1.0)
 Returns the squared norm of the frequency-domain data vector. More...
 
double GetNorm2_f (const unsigned int aIndex)
 Returns the squared norm of a frequency-domain data vector element. More...
 
double * GetNorm2_t (const double aNorm=1.0)
 Returns the squared norm of the time-domain data vector. More...
 
double GetNorm2_t (const unsigned int aIndex)
 Returns the squared norm of a time-domain data vector element. More...
 
double * GetNorm_f (const double aNorm=1.0)
 Returns the norm of the frequency-domain data vector. More...
 
double GetNorm_f (const unsigned int aIndex)
 Returns the norm of a frequency-domain data vector element. More...
 
double * GetNorm_t (const double aNorm=1.0)
 Returns the norm of the time-domain data vector. More...
 
double GetNorm_t (const unsigned int aIndex)
 Returns the norm of a time-domain data vector element. More...
 
double GetPhase_f (const unsigned int aIndex)
 Returns the phase of a frequency-domain data vector element. More...
 
double * GetPhase_f (void)
 Returns the phase of the frequency-domain data vector. More...
 
double GetPhase_t (const unsigned int aIndex)
 Returns the phase of a time-domain data vector element. More...
 
double * GetPhase_t (void)
 Returns the phase of the time-domain data vector. More...
 
double * GetRe_f (const double aNorm=1.0)
 Returns the real part of the frequency-domain data vector. More...
 
double GetRe_f (const unsigned int aIndex)
 Returns the real part of a frequency-domain data vector element. More...
 
double * GetRe_t (const double aNorm=1.0)
 Returns the real part of the time-domain data vector. More...
 
double GetRe_t (const unsigned int aIndex)
 Returns the real part of a time-domain data vector element. More...
 
unsigned int GetSize_f (void)
 Returns the size of the frequency-domain vector. More...
 
unsigned int GetSize_t (void)
 Returns the size of the time-domain vector. More...
 
void SetIm_f (const unsigned int aIndex, const double aValue)
 Sets a frequency-domain vector element (imaginary part) to a given value. More...
 
void SetIm_t (const unsigned int aIndex, const double aValue)
 Sets a time-domain vector element (imaginary part) to a given value. More...
 
void SetRe_f (const unsigned int aIndex, const double aValue)
 Sets a frequency-domain vector element (real part) to a given value. More...
 
void SetRe_t (const unsigned int aIndex, const double aValue)
 Sets a time-domain vector element (real part) to a given value. More...
 
Constructors and destructors
 fft (const unsigned int aSize_t, const string aPlan="FFTW_ESTIMATE", const string aType="c2c")
 Constructor of the fft class. More...
 
virtual ~fft (void)
 Destructor of the fft class. More...
 

Private Attributes

fftw_plan FFT_Backward
 Backward plan. More...
 
fftw_plan FFT_Forward
 Forward plan. More...
 
string fPlan
 FFTW plan. More...
 
unsigned int fSize_f
 Frequency-domain vector size. More...
 
unsigned int fSize_t
 Time-domain vector size. More...
 
bool fType
 FFT type true=c2c, false=r2c. More...
 
fftw_complex * x_f
 Frequency-domain data. More...
 
fftw_complex * xc_t
 Time-domain data - complex. More...
 
double * xr_t
 Time-domain data - real. More...
 

Detailed Description

Wrap and optimize FFTW.

See also
FFTW.
Author
Florent Robinet

Constructor & Destructor Documentation

◆ fft()

fft::fft ( const unsigned int  aSize_t,
const string  aPlan = "FFTW_ESTIMATE",
const string  aType = "c2c" 
)

Constructor of the fft class.

A FFT plan is contructed for a given time-domain vector size. The input size must be an even number or else it will be forced as such.

A FFT type must be provided. The default "c2c" type provides a complex-to-complex FFT. In this case, the time-domain and the frequency domain vectors are both complex and have the same size. If a "r2c" type is selected, the time-domain vector is purely real. The frequency-domain vector is complex and its size is "aSize_t/2+1".

Warning
When using the "r2c" type, the backward transform will alter the frequency-domain vector. Do not use it after the transform.

Several FFTW plans are supported: "FFTW_ESTIMATE", "FFTW_MEASURE" and "FFTW_PATIENT".

Parameters
[in]aSize_tTime-domain vector size: must be an even number.
[in]aPlanFFT plan.
[in]aTypeFFT type "r2c" or "c2c".

◆ ~fft()

fft::~fft ( void  )
virtual

Destructor of the fft class.

Member Function Documentation

◆ Backward() [1/3]

bool fft::Backward ( double *  aRe_f)

Performs the backward FFT of a frequency-domain purely real data vector.

The input data vector is copied before being transformed.

Parameters
[in]aRe_fFrequency-domain real data vector.

◆ Backward() [2/3]

bool fft::Backward ( double *  aRe_f,
double *  aIm_f 
)

Performs the backward FFT of a frequency-domain complex data vector.

The input data vector is defined by 2 arrays for the real and imaginary parts. The input data vector is copied before being transformed.

Parameters
[in]aRe_fFrequency-domain data vector (real part).
[in]aIm_fFrequency-domain data vector (imaginary part).
Precondition
The input pointers must point to arrays the size of which must match the definition given in the constructor.

◆ Backward() [3/3]

void fft::Backward ( void  )
inline

Performs the backward FFT.

◆ Forward() [1/3]

bool fft::Forward ( double *  aRe_t)

Performs the forward FFT of a time-domain purely real data vector.

The input data vector is copied before being transformed.

Note
This function can be called both in a "c2c" and a "r2c" mode. In a "c2c" mode, the imaginary part is set to 0.
Parameters
[in]aRe_tTime-domain real data vector.
Precondition
The input pointer must point to an array the size of which must match the definition given in the constructor.

◆ Forward() [2/3]

bool fft::Forward ( double *  aRe_t,
double *  aIm_t 
)

Performs the forward FFT of a time-domain complex data vector.

The input data vector is defined by 2 arrays; the real and imaginary parts. The input data vector is copied before being transformed.

Note
This function can only be called in a "c2c" mode.
Returns
false if the "c2c" mode was not selected in the constructor; true if success.
Parameters
[in]aRe_tTime-domain data vector (real part).
[in]aIm_tTime-domain data vector (imaginary part).
Precondition
The input pointers must point to arrays the size of which must match the definition given in the constructor.

◆ Forward() [3/3]

void fft::Forward ( void  )
inline

Performs the forward FFT.

◆ GetIm_f() [1/2]

double * fft::GetIm_f ( const double  aNorm = 1.0)

Returns the imaginary part of the frequency-domain data vector.

The user is in charge of deleting the returned array.

Parameters
[in]aNormThe returned vector is normalized by this parameter.

◆ GetIm_f() [2/2]

double fft::GetIm_f ( const unsigned int  aIndex)
inline

Returns the imaginary part of a frequency-domain data vector element.

Parameters
[in]aIndexVector index.

◆ GetIm_t() [1/2]

double * fft::GetIm_t ( const double  aNorm = 1.0)

Returns the imaginary part of the time-domain data vector.

The user is in charge of deleting the returned array.

Parameters
[in]aNormThe returned vector is normalized by this parameter.

◆ GetIm_t() [2/2]

double fft::GetIm_t ( const unsigned int  aIndex)
inline

Returns the imaginary part of a time-domain data vector element.

Parameters
[in]aIndexvector index

◆ GetNorm2_f() [1/2]

double * fft::GetNorm2_f ( const double  aNorm = 1.0)

Returns the squared norm of the frequency-domain data vector.

The user is in charge of deleting the returned array.

Parameters
[in]aNormThe returned vector is normalized by this parameter.

◆ GetNorm2_f() [2/2]

double fft::GetNorm2_f ( const unsigned int  aIndex)
inline

Returns the squared norm of a frequency-domain data vector element.

Parameters
[in]aIndexVector index.

◆ GetNorm2_t() [1/2]

double * fft::GetNorm2_t ( const double  aNorm = 1.0)

Returns the squared norm of the time-domain data vector.

The user is in charge of deleting the returned array.

Parameters
[in]aNormThe returned vector is normalized by this parameter.

◆ GetNorm2_t() [2/2]

double fft::GetNorm2_t ( const unsigned int  aIndex)
inline

Returns the squared norm of a time-domain data vector element.

Parameters
[in]aIndexVector index.

◆ GetNorm_f() [1/2]

double * fft::GetNorm_f ( const double  aNorm = 1.0)

Returns the norm of the frequency-domain data vector.

The user is in charge of deleting the returned array.

Parameters
[in]aNormThe returned vector is normalized by this parameter.

◆ GetNorm_f() [2/2]

double fft::GetNorm_f ( const unsigned int  aIndex)
inline

Returns the norm of a frequency-domain data vector element.

Parameters
[in]aIndexVector index.

◆ GetNorm_t() [1/2]

double * fft::GetNorm_t ( const double  aNorm = 1.0)

Returns the norm of the time-domain data vector.

The user is in charge of deleting the returned array.

Parameters
[in]aNormThe returned vector is normalized by this parameter.

◆ GetNorm_t() [2/2]

double fft::GetNorm_t ( const unsigned int  aIndex)
inline

Returns the norm of a time-domain data vector element.

Note
For "r2c", prefer the use of GetRe_t() (faster).
Parameters
[in]aIndexVector index.

◆ GetPhase_f() [1/2]

double fft::GetPhase_f ( const unsigned int  aIndex)
inline

Returns the phase of a frequency-domain data vector element.

The returned phase is between \(-\pi\) and \(-\pi\).

Parameters
[in]aIndexVector index.

◆ GetPhase_f() [2/2]

double * fft::GetPhase_f ( void  )

Returns the phase of the frequency-domain data vector.

The user is in charge of deleting the returned array.

◆ GetPhase_t() [1/2]

double fft::GetPhase_t ( const unsigned int  aIndex)
inline

Returns the phase of a time-domain data vector element.

The returned phase is between \(-\pi\) and \(-\pi\).

Parameters
[in]aIndexvector index

◆ GetPhase_t() [2/2]

double * fft::GetPhase_t ( void  )

Returns the phase of the time-domain data vector.

The user is in charge of deleting the returned array.

◆ GetRe_f() [1/2]

double * fft::GetRe_f ( const double  aNorm = 1.0)

Returns the real part of the frequency-domain data vector.

The user is in charge of deleting the returned array.

Parameters
[in]aNormThe returned vector is normalized by this parameter.

◆ GetRe_f() [2/2]

double fft::GetRe_f ( const unsigned int  aIndex)
inline

Returns the real part of a frequency-domain data vector element.

Parameters
[in]aIndexVector index.

◆ GetRe_t() [1/2]

double * fft::GetRe_t ( const double  aNorm = 1.0)

Returns the real part of the time-domain data vector.

The user is in charge of deleting the returned array.

Parameters
[in]aNormThe returned vector is normalized by this parameter.

◆ GetRe_t() [2/2]

double fft::GetRe_t ( const unsigned int  aIndex)
inline

Returns the real part of a time-domain data vector element.

Parameters
[in]aIndexVector index.

◆ GetSize_f()

unsigned int fft::GetSize_f ( void  )
inline

Returns the size of the frequency-domain vector.

◆ GetSize_t()

unsigned int fft::GetSize_t ( void  )
inline

Returns the size of the time-domain vector.

◆ SetIm_f()

void fft::SetIm_f ( const unsigned int  aIndex,
const double  aValue 
)
inline

Sets a frequency-domain vector element (imaginary part) to a given value.

Parameters
[in]aIndexVector index.
[in]aValueNew value.

◆ SetIm_t()

void fft::SetIm_t ( const unsigned int  aIndex,
const double  aValue 
)
inline

Sets a time-domain vector element (imaginary part) to a given value.

Parameters
[in]aIndexVector index.
[in]aValueNew value.

◆ SetRe_f()

void fft::SetRe_f ( const unsigned int  aIndex,
const double  aValue 
)
inline

Sets a frequency-domain vector element (real part) to a given value.

Parameters
[in]aIndexVector index.
[in]aValueNew value.

◆ SetRe_t()

void fft::SetRe_t ( const unsigned int  aIndex,
const double  aValue 
)
inline

Sets a time-domain vector element (real part) to a given value.

Parameters
[in]aIndexVector index.
[in]aValueNew value.

Member Data Documentation

◆ FFT_Backward

fftw_plan fft::FFT_Backward
private

Backward plan.

◆ FFT_Forward

fftw_plan fft::FFT_Forward
private

Forward plan.

◆ fPlan

string fft::fPlan
private

FFTW plan.

◆ fSize_f

unsigned int fft::fSize_f
private

Frequency-domain vector size.

◆ fSize_t

unsigned int fft::fSize_t
private

Time-domain vector size.

◆ fType

bool fft::fType
private

FFT type true=c2c, false=r2c.

◆ x_f

fftw_complex* fft::x_f
private

Frequency-domain data.

◆ xc_t

fftw_complex* fft::xc_t
private

Time-domain data - complex.

◆ xr_t

double* fft::xr_t
private

Time-domain data - real.


The documentation for this class was generated from the following files: