GWOLLUM 4.2.0
Tools for gravitational-wave analyses
Loading...
Searching...
No Matches
FFT.h
Go to the documentation of this file.
1
6#ifndef __fft__
7#define __fft__
8
9#include "CUtils.h"
10#include <fftw3.h>
11
12using namespace std;
13
19class fft {
20
21 public:
22
45 fft(const unsigned int aSize_t, const string aPlan = "FFTW_ESTIMATE", const string aType="c2c");
46
50 virtual ~fft(void);// destructor
58 inline void Forward(void){ fftw_execute(FFT_Forward); };
59
63 inline void Backward(void){ fftw_execute(FFT_Backward); };
64
75 bool Forward(double *aRe_t, double *aIm_t);
76
84 bool Forward(double *aRe_t);
85
94 bool Backward(double *aRe_f, double *aIm_f);
95
101 bool Backward(double *aRe_f);
102
107 inline double GetNorm2_f(const unsigned int aIndex){
108 return x_f[aIndex][0]*x_f[aIndex][0] + x_f[aIndex][1]*x_f[aIndex][1];
109 };
110
115 inline double GetNorm2_t(const unsigned int aIndex){
116 if(fType) return xc_t[aIndex][0]*xc_t[aIndex][0] + xc_t[aIndex][1]*xc_t[aIndex][1];
117 return xr_t[aIndex]*xr_t[aIndex];
118 };
119
124 inline double GetNorm_f(const unsigned int aIndex){
125 return sqrt(GetNorm2_f(aIndex));
126 };
127
133 inline double GetNorm_t(const unsigned int aIndex){
134 return sqrt(GetNorm2_t(aIndex));
135 };
136
141 inline double GetRe_f(const unsigned int aIndex){
142 return x_f[aIndex][0];
143 };
144
149 inline double GetRe_t(const unsigned int aIndex){
150 if(fType) return xc_t[aIndex][0];
151 return xr_t[aIndex];
152 };
153
158 inline double GetIm_f(const unsigned int aIndex){
159 return x_f[aIndex][1];
160 };
161
166 inline double GetIm_t(const unsigned int aIndex){
167 if(fType) return xc_t[aIndex][1];
168 return 0.0;
169 };
170
176 inline double GetPhase_f(const unsigned int aIndex){
177 return atan2(x_f[aIndex][1], x_f[aIndex][0]);
178 };
179
185 inline double GetPhase_t(const unsigned int aIndex){
186 if(fType) return atan2(xc_t[aIndex][1], xc_t[aIndex][0]);
187 else return 0.0;
188 };
189
195 double* GetNorm2_f(const double aNorm=1.0);
196
202 double* GetNorm2_t(const double aNorm=1.0);
203
209 double* GetNorm_f(const double aNorm=1.0);
210
216 double* GetNorm_t(const double aNorm=1.0);
217
223 double* GetRe_f(const double aNorm=1.0);
224
230 double* GetRe_t(const double aNorm=1.0);
231
237 double* GetIm_f(const double aNorm=1.0);
238
244 double* GetIm_t(const double aNorm=1.0);
245
250 double* GetPhase_f(void);
251
256 double* GetPhase_t(void);
257
261 inline unsigned int GetSize_t(void){ return fSize_t; };
262
266 inline unsigned int GetSize_f(void){ return fSize_f; };
267
273 inline void SetRe_t(const unsigned int aIndex, const double aValue){
274 if(fType) xc_t[aIndex][0]=aValue;
275 else xr_t[aIndex]=aValue;
276 };
277
283 inline void SetIm_t(const unsigned int aIndex, const double aValue){
284 xc_t[aIndex][1]=aValue;
285 };
286
292 inline void SetRe_f(const unsigned int aIndex, const double aValue){
293 x_f[aIndex][0]=aValue;
294 };
295
301 inline void SetIm_f(const unsigned int aIndex, const double aValue){
302 x_f[aIndex][1]=aValue;
303 };
304
305
306 private:
307 unsigned int fSize_t;
308 unsigned int fSize_f;
309 string fPlan;
310 bool fType;
311
312 fftw_complex* xc_t;
313 double* xr_t;
314 fftw_complex* x_f;
315 fftw_plan FFT_Forward;
316 fftw_plan FFT_Backward;
317};
318
319
320#endif
Generic C utility functions.
Wrap and optimize FFTW.
Definition FFT.h:19
double GetIm_f(const unsigned int aIndex)
Returns the imaginary part of a frequency-domain data vector element.
Definition FFT.h:158
fftw_complex * x_f
Frequency-domain data.
Definition FFT.h:314
fftw_plan FFT_Forward
Forward plan.
Definition FFT.h:315
void SetIm_t(const unsigned int aIndex, const double aValue)
Sets a time-domain vector element (imaginary part) to a given value.
Definition FFT.h:283
double * GetPhase_t(void)
Returns the phase of the time-domain data vector.
Definition FFT.cc:250
void Forward(void)
Performs the forward FFT.
Definition FFT.h:58
fftw_plan FFT_Backward
Backward plan.
Definition FFT.h:316
bool fType
FFT type true=c2c, false=r2c.
Definition FFT.h:310
unsigned int fSize_t
Time-domain vector size.
Definition FFT.h:307
double GetNorm2_f(const unsigned int aIndex)
Returns the squared norm of a frequency-domain data vector element.
Definition FFT.h:107
void SetRe_f(const unsigned int aIndex, const double aValue)
Sets a frequency-domain vector element (real part) to a given value.
Definition FFT.h:292
string fPlan
FFTW plan.
Definition FFT.h:309
double * xr_t
Time-domain data - real.
Definition FFT.h:313
unsigned int GetSize_f(void)
Returns the size of the frequency-domain vector.
Definition FFT.h:266
double GetRe_f(const unsigned int aIndex)
Returns the real part of a frequency-domain data vector element.
Definition FFT.h:141
double GetPhase_t(const unsigned int aIndex)
Returns the phase of a time-domain data vector element.
Definition FFT.h:185
double * GetPhase_f(void)
Returns the phase of the frequency-domain data vector.
Definition FFT.cc:242
void SetIm_f(const unsigned int aIndex, const double aValue)
Sets a frequency-domain vector element (imaginary part) to a given value.
Definition FFT.h:301
virtual ~fft(void)
Destructor of the fft class.
Definition FFT.cc:81
double GetRe_t(const unsigned int aIndex)
Returns the real part of a time-domain data vector element.
Definition FFT.h:149
double GetPhase_f(const unsigned int aIndex)
Returns the phase of a frequency-domain data vector element.
Definition FFT.h:176
double GetIm_t(const unsigned int aIndex)
Returns the imaginary part of a time-domain data vector element.
Definition FFT.h:166
void SetRe_t(const unsigned int aIndex, const double aValue)
Sets a time-domain vector element (real part) to a given value.
Definition FFT.h:273
double GetNorm2_t(const unsigned int aIndex)
Returns the squared norm of a time-domain data vector element.
Definition FFT.h:115
fftw_complex * xc_t
Time-domain data - complex.
Definition FFT.h:312
unsigned int GetSize_t(void)
Returns the size of the time-domain vector.
Definition FFT.h:261
unsigned int fSize_f
Frequency-domain vector size.
Definition FFT.h:308
double GetNorm_t(const unsigned int aIndex)
Returns the norm of a time-domain data vector element.
Definition FFT.h:133
double GetNorm_f(const unsigned int aIndex)
Returns the norm of a frequency-domain data vector element.
Definition FFT.h:124
void Backward(void)
Performs the backward FFT.
Definition FFT.h:63