Diana Software
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
QRealComplexFFTW3 Class Reference

Real to complex FFT using FFTW3. More...

Inheritance diagram for QRealComplexFFTW3:
QFFT QFFTW3

Public Member Functions

 QRealComplexFFTW3 (const size_t size)
 constructor More...
 
 QRealComplexFFTW3 ()
 empty constructor More...
 
virtual ~QRealComplexFFTW3 ()
 destructor More...
 
virtual int TransformToFreq (const QVector &data, QVectorC &FFT, bool compress=false)
 transform from the times to the frequencies More...
 
virtual int TransformFromFreq (const QVectorC &FT, QVector &spectrum, bool compress=false)
 transform from the frequencies to the times More...
 
virtual void SetWindowType (WindowType wt, int coherent=0)
 resize working table and space More...
 
- Public Member Functions inherited from QFFT
 QFFT (const size_t size)
 constructor with its size More...
 
 QFFT ()
 empty constructor with its size More...
 
virtual ~QFFT ()
 destructor More...
 
int GetWindowType () const
 Return the window type. More...
 
void SetForward (bool isForward)
 virtual method. Must be implemented by daughter classes set the direction of the transform. More...
 
void SetNormalized (bool isNormalized)
 virtual method. Must be implemented by daughter classes set the normalization of the transform. More...
 
- Public Member Functions inherited from QFFTW3
 QFFTW3 ()
 ctor loads the wisdom singleton. More...
 
virtual ~QFFTW3 ()
 dtor does nothing. More...
 
bool ExportWisdomToFile (const char *outFile)
 Save the accumulated wisdom to file. More...
 
bool ImportWisdomFromFile (const char *inFile)
 Load wisdom from file. More...
 
bool ImportSystemWisdom ()
 Load wisdom from a system file. More...
 
bool ForgetWisdom (void)
 Forget any accumulated wisdom. More...
 
 ClassDef (QFFTW3, 0)
 

Protected Member Functions

virtual bool CreatePlan (bool Forward)
 Create a plan for this FFT (or any of this size) More...
 
virtual bool ClearCache (bool clearMem=true)
 Free the allocated memory. More...
 
virtual int Transform (const QVector &data, QVector &result)
 Perform the transformation. More...
 
bool CheckPlanCompatible (const QVector &data, const QVector &result) const
 Check to see if the fftw can be done without copying the data. More...
 
bool Resize (size_t s)
 Resize the allocated memory. More...
 
 ClassDef (QRealComplexFFTW3, 0)
 
- Protected Member Functions inherited from QFFT
void IsPowerOf2 ()
 
 ClassDef (QFFT, 0)
 

Protected Attributes

bool fSizeLock
 Lock the FFT size. More...
 
const int fStride
 The stride of the arrays. More...
 
size_t fAllocatedSize
 Size of the allocated array. More...
 
- Protected Attributes inherited from QFFT
size_t fSize
 
FFTDirection fDirection
 
bool fNormalized
 
bool fIsPowerOf2
 
double * fData
 
double * fWindow
 
size_t fWindowSize
 
WindowType fWindowType
 

Additional Inherited Members

- Public Types inherited from QFFT
enum  WindowType {
  WT_None = 0 , WT_Welch = 1 , WT_Rect = 2 , WT_Hann = 3 ,
  WT_Hamming = 4 , WT_Cosinus = 5 , WT_BlackmanHarris = 6 , WT_Kaiser3 = 7
}
 window type More...
 
enum  ZeroPadSide { kMiddle = 0 , kSymmetric = 1 , kLeft = 2 , kRight = 3 }
 
enum  FFTDirection { kForward = 1 , kBackward = 0 }
 
- Static Public Member Functions inherited from QFFT
static QVector FFTAntiSym (const QVector &input)
 antisimmetrize time domain vector More...
 
static QVector FFTSym (const QVector &input)
 simmetrize time domain vector More...
 
static QVector ZeroPad (const QVector &input, size_t n_zeros, int Side=kMiddle, double zeroVal=0.)
 Add zeros to the input. More...
 
static QVectorC ZeroPad (const QVectorC &input, size_t n_zeros, int Side=kMiddle, Diana::QComplex zeroVal=Diana::QComplex(0, 0))
 Add zeros to the input. More...
 
