Fourier library: Operations in the Fourier domain
fftwaff.h
Go to the documentation of this file.
1 
72 // include guard
73 #ifndef TF_FFTWAFF_H_VERSION
74 
75 #define TF_FFTWAFF_H_VERSION \
76  "TF_FFTWAFF_H V1.3"
77 
78 #include<complex>
79 #ifdef FFTWFALLBACK
80 #include<drfftw.h>
81 #else
82 #include<fftw3.h>
83 #endif
84 #include<aff/series.h>
85 
86 namespace fourier {
87 
90  namespace fft {
91 
129  class DRFFTWAFF {
130  public:
131  typedef double Tsample;
132  typedef std::complex<Tsample> Tcoeff;
133  typedef aff::Series<Tsample> Tseries;
134  typedef aff::Series<Tcoeff> Tspectrum;
135 #ifdef FFTWFALLBACK
136  DRFFTWAFF(const unsigned int& n):
137  Msize(n), Mplan_forward(0), Mplan_backward(0) { }
138  DRFFTWAFF():
139  Msize(0), Mplan_forward(0), Mplan_backward(0) { }
140 #else
141  DRFFTWAFF(const unsigned int& n, const bool& deletearrays=false):
144  Mdeletearrays(deletearrays) { }
145  DRFFTWAFF(const bool& deletearrays=false):
148  Mdeletearrays(deletearrays) { }
149 #endif
150  ~DRFFTWAFF();
151  Tspectrum operator()(const Tseries::Tcoc& s,
152  const bool& debug=false) const;
153  Tseries operator()(const Tspectrum::Tcoc& s,
154  const bool& debug=false) const;
155  Tspectrum operator()(const Tseries::Tcoc& s,
156  const Tsample& dt,
157  const bool& debug=false) const;
158  Tseries operator()(const Tspectrum::Tcoc& s,
159  const Tsample& dt,
160  const bool& debug=false) const;
161  Tsample scale_series(const Tsample& dt) const;
162  Tsample scale_spectrum(const Tsample& dt) const;
163  unsigned int size() const { return(Msize); }
164  void size(const unsigned int& s) const { this->set_size(s); }
165 
167  inline
168  static unsigned int spectrumsize(const unsigned int& n)
169  { return(n/2+1); }
170 
172  inline
173  static unsigned int seriessize(const unsigned int& n)
174  { return(n*2-1); }
175 
176  private:
177  void create_plan_forward() const;
178  void create_plan_backward() const;
179  void delete_plans() const;
180 #ifndef FFTWFALLBACK
181  void create_arrays() const;
182  void delete_arrays() const;
183  unsigned int ssize() const
184  { return(DRFFTWAFF::spectrumsize(this->size())); }
185 #endif
186  void set_size(const unsigned int& n) const;
187  mutable unsigned int Msize;
188 #ifdef FFTWFALLBACK
189  mutable rfftw_plan Mplan_forward;
190  mutable rfftw_plan Mplan_backward;
191 #else
192  mutable fftw_plan Mplan_forward;
193  mutable fftw_plan Mplan_backward;
194  mutable double *Mseriesarray;
195  mutable fftw_complex *Mspectrumarray;
197  mutable Tseries Mseries;
199 #endif
200  }; // class DRFFTWAFF
201 
202  } // namespace ftt
203 
204 } // namespace fourier
205 
206 #endif // TF_FFTWAFF_H_VERSION (includeguard)
207 
208 /* ----- END OF fftwaff.h ----- */
DRFFTWAFF(const unsigned int &n, const bool &deletearrays=false)
Definition: fftwaff.h:141
Tsample scale_series(const Tsample &dt) const
Return appropriate scaling factor for sampling interval dt.
Definition: fftwaff.cc:317
aff::Series< Tsample > Tseries
Definition: fftwaff.h:133
void size(const unsigned int &s) const
Definition: fftwaff.h:164
DRFFTWAFF::Tspectrum Tspectrum
Definition: cxxfftwtest.cc:49
fftw_plan Mplan_backward
Definition: fftwaff.h:193
fftw_complex * Mspectrumarray
Definition: fftwaff.h:195
Tspectrum operator()(const Tseries::Tcoc &s, const bool &debug=false) const
Transform time series to Fourier coefficients.
Definition: fftwaff.cc:186
void create_plan_backward() const
create plan.
Definition: fftwaff.cc:83
unsigned int size() const
Definition: fftwaff.h:163
aff::Series< Tcoeff > Tspectrum
Definition: fftwaff.h:134
std::complex< Tsample > Tcoeff
Definition: fftwaff.h:132
void delete_plans() const
delete plans.
Definition: fftwaff.cc:130
static unsigned int seriessize(const unsigned int &n)
return number of samples for given number of coefficients
Definition: fftwaff.h:173
Definition: error.cc:44
void create_plan_forward() const
create plan.
Definition: fftwaff.cc:60
DRFFTWAFF(const bool &deletearrays=false)
Definition: fftwaff.h:145
fftw_plan Mplan_forward
Definition: fftwaff.h:192
void set_size(const unsigned int &n) const
prepare FFT settings for size n.
Definition: fftwaff.cc:115
aff::Series< Tvalue > Tseries
~DRFFTWAFF()
delete plan.
Definition: fftwaff.cc:104
Tsample scale_spectrum(const Tsample &dt) const
Return appropriate scaling factor for sampling interval dt.
Definition: fftwaff.cc:333
unsigned int Msize
Definition: fftwaff.h:187
unsigned int ssize() const
Definition: fftwaff.h:183
void create_arrays() const
create plans.
Definition: fftwaff.cc:147
void delete_arrays() const
delete arrays.
Definition: fftwaff.cc:168
static unsigned int spectrumsize(const unsigned int &n)
return number of coefficients for given number of samples
Definition: fftwaff.h:168