Diana Software
QRealComplexFFTW3.hh
Go to the documentation of this file.
1 
2 #ifndef __Q_REAL_COMPLEX_FFTW3_
3 #define __Q_REAL_COMPLEX_FFTW3_
4 
5 #include "QFFT.hh"
6 #include "QFFTW3.hh"
7 
8 // Portability Issue: If the environment cannot find the FFTW3
9 // headers, QRealComplexFFTW3 just becomes a wrapper for the
10 // QRealComplexFFTGSL class which uses gsl.
11 #if defined(_HAVE_FFTW3_) && !defined(__CINT__) && !defined(__CLING__)
12 #include <fftw3.h>
13 #endif
14 
16 
17 class QVectorC;
18 
54 class QRealComplexFFTW3: public QFFT, virtual public QFFTW3 {
55 
56 public:
61  QRealComplexFFTW3(const size_t size);
69  virtual ~QRealComplexFFTW3();
76  virtual int TransformToFreq(const QVector& data, QVectorC& FFT,
77  bool compress=false);
84  virtual int TransformFromFreq(const QVectorC& FT, QVector& spectrum,
85  bool compress=false);
86 
87 
92  virtual void SetWindowType(WindowType wt, int coherent = 0);
93 
94 protected:
95 
97  virtual bool CreatePlan(bool Forward);
105  virtual bool ClearCache(bool clearMem=true);
106 
125  virtual int Transform(const QVector& data, QVector& result);
126 
149  bool CheckPlanCompatible(const QVector& data,const QVector& result) const
150  { return false; }
151 
157  bool Resize(size_t s);
158 
159 
160 protected:
161 
165  bool fSizeLock;
166 
170  const int fStride;
173 
174 #if defined(_HAVE_FFTW3_) && !defined(__CINT__) && !defined(__CLING__)
183  fftw_plan fFFTW_Forward_Plan;
187  fftw_plan fFFTW_Backward_Plan;
188 
195  fftw_complex *fComplexData;
196 
197 #endif
198 
200 };
202 
203 
204 #endif
#define Q_END_NAMESPACE
Definition: QDiana.hh:22
#define Q_BEGIN_NAMESPACE
Definition: QDiana.hh:20
Base class for FFTW3 classes.
Definition: QFFTW3.hh:130
Interface for ffts in Diana analysis.
Definition: QFFT.hh:26
WindowType
window type
Definition: QFFT.hh:30
Real to complex FFT using FFTW3.
virtual bool ClearCache(bool clearMem=true)
Free the allocated memory.
ClassDef(QRealComplexFFTW3, 0)
bool Resize(size_t s)
Resize the allocated memory.
size_t fAllocatedSize
Size of the allocated array.
const int fStride
The stride of the arrays.
virtual ~QRealComplexFFTW3()
destructor
virtual void SetWindowType(WindowType wt, int coherent=0)
resize working table and space
virtual bool CreatePlan(bool Forward)
Create a plan for this FFT (or any of this size)
virtual int Transform(const QVector &data, QVector &result)
Perform the transformation.
virtual int TransformFromFreq(const QVectorC &FT, QVector &spectrum, bool compress=false)
transform from the frequencies to the times
bool CheckPlanCompatible(const QVector &data, const QVector &result) const
Check to see if the fftw can be done without copying the data.
QRealComplexFFTW3()
empty constructor
bool fSizeLock
Lock the FFT size.
virtual int TransformToFreq(const QVector &data, QVectorC &FFT, bool compress=false)
transform from the times to the frequencies
Interface for complex vectors in Diana analysis.
Definition: QVectorC.hh:25
Interface for vectors in Diana analysis.
Definition: QVector.hh:30