static WindowType StrToWindowType (const std::string &winName)
 Convert string to window type. More...
 
static QVector GetWindow (WindowType wt, size_t size, size_t zeros=0, int coherent=0, int param=0)
 create window and add zeros/2 to the left and zeros/2 to the right More...
 
static QVector CutSides (const QVector &input, size_t ncut, bool isSym)
 cut left and right sides by ncut/2 More...
 
static bool IsPowerOf2 (UInt_t Value)
 Is the input a power of 2. More...
 
static bool IsPowerOf2 (UInt_t Value, UInt_t &Power, UInt_t &Remainder)
 Is the input a power of 2? More...
 

Detailed Description

Real to complex FFT using FFTW3.

Author
J. Ouellet

This is an adaptation of QRealComplexFFT to use a FFTW3 library. This carries with it several advantages.

First, if an FFT of on size is to be performed many many times, this algorithm does not need to reallocate memory multiple times (so long as the instance is persistant, the memory is deallocated upon deletion of the object).

The bigger advantage is the use of FFTW3, since FFTW3 makes a plan for how it will perform and FFT and stores that plan the time to complete the FFT is much much faster (once the plan has been estabilished). FFTW3 needs to make two plans, one plan for a forward FFT, one plan for a backward FFT, these are made only when needed. Once the plans are made, they persist for the lifetime of the program (not just the object). This is the FFTW3 wisdom. The biggest speed up is for arrays whose size is a multiple of a large prime (like say 626=2*313). Here FFTW3 chooses the fastest algorithm to perform the FFT and can result in a factor of 10 speedup of the code.

Limitations: For a single instance of the object, the size of the array is fixed. If you wish to perform an FFT of another size, you must create a new instance of QRealComplexFFTW3. (Reallocate the memory, and replan the FFT).

For this class to be implemented, the os must be able to find the FFTW3 headers. If those are not found, this class will just be an alias of QRealComplexFFT using the gsl headers.

Definition at line 54 of file QRealComplexFFTW3.hh.

Constructor & Destructor Documentation

◆ QRealComplexFFTW3() [1/2]

QRealComplexFFTW3::QRealComplexFFTW3 ( const size_t  size)

constructor

Parameters
sizeis the length of the real vector to be transformed

Definition at line 11 of file QRealComplexFFTW3.cc.

References fAllocatedSize, QFFT::fDirection, QFFT::fNormalized, QFFT::fSize, fSizeLock, QFFT::fWindow, QFFT::fWindowType, QFFT::kForward, QERR_UNKNOWN_ERR, Resize(), and QFFT::WT_None.

◆ QRealComplexFFTW3() [2/2]

QRealComplexFFTW3::QRealComplexFFTW3 ( )

◆ ~QRealComplexFFTW3()

QRealComplexFFTW3::~QRealComplexFFTW3 ( )
virtual

destructor

Definition at line 52 of file QRealComplexFFTW3.cc.

References ClearCache().

Member Function Documentation

◆ CheckPlanCompatible()

bool QRealComplexFFTW3::CheckPlanCompatible ( const QVector data,
const QVector result 
) const
inlineprotected

Check to see if the fftw can be done without copying the data.

Normally, this class creates a plan which FFTs data back and forth between fData<->fComplexData. However, this means that data must be copied between the original vector and the local memory slot. This check determines if the fft can be done without moving the data.

There are a few requirements:

  1. The real/complex vectors must have the same stride as fData/fComplexData (in this case fStride=1).
  2. The real/complex vectors must have the same alignment in memory as fData/fComplex data. The easiest way to achieve this is using fftw_malloc as is implemented in QVector.

However, as of FFTW 3.3.3, this code actually caused significant slow down. Currently, this function returns false, forcing the memory copy.

Definition at line 149 of file QRealComplexFFTW3.hh.

◆ ClassDef()

QRealComplexFFTW3::ClassDef ( QRealComplexFFTW3  ,
 
)
protected

◆ ClearCache()

bool QRealComplexFFTW3::ClearCache ( bool  clearMem = true)
protectedvirtual

Free the allocated memory.

This occurs on destruction of the object and resizing

Parameters
clearMemIf true, also free the memory

Definition at line 238 of file QRealComplexFFTW3.cc.

References QVector::ArrayFree(), and QFFT::fData.

Referenced by Resize(), and ~QRealComplexFFTW3().

◆ CreatePlan()

bool QRealComplexFFTW3::CreatePlan ( bool  Forward)
protectedvirtual

Create a plan for this FFT (or any of this size)

Definition at line 165 of file QRealComplexFFTW3.cc.

References QFFT::fData, and QFFT::fSize.

Referenced by TransformFromFreq(), and TransformToFreq().

◆ Resize()

bool QRealComplexFFTW3::Resize ( size_t  s)
protectedvirtual

Resize the allocated memory.

If fSizeLock == true this returns false.

Implements QFFT.

Definition at line 207 of file QRealComplexFFTW3.cc.

References QVector::ArrayAlloc(), ClearCache(), fAllocatedSize, QFFT::fData, QFFT::fSize, and fSizeLock.

Referenced by QRealComplexFFTW3(), TransformFromFreq(), and TransformToFreq().

◆ SetWindowType()

void QRealComplexFFTW3::SetWindowType ( WindowType  wt,
int  coherent = 0 
)
virtual

resize working table and space

Parameters
sis the size of the real array to be transformed

Implements QFFT.

Definition at line 259 of file QRealComplexFFTW3.cc.

References QFFT::fSize, QFFT::fWindow, QFFT::fWindowSize, QFFT::fWindowType, and QFFT::GetWindow().

Referenced by QRealComplexFFT::SetWindowType(), TransformFromFreq(), and TransformToFreq().

◆ Transform()

int QRealComplexFFTW3::Transform ( const QVector data,
QVector result 
)
protectedvirtual

Perform the transformation.

This method takes the input data vector and performs a forward or backward transform to the result vector. If this is a forward transformation result must have size 2*fSize. If this is a reverse transformation, data needs to be (at least) size fSize+(fSize % 2 ? 1 : 2) i.e. the next even number larger than fSize.

This method is different from QRealComplexFFTGSL, in that the latter class requires the complex vector to be of size 2*fSize, and be formatted with the real terms first, followed by the complex terms. On the other hand, this class requires that the complex vector be interleaved. (This is more efficient since, that is the way the data are arranged after both gsl and fftw algorithms.)

Implements QFFT.

Definition at line 162 of file QRealComplexFFTW3.cc.

◆ TransformFromFreq()

int QRealComplexFFTW3::TransformFromFreq ( const QVectorC FT,
QVector spectrum,
bool  compress = false 
)
virtual

transform from the frequencies to the times

Parameters
resultstores the spectra size of the FFT is size of @FFF
Returns
0 upon success

Definition at line 110 of file QRealComplexFFTW3.cc.

References CreatePlan(), QFFT::fData, QFFT::fSize, QFFT::fWindow, QFFT::fWindowType, Im(), QERR_UNKNOWN_ERR, Re(), Resize(), QFFT::SetForward(), and SetWindowType().

Referenced by QRealComplexFFT::TransformFromFreq().

◆ TransformToFreq()

int QRealComplexFFTW3::TransformToFreq ( const QVector data,
QVectorC FFT,
bool  compress = false 
)
virtual

transform from the times to the frequencies

Parameters
resultstores the fft size of the FFT is size of @data
Returns
0 upon success

Definition at line 56 of file QRealComplexFFTW3.cc.

References CreatePlan(), QFFT::fData, QFFT::fSize, QFFT::fWindow, QFFT::fWindowType, Im(), QERR_UNKNOWN_ERR, Re(), Resize(), QFFT::SetForward(), and SetWindowType().

Referenced by QRealComplexFFT::TransformToFreq().

Member Data Documentation

◆ fAllocatedSize

size_t QRealComplexFFTW3::fAllocatedSize
protected

Size of the allocated array.

Definition at line 172 of file QRealComplexFFTW3.hh.

Referenced by QRealComplexFFTW3(), and Resize().

◆ fSizeLock

bool QRealComplexFFTW3::fSizeLock
protected

Lock the FFT size.

Definition at line 165 of file QRealComplexFFTW3.hh.

Referenced by QRealComplexFFTW3(), and Resize().

◆ fStride

const int QRealComplexFFTW3::fStride
protected

The stride of the arrays.

Definition at line 170 of file QRealComplexFFTW3.hh.


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