fftw3: Real-data DFTs

 
 4.3.3 Real-data DFTs
 --------------------
 
      fftw_plan fftw_plan_dft_r2c_1d(int n0,
                                     double *in, fftw_complex *out,
                                     unsigned flags);
      fftw_plan fftw_plan_dft_r2c_2d(int n0, int n1,
                                     double *in, fftw_complex *out,
                                     unsigned flags);
      fftw_plan fftw_plan_dft_r2c_3d(int n0, int n1, int n2,
                                     double *in, fftw_complex *out,
                                     unsigned flags);
      fftw_plan fftw_plan_dft_r2c(int rank, const int *n,
                                  double *in, fftw_complex *out,
                                  unsigned flags);
 
    Plan a real-input/complex-output discrete Fourier transform (DFT) in
 zero or more dimensions, returning an 'fftw_plan' (SeeUsing Plans).
 
    Once you have created a plan for a certain transform type and
 parameters, then creating another plan of the same type and parameters,
 but for different arrays, is fast and shares constant data with the
 first plan (if it still exists).
 
    The planner returns 'NULL' if the plan cannot be created.  A
 non-'NULL' plan is always returned by the basic interface unless you are
 using a customized FFTW configuration supporting a restricted set of
 transforms, or if you use the 'FFTW_PRESERVE_INPUT' flag with a
 multi-dimensional out-of-place c2r transform (see below).
 
 Arguments
 .........
 
    * 'rank' is the rank of the transform (it should be the size of the
      array '*n'), and can be any non-negative integer.  (SeeComplex
      Multi-Dimensional DFTs, for the definition of "rank".)  The
      '_1d', '_2d', and '_3d' planners correspond to a 'rank' of '1',
      '2', and '3', respectively.  The rank may be zero, which is
      equivalent to a rank-1 transform of size 1, i.e.  a copy of one
      real number (with zero imaginary part) from input to output.
 
    * 'n0', 'n1', 'n2', or 'n[0..rank-1]', (as appropriate for each
      routine) specify the size of the transform dimensions.  They can be
      any positive integer.  This is different in general from the
      _physical_ array dimensions, which are described in SeeReal-data
      DFT Array Format.
 
         - FFTW is best at handling sizes of the form 2^a 3^b 5^c 7^d
           11^e 13^f, where e+f is either 0 or 1, and the other exponents
           are arbitrary.  Other sizes are computed by means of a slow,
           general-purpose algorithm (which nevertheless retains O(n log
           n) performance even for prime sizes).  (It is possible to
           customize FFTW for different array sizes; see See
           Installation and Customization.)  Transforms whose sizes are
           powers of 2 are especially fast, and it is generally
           beneficial for the _last_ dimension of an r2c/c2r transform to
           be _even_.
 
    * 'in' and 'out' point to the input and output arrays of the
      transform, which may be the same (yielding an in-place transform).
      These arrays are overwritten during planning, unless
      'FFTW_ESTIMATE' is used in the flags.  (The arrays need not be
      initialized, but they must be allocated.)  For an in-place
      transform, it is important to remember that the real array will
      require padding, described in SeeReal-data DFT Array Format.
 
    * 'flags' is a bitwise OR ('|') of zero or more planner flags, as
      defined in SeePlanner Flags.
 
    The inverse transforms, taking complex input (storing the
 non-redundant half of a logically Hermitian array) to real output, are
 given by:
 
      fftw_plan fftw_plan_dft_c2r_1d(int n0,
                                     fftw_complex *in, double *out,
                                     unsigned flags);
      fftw_plan fftw_plan_dft_c2r_2d(int n0, int n1,
                                     fftw_complex *in, double *out,
                                     unsigned flags);
      fftw_plan fftw_plan_dft_c2r_3d(int n0, int n1, int n2,
                                     fftw_complex *in, double *out,
                                     unsigned flags);
      fftw_plan fftw_plan_dft_c2r(int rank, const int *n,
                                  fftw_complex *in, double *out,
                                  unsigned flags);
 
    The arguments are the same as for the r2c transforms, except that the
 input and output data formats are reversed.
 
    FFTW computes an unnormalized transform: computing an r2c followed by
 a c2r transform (or vice versa) will result in the original data
 multiplied by the size of the transform (the product of the logical
 dimensions).  An r2c transform produces the same output as a
 'FFTW_FORWARD' complex DFT of the same input, and a c2r transform is
 correspondingly equivalent to 'FFTW_BACKWARD'.  For more information,
 see SeeWhat FFTW Really Computes.