d@0: d@0: d@0: Real-data DFTs - FFTW 3.2.1 d@0: d@0: d@0: d@0: d@0: d@0: d@0: d@0: d@0: d@0: d@0: d@0: d@0: d@0:
d@0:

d@0: d@0: d@0: Next: , d@0: Previous: Planner Flags, d@0: Up: Basic Interface d@0:


d@0:
d@0: d@0:

4.3.3 Real-data DFTs

d@0: d@0:
     fftw_plan fftw_plan_dft_r2c_1d(int n,
d@0:                                     double *in, fftw_complex *out,
d@0:                                     unsigned flags);
d@0:      fftw_plan fftw_plan_dft_r2c_2d(int n0, int n1,
d@0:                                     double *in, fftw_complex *out,
d@0:                                     unsigned flags);
d@0:      fftw_plan fftw_plan_dft_r2c_3d(int n0, int n1, int n2,
d@0:                                     double *in, fftw_complex *out,
d@0:                                     unsigned flags);
d@0:      fftw_plan fftw_plan_dft_r2c(int rank, const int *n,
d@0:                                  double *in, fftw_complex *out,
d@0:                                  unsigned flags);
d@0: 
d@0:

d@0: Plan a real-input/complex-output discrete Fourier transform (DFT) in d@0: zero or more dimensions, returning an fftw_plan (see Using Plans). d@0: d@0:

Once you have created a plan for a certain transform type and d@0: parameters, then creating another plan of the same type and parameters, d@0: but for different arrays, is fast and shares constant data with the d@0: first plan (if it still exists). d@0: d@0:

The planner returns NULL if the plan cannot be created. A d@0: non-NULL plan is always returned by the basic interface unless d@0: you are using a customized FFTW configuration supporting a restricted d@0: set of transforms, or if you use the FFTW_PRESERVE_INPUT flag d@0: with a multi-dimensional out-of-place c2r transform (see below). d@0: d@0:

Arguments
d@0: d@0: d@0: d@0:

The inverse transforms, taking complex input (storing the non-redundant d@0: half of a logically Hermitian array) to real output, are given by: d@0: d@0:

     fftw_plan fftw_plan_dft_c2r_1d(int n,
d@0:                                     fftw_complex *in, double *out,
d@0:                                     unsigned flags);
d@0:      fftw_plan fftw_plan_dft_c2r_2d(int n0, int n1,
d@0:                                     fftw_complex *in, double *out,
d@0:                                     unsigned flags);
d@0:      fftw_plan fftw_plan_dft_c2r_3d(int n0, int n1, int n2,
d@0:                                     fftw_complex *in, double *out,
d@0:                                     unsigned flags);
d@0:      fftw_plan fftw_plan_dft_c2r(int rank, const int *n,
d@0:                                  fftw_complex *in, double *out,
d@0:                                  unsigned flags);
d@0: 
d@0:

d@0: The arguments are the same as for the r2c transforms, except that the d@0: input and output data formats are reversed. d@0: d@0:

FFTW computes an unnormalized transform: computing an r2c followed by a d@0: c2r transform (or vice versa) will result in the original data d@0: multiplied by the size of the transform (the product of the logical d@0: dimensions). d@0: An r2c transform produces the same output as a FFTW_FORWARD d@0: complex DFT of the same input, and a c2r transform is correspondingly d@0: equivalent to FFTW_BACKWARD. For more information, see What FFTW Really Computes. d@0: d@0: d@0: d@0: