GWOLLUM
4.2.0
Tools for gravitational-wave analyses
|
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... | |
Wrap and optimize FFTW.
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".
Several FFTW plans are supported: "FFTW_ESTIMATE", "FFTW_MEASURE" and "FFTW_PATIENT".
[in] | aSize_t | Time-domain vector size: must be an even number. |
[in] | aPlan | FFT plan. |
[in] | aType | FFT type "r2c" or "c2c". |
|
virtual |
Destructor of the fft class.
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.
[in] | aRe_f | Frequency-domain real data vector. |
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.
[in] | aRe_f | Frequency-domain data vector (real part). |
[in] | aIm_f | Frequency-domain data vector (imaginary part). |
|
inline |
Performs the backward FFT.
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.
[in] | aRe_t | Time-domain real data vector. |
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.
[in] | aRe_t | Time-domain data vector (real part). |
[in] | aIm_t | Time-domain data vector (imaginary part). |
|
inline |
Performs the forward FFT.
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.
[in] | aNorm | The returned vector is normalized by this parameter. |
|
inline |
Returns the imaginary part of a frequency-domain data vector element.
[in] | aIndex | Vector index. |
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.
[in] | aNorm | The returned vector is normalized by this parameter. |
|
inline |
Returns the imaginary part of a time-domain data vector element.
[in] | aIndex | vector index |
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.
[in] | aNorm | The returned vector is normalized by this parameter. |
|
inline |
Returns the squared norm of a frequency-domain data vector element.
[in] | aIndex | Vector index. |
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.
[in] | aNorm | The returned vector is normalized by this parameter. |
|
inline |
Returns the squared norm of a time-domain data vector element.
[in] | aIndex | Vector index. |
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.
[in] | aNorm | The returned vector is normalized by this parameter. |
|
inline |
Returns the norm of a frequency-domain data vector element.
[in] | aIndex | Vector index. |
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.
[in] | aNorm | The returned vector is normalized by this parameter. |
|
inline |
Returns the norm of a time-domain data vector element.
[in] | aIndex | Vector index. |
|
inline |
Returns the phase of a frequency-domain data vector element.
The returned phase is between \(-\pi\) and \(-\pi\).
[in] | aIndex | Vector index. |
double * fft::GetPhase_f | ( | void | ) |
Returns the phase of the frequency-domain data vector.
The user is in charge of deleting the returned array.
|
inline |
Returns the phase of a time-domain data vector element.
The returned phase is between \(-\pi\) and \(-\pi\).
[in] | aIndex | vector index |
double * fft::GetPhase_t | ( | void | ) |
Returns the phase of the time-domain data vector.
The user is in charge of deleting the returned array.
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.
[in] | aNorm | The returned vector is normalized by this parameter. |
|
inline |
Returns the real part of a frequency-domain data vector element.
[in] | aIndex | Vector index. |
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.
[in] | aNorm | The returned vector is normalized by this parameter. |
|
inline |
Returns the real part of a time-domain data vector element.
[in] | aIndex | Vector index. |
|
inline |
Returns the size of the frequency-domain vector.
|
inline |
Returns the size of the time-domain vector.
|
inline |
Sets a frequency-domain vector element (imaginary part) to a given value.
[in] | aIndex | Vector index. |
[in] | aValue | New value. |
|
inline |
Sets a time-domain vector element (imaginary part) to a given value.
[in] | aIndex | Vector index. |
[in] | aValue | New value. |
|
inline |
Sets a frequency-domain vector element (real part) to a given value.
[in] | aIndex | Vector index. |
[in] | aValue | New value. |
|
inline |
Sets a time-domain vector element (real part) to a given value.
[in] | aIndex | Vector index. |
[in] | aValue | New value. |
|
private |
Backward plan.
|
private |
Forward plan.
|
private |
FFTW plan.
|
private |
Frequency-domain vector size.
|
private |
Time-domain vector size.
|
private |
FFT type true=c2c, false=r2c.
|
private |
Frequency-domain data.
|
private |
Time-domain data - complex.
|
private |
Time-domain data - real.