Fourier library: Operations in the Fourier domain

◆ operator()() [1/4]

DRFFTWAFF::Tspectrum fourier::fft::DRFFTWAFF::operator() ( const Tseries::Tcoc &  s,
const bool &  debug = false 
) const

Transform time series to Fourier coefficients.

No scaling is applied.

Definition at line 186 of file fftwaff.cc.

References create_plan_forward(), FOURIER_debug, Mplan_forward, Mseries, Msize, Mspectrum, set_size(), and ssize().

Referenced by operator()().

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
DRFFTWAFF::Tspectrum Tspectrum
Definition: cxxfftwtest.cc:49
std::complex< Tsample > Tcoeff
Definition: fftwaff.h:132
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
unsigned int Msize
Definition: fftwaff.h:187
unsigned int ssize() const
Definition: fftwaff.h:183
Here is the call graph for this function:
Here is the caller graph for this function: