Next: , Previous: Guru vector and transform sizes, Up: Guru Interface


4.5.3 Guru Complex DFTs

     fftw_plan fftw_plan_guru_dft(
          int rank, const fftw_iodim *dims,
          int howmany_rank, const fftw_iodim *howmany_dims,
          fftw_complex *in, fftw_complex *out,
          int sign, unsigned flags);
     
     fftw_plan fftw_plan_guru_split_dft(
          int rank, const fftw_iodim *dims,
          int howmany_rank, const fftw_iodim *howmany_dims,
          double *ri, double *ii, double *ro, double *io,
          unsigned flags);

These two functions plan a complex-data, multi-dimensional DFT for the interleaved and split format, respectively. Transform dimensions are given by (rank, dims) over a multi-dimensional vector (loop) of dimensions (howmany_rank, howmany_dims). dims and howmany_dims should point to fftw_iodim arrays of length rank and howmany_rank, respectively.

flags is a bitwise OR (‘|’) of zero or more planner flags, as defined in Planner Flags.

In the fftw_plan_guru_dft function, the pointers in and out point to the interleaved input and output arrays, respectively. The sign can be either -1 (= FFTW_FORWARD) or +1 (= FFTW_BACKWARD). If the pointers are equal, the transform is in-place.

In the fftw_plan_guru_split_dft function, ri and ii point to the real and imaginary input arrays, and ro and io point to the real and imaginary output arrays. The input and output pointers may be the same, indicating an in-place transform. For example, for fftw_complex pointers in and out, the corresponding parameters are:

     ri = (double *) in;
     ii = (double *) in + 1;
     ro = (double *) out;
     io = (double *) out + 1;

Because fftw_plan_guru_split_dft accepts split arrays, strides are expressed in units of double. For a contiguous fftw_complex array, the overall stride of the transform should be 2, the distance between consecutive real parts or between consecutive imaginary parts; see Guru vector and transform sizes. Note that the dimension strides are applied equally to the real and imaginary parts; real and imaginary arrays with different strides are not supported.

There is no sign parameter in fftw_plan_guru_split_dft. This function always plans for an FFTW_FORWARD transform. To plan for an FFTW_BACKWARD transform, you can exploit the identity that the backwards DFT is equal to the forwards DFT with the real and imaginary parts swapped. For example, in the case of the fftw_complex arrays above, the FFTW_BACKWARD transform is computed by the parameters:

     ri = (double *) in + 1;
     ii = (double *) in;
     ro = (double *) out + 1;
     io = (double *) out;