Fourier library: Operations in the Fourier domain
fftwaffar.h
Go to the documentation of this file.
1 
41 // include guard
42 #ifndef TF_FFTWAFFAR_H_VERSION
43 
44 #define TF_FFTWAFFAR_H_VERSION \
45  "TF_FFTWAFFAR_H V1.2"
46 
47 #include<complex>
48 #include<fftw3.h>
49 #include<aff/array.h>
50 #include<aff/series.h>
51 
52 namespace fourier {
53 
54  namespace fft {
55 
88  public:
89  typedef double Tsample;
90  typedef std::complex<Tsample> Tcoeff;
91  typedef aff::Array<Tsample> TAseries;
92  typedef aff::Array<Tcoeff> TAspectrum;
93  DRFFTWAFFArrayEngine(const int& nseis,
94  const int& nsamp);
96  const TAspectrum& spec);
101  void r2c();
102  void c2r();
103  Tsample scale_series(const Tsample& dt) const;
104  Tsample scale_spectrum(const Tsample& dt) const;
106  TAseries series() const { return Mseriesarray; }
108  TAspectrum spectrum() const { return Mspectrumarray; }
110  unsigned int nseries() const;
112  unsigned int nsamples() const;
114  unsigned int nfrequencies() const;
116  TAseries series(const unsigned int& i) const;
118  TAspectrum spectrum(const unsigned int& i) const;
119 
121  inline
122  static unsigned int ncoeff(const unsigned int& nsamples)
123  { return(nsamples/2+1); }
124 
126  inline
127  static unsigned int nsamples(const unsigned int& n)
128  { return(n*2-1); }
129 
130  private:
131  void checkconsistency();
132  void createplanr2c();
133  void createplanc2r();
134  void delete_plans();
137  fftw_plan Mplanr2c;
138  fftw_plan Mplanc2r;
139  }; // class DRFFTWAFFArrayEngine
140 
141  } // namespace fft
142 
143 } // namespace fourier
144 
145 #endif // TF_FFTWAFFAR_H_VERSION (includeguard)
146 
147 /* ----- END OF fftwaffar.h ----- */
unsigned int nsamples() const
return the number of samples in time series
Definition: fftwaffar.cc:287
Tsample scale_spectrum(const Tsample &dt) const
Return appropriate scaling factor for sampling interval dt.
Definition: fftwaffar.cc:280
aff::Array< Tsample > TAseries
Definition: fftwaffar.h:91
aff::Array< Tcoeff > TAspectrum
Definition: fftwaffar.h:92
~DRFFTWAFFArrayEngine()
delete plan.
Definition: fftwaffar.cc:96
TAseries series() const
return a reference to the time series arrays
Definition: fftwaffar.h:106
void r2c()
execute r2c plan
Definition: fftwaffar.cc:127
unsigned int nfrequencies() const
return the number of positive frequencies used
Definition: fftwaffar.cc:294
static unsigned int nsamples(const unsigned int &n)
return number of samples for given number of coefficients
Definition: fftwaffar.h:127
Tsample scale_series(const Tsample &dt) const
Return appropriate scaling factor for sampling interval dt.
Definition: fftwaffar.cc:263
Definition: error.cc:44
std::complex< Tsample > Tcoeff
Definition: fftwaffar.h:90
DRFFTWAFFArrayEngine & operator=(const DRFFTWAFFArrayEngine &e)
Definition: fftwaffar.cc:85
unsigned int nseries() const
return the number of series in the arrays
Definition: fftwaffar.cc:301
TAspectrum spectrum() const
return a reference to the Fourier coefficient arrays
Definition: fftwaffar.h:108
static unsigned int ncoeff(const unsigned int &nsamples)
return number of coefficients for given number of samples
Definition: fftwaffar.h:122
Engine to transform several signals at once. References to workspace are passed to the contstructor...
Definition: fftwaffar.h:87
void c2r()
execute c2r plan
Definition: fftwaffar.cc:115
void delete_plans()
delete plans.
Definition: fftwaffar.cc:104