Diana Software
QFFT.hh
Go to the documentation of this file.
1 
11 #ifndef __QFFT_HH_
12 #define __QFFT_HH_
13 
14 #include "QVector.hh"
15 #include "QVectorC.hh"
16 #include "QComplex.hh"
17 
18 #include <gsl/gsl_fft_halfcomplex.h>
19 #include <gsl/gsl_fft_complex.h>
20 #include <gsl/gsl_errno.h>
21 #include <string>
22 #include <gsl/gsl_sf.h>
23 
25 
26 class QFFT {
27 
28 public:
30  enum WindowType {
31  WT_None = 0,
32  WT_Welch = 1,
33  WT_Rect = 2,
34  WT_Hann = 3,
39  };
40  enum ZeroPadSide {
41  kMiddle = 0,
43  kLeft = 2,
44  kRight = 3
45  };
46  enum FFTDirection {
47  kForward = 1,
48  kBackward = 0,
49  };
50 
51 
53  static QVector FFTAntiSym(const QVector& input);
55  static QVector FFTSym(const QVector& input);
73  static QVector ZeroPad(const QVector& input,size_t n_zeros, int Side = kMiddle, double zeroVal = 0.);
74 
92  static QVectorC ZeroPad(const QVectorC& input,size_t n_zeros, int Side = kMiddle, Diana::QComplex zeroVal = Diana::QComplex(0,0));
93 
95  static WindowType StrToWindowType(const std::string& winName);
97  static QVector GetWindow(WindowType wt, size_t size, size_t zeros = 0,int coherent = 0, int param = 0);
99  static QVector CutSides(const QVector& input, size_t ncut, bool isSym);
100 
105  QFFT( const size_t size);
109  QFFT();
113  virtual ~QFFT();
119  virtual int Transform(const QVector& data, QVector& result) = 0;
124  virtual bool Resize(size_t s) = 0;
130  virtual void SetWindowType(WindowType wt, int coherent = false) = 0;
135  int GetWindowType() const { return fWindowType; }
136 
142  void SetForward(bool isForward)
143  { fDirection = (isForward ? kForward : kBackward); }
149  void SetNormalized(bool isNormalized){fNormalized = isNormalized;}
150 
156  static bool IsPowerOf2(UInt_t Value);
166  static bool IsPowerOf2(UInt_t Value,UInt_t &Power,
167  UInt_t &Remainder);
168 
169 protected:
170 
171  void IsPowerOf2();
172  size_t fSize;
176  double* fData;
177  double* fWindow;
178  size_t fWindowSize;
180 
182 
183 };
184 
186 #endif
#define Q_END_NAMESPACE
Definition: QDiana.hh:22
#define Q_BEGIN_NAMESPACE
Definition: QDiana.hh:20
Interface for ffts in Diana analysis.
Definition: QFFT.hh:26
bool fIsPowerOf2
Definition: QFFT.hh:175
virtual int Transform(const QVector &data, QVector &result)=0
virtual method. Must be implemented by daughter classes.
static WindowType StrToWindowType(const std::string &winName)
Convert string to window type.
Definition: QFFT.cc:64
virtual void SetWindowType(WindowType wt, int coherent=false)=0
virtual method. Must be implemented by daughter classes set the window type.
size_t fWindowSize
Definition: QFFT.hh:178
void SetNormalized(bool isNormalized)
virtual method. Must be implemented by daughter classes set the normalization of the transform.
Definition: QFFT.hh:149
virtual bool Resize(size_t s)=0
virtual method. Must be implemented by daughter classes
virtual ~QFFT()
destructor
Definition: QFFT.cc:35
void SetForward(bool isForward)
virtual method. Must be implemented by daughter classes set the direction of the transform.
Definition: QFFT.hh:142
double * fData
Definition: QFFT.hh:176
WindowType
window type
Definition: QFFT.hh:30
@ WT_Rect
Definition: QFFT.hh:33
@ WT_Hann
Definition: QFFT.hh:34
@ WT_Welch
Definition: QFFT.hh:32
@ WT_Kaiser3
Definition: QFFT.hh:38
@ WT_Cosinus
Definition: QFFT.hh:36
@ WT_Hamming
Definition: QFFT.hh:35
@ WT_None
Definition: QFFT.hh:31
@ WT_BlackmanHarris
Definition: QFFT.hh:37
FFTDirection fDirection
Definition: QFFT.hh:173
WindowType fWindowType
Definition: QFFT.hh:179
static QVector ZeroPad(const QVector &input, size_t n_zeros, int Side=kMiddle, double zeroVal=0.)
Add zeros to the input.
Definition: QFFT.cc:216
static QVector CutSides(const QVector &input, size_t ncut, bool isSym)
cut left and right sides by ncut/2
Definition: QFFT.cc:290
ClassDef(QFFT, 0)
static QVector FFTAntiSym(const QVector &input)
antisimmetrize time domain vector
Definition: QFFT.cc:201
void IsPowerOf2()
Definition: QFFT.cc:61
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
Definition: QFFT.cc:82
QFFT()
empty constructor with its size
Definition: QFFT.cc:24
ZeroPadSide
Definition: QFFT.hh:40
@ kMiddle
Definition: QFFT.hh:41
@ kLeft
Definition: QFFT.hh:43
@ kSymmetric
Definition: QFFT.hh:42
@ kRight
Definition: QFFT.hh:44
static QVector FFTSym(const QVector &input)
simmetrize time domain vector
Definition: QFFT.cc:208
double * fWindow
Definition: QFFT.hh:177
size_t fSize
Definition: QFFT.hh:172
FFTDirection
Definition: QFFT.hh:46
@ kBackward
Definition: QFFT.hh:48
@ kForward
Definition: QFFT.hh:47
int GetWindowType() const
Return the window type.
Definition: QFFT.hh:135
bool fNormalized
Definition: QFFT.hh:174
Interface for complex vectors in Diana analysis.
Definition: QVectorC.hh:25
Interface for vectors in Diana analysis.
Definition: QVector.hh:30