Fourier library: Operations in the Fourier domain
fftwaff.cc
Go to the documentation of this file.
1 
42 #define TF_FFTWAFF_CC_VERSION \
43  "TF_FFTWAFF_CC V1.4"
44 
45 #include <iostream>
46 #include <fourier/fftwaff.h>
47 #include <fourier/error.h>
48 #include <fourier/error.h>
49 #include <aff/seriesoperators.h>
50 #include <aff/dump.h>
51 #include <cmath>
52 
53 namespace fourier {
54 
57  namespace fft {
58 
61  {
62  if (Mplan_forward==0)
63  {
64  //Mplan_forward=rfftw_create_plan(Msize, FFTW_FORWARD, 0);
65 #ifdef FFTWFALLBACK
66  Mplan_forward=rfftw_create_plan(Msize, FFTW_REAL_TO_COMPLEX,
67  FFTW_ESTIMATE);
68 #else
69  this->create_arrays();
70  Mplan_forward=fftw_plan_dft_r2c_1d(Msize, Mseriesarray,
72  FFTW_ESTIMATE);
73 #endif
75  "Error (DRFFTWAFF::create_plan_forward): "
76  "could not create plan!")
77  }
78  } // void DRFFTWAFF::create_plan_forward() const
79 
80  /*----------------------------------------------------------------------*/
81 
84  {
85  if (Mplan_backward==0)
86  {
87 #ifdef FFTWFALLBACK
88  Mplan_backward=rfftw_create_plan(Msize, FFTW_BACKWARD, 0);
89 #else
90  this->create_arrays();
91  Mplan_backward=fftw_plan_dft_c2r_1d(Msize, Mspectrumarray,
93  FFTW_ESTIMATE);
94 #endif
96  "Error (DRFFTWAFF::create_plan_backward): "
97  "could not create plan!")
98  }
99  } // void DRFFTWAFF::create_plan_backward() const
100 
101  /*----------------------------------------------------------------------*/
102 
105  {
106  this->delete_plans();
107 #ifndef FFTWFALLBACK
108  fftw_cleanup();
109 #endif
110  } // DRFFTWAFF::~DRFFTWAFF()
111 
112  /*----------------------------------------------------------------------*/
113 
115  void DRFFTWAFF::set_size(const unsigned int& n) const
116  {
117  if (n != Msize)
118  {
119  Msize=n;
120  this->delete_plans();
121 #ifndef FFTWFALLBACK
122  this->delete_arrays();
123 #endif
124  }
125  } // DRFFTWAFF::~DRFFTWAFF()
126 
127  /*----------------------------------------------------------------------*/
128 
131  {
132 #ifdef FFTWFALLBACK
133  if (Mplan_forward != 0) { rfftw_destroy_plan(Mplan_forward); }
134  if (Mplan_backward != 0) { rfftw_destroy_plan(Mplan_backward); }
135 #else
136  if (Mplan_forward != 0) { fftw_destroy_plan(Mplan_forward); }
137  if (Mplan_backward != 0) { fftw_destroy_plan(Mplan_backward); }
138 #endif
139  Mplan_forward=0;
140  Mplan_backward=0;
141  } // DRFFTWAFF::delete_plans()
142 
143  /*----------------------------------------------------------------------*/
144 
145 #ifndef FFTWFALLBACK
146  void DRFFTWAFF::create_arrays() const
148  {
149  this->delete_arrays();
150  Mseriesarray = (double *) fftw_malloc(sizeof(double)*Msize);
152  "Error (DRFFTWAFF::create_plan_forward): "
153  "could not create series array!")
154  Mseries=Tseries(Tseries::Trepresentation(Mseriesarray, Msize));
155  Mspectrumarray = (fftw_complex *)
156  fftw_malloc(sizeof(fftw_complex)*this->ssize());
158  "Error (DRFFTWAFF::create_plan_forward): "
159  "could not create spectrum array!")
160  Mspectrum=Tspectrum(Tspectrum::Trepresentation(
161  reinterpret_cast<Tcoeff*>(Mspectrumarray),
162  this->ssize()));
163  } // DRFFTWAFF::create_arrays()
164 
165  /*----------------------------------------------------------------------*/
166 
168  void DRFFTWAFF::delete_arrays() const
169  {
170  this->delete_plans();
171  if (Mseriesarray != 0) { fftw_free(Mseriesarray); }
172  if (Mspectrumarray != 0) { fftw_free(Mspectrumarray); }
173  Mseries=Tseries(0);
174  Mspectrum=Tspectrum(0);
175  // Mplan_forward=0;
176  // Mplan_backward=0;
177  } // DRFFTWAFF::delete_arrays()
178 #endif
179 
180  /*----------------------------------------------------------------------*/
181 
187  const bool& debug) const
188  {
189  this->set_size(s.size());
190 #ifdef FFTWFALLBACK
191  FOURIER_debug(debug, "DRFFTWAFF::operator()",
192  "use fallback code (FFTW2)");
193  Tspectrum retval(this->Msize/2+1);
194  aff::Series<fftw_real> out(this->Msize);
195  aff::Series<fftw_real> in(this->Msize);
196  fftw_real* pout=out.pointer();
197  fftw_real* pin=in.pointer();
198  FOURIER_debug(debug, "DRFFTWAFF::operator()",
199  "processing arrays are created; copy in series");
200  in.copyin(s);
201 
202  FOURIER_debug(debug, "DRFFTWAFF::operator()",
203  "create plan forward");
204  this->create_plan_forward();
205  FOURIER_debug(debug, "DRFFTWAFF::operator()",
206  "call rfftw_one");
207  rfftw_one(Mplan_forward, pin, pout);
208  FOURIER_debug(debug, "DRFFTWAFF::operator()",
209  "copy results to output");
210  retval(0)=out(0);
211  for (int i=1; i<((Msize+1)/2); ++i)
212  {
213  retval(i)=Tcoeff(out(i),out(Msize-i));
214  }
215  if ((Msize % 2) == 0)
216  {
217  retval(Msize/2)=out(Msize/2);
218  }
219  /*
220  if (debug)
221  {
222  DUMP(in);
223  DUMP(out);
224  }
225  */
226 #else
227  FOURIER_debug(debug, "DRFFTWAFF::operator()",
228  "use recent code (FFTW3)");
229  Tspectrum retval(this->ssize());
230  this->create_plan_forward();
231  Mseries.copyin(s);
232  /*
233  for (int i=0; i<this->size(); ++i)
234  { Mseries(Mseries.first()+i)=s(s.first()+i); }
235  */
236  fftw_execute(Mplan_forward);
237  retval.copyin(Mspectrum);
238  /*
239  for (int i=0; <this->ssize(); ++i)
240  { retval(retval.first()+i)=Mspectrum(Mspectrum.first()+i); }
241  */
242 #endif
243  /*
244  if (debug)
245  {
246  DUMP(s);
247  DUMP(retval);
248  }
249  */
250  FOURIER_debug(debug, "DRFFTWAFF::operator()",
251  "finished; return");
252  return(retval);
253  } // Tspectrum DRFFTWAFF::operator()(const Tseries::Tcoc& s) const
254 
255  /*----------------------------------------------------------------------*/
256 
261  DRFFTWAFF::Tseries DRFFTWAFF::operator()(const Tspectrum::Tcoc& s,
262  const bool& debug) const
263  {
264  // check number of expected Fourier coefficients
265  if (this->ssize() != s.size())
266  {
267  // adjust FFT size
268  int seriessize=DRFFTWAFF::seriessize(s.size());
269  // is Nyquits coefficients real?
270  // this measure can only be effective if the Nyquist coefficient is
271  // finite, which will not be the case for most signals
272  if (std::abs(s(s.size()).imag()) < 1.e-8*std::abs(s(s.size()).real()))
273  { --seriessize; }
274  this->set_size(seriessize);
275  }
276 
277 #ifdef FFTWFALLBACK
278  Tseries retval(Msize);
279  aff::Series<fftw_real> out(Msize);
280  aff::Series<fftw_real> in(Msize);
281  fftw_real* pout=out.pointer();
282  fftw_real* pin=in.pointer();
283  in(0)=s(0).real();
284  for (int i=1; i<((Msize+1)/2); ++i)
285  {
286  in(i)=s(i).real();
287  in(Msize-i)=s(i).imag();
288  }
289  if ((Msize % 2) == 0)
290  {
291  in(Msize/2)=s(Msize/2).real();
292  }
293  this->create_plan_backward();
294  rfftw_one(Mplan_backward, pin, pout);
295  retval.copyin(out);
296 #else
297  Tseries retval(this->size());
298  this->create_plan_backward();
299  Mspectrum.copyin(s);
300  fftw_execute(Mplan_backward);
301  retval.copyin(Mseries);
302 #endif
303  return(retval);
304  } // Tseries DRFFTWAFF::operator()(const Tspectrum::Tcoc& s) const
305 
306  /*----------------------------------------------------------------------*/
307 
318  {
319  return(1./(Msize*dt));
320  } // Tsample DRFFTWAFF::scale_series(const Tsample& dt) const
321 
322  /*----------------------------------------------------------------------*/
323 
334  {
335  return(dt);
336  } // Tsample DRFFTWAFF::scale_spectrum(const Tsample& dt) const
337 
338  /*----------------------------------------------------------------------*/
339 
345  const double& dt,
346  const bool& debug) const
347  {
348  Tspectrum retval=this->operator()(s, debug);
349  retval *= this->scale_spectrum(dt);
350  return(retval);
351  }
352 
353  /*----------------------------------------------------------------------*/
354 
359  DRFFTWAFF::Tseries DRFFTWAFF::operator()(const Tspectrum::Tcoc& s,
360  const double& dt,
361  const bool& debug) const
362  {
363  return(this->scale_series(dt)*
364  this->operator()(s, debug));
365  }
366 
367  } // namespace ftt
368 
369 } // namespace fourier
370 
371 /* ----- END OF fftwaff.cc ----- */
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
#define FOURIER_assert(C, M)
Check an assertion and report by throwing an exception.
Definition: error.h:141
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
error handling for libfourier (prototypes)
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
#define FOURIER_debug(C, N, M)
produce debug output
Definition: error.h:194
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
use fftw together with aff containers (prototypes)
~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