GNU Info

Info Node: (fftw.info)Complex One-dimensional Transforms Tutorial

(fftw.info)Complex One-dimensional Transforms Tutorial


Next: Complex Multi-dimensional Transforms Tutorial Prev: Tutorial Up: Tutorial
Enter node , (file) or (file)node

Complex One-dimensional Transforms Tutorial
===========================================

   The basic usage of FFTW is simple.  A typical call to FFTW looks
like:

     #include <fftw.h>
     ...
     {
          fftw_complex in[N], out[N];
          fftw_plan p;
          ...
          p = fftw_create_plan(N, FFTW_FORWARD, FFTW_ESTIMATE);
          ...
          fftw_one(p, in, out);
          ...
          fftw_destroy_plan(p);
     }

   The first thing we do is to create a "plan", which is an object that
contains all the data that FFTW needs to compute the FFT, using the
following function:

     fftw_plan fftw_create_plan(int n, fftw_direction dir, int flags);

   The first argument, `n', is the size of the transform you are trying
to compute.  The size `n' can be any positive integer, but sizes that
are products of small factors are transformed most efficiently.  The
second argument, `dir', can be either `FFTW_FORWARD' or
`FFTW_BACKWARD', and indicates the direction of the transform you are
interested in.  Alternatively, you can use the sign of the exponent in
the transform, -1 or +1, which corresponds to `FFTW_FORWARD' or
`FFTW_BACKWARD' respectively.  The `flags' argument is either
`FFTW_MEASURE' or `FFTW_ESTIMATE'.  `FFTW_MEASURE' means that FFTW
actually runs and measures the execution time of several FFTs in order
to find the best way to compute the transform of size `n'.  This may
take some time, depending on your installation and on the precision of
the timer in your machine.  `FFTW_ESTIMATE', on the contrary, does not
run any computation, and just builds a reasonable plan, which may be
sub-optimal.  In other words, if your program performs many transforms
of the same size and initialization time is not important, use
`FFTW_MEASURE'; otherwise use the estimate.  (A compromise between
these two extremes exists.  Note: Words of Wisdom.)

   Once the plan has been created, you can use it as many times as you
like for transforms on arrays of the same size.  When you are done with
the plan, you deallocate it by calling `fftw_destroy_plan(plan)'.

   The transform itself is computed by passing the plan along with the
input and output arrays to `fftw_one':

     void fftw_one(fftw_plan plan, fftw_complex *in, fftw_complex *out);

   Note that the transform is out of place: `in' and `out' must point
to distinct arrays. It operates on data of type `fftw_complex', a data
structure with real (`in[i].re') and imaginary (`in[i].im')
floating-point components.  The `in' and `out' arrays should have the
length specified when the plan was created.  An alternative function,
`fftw', allows you to efficiently perform multiple and/or strided
transforms (Note: FFTW Reference.).

   The DFT results are stored in-order in the array `out', with the
zero-frequency (DC) component in `out[0]'.  The array `in' is not
modified.  Users should note that FFTW computes an unnormalized DFT,
the sign of whose exponent is given by the `dir' parameter of
`fftw_create_plan'.  Thus, computing a forward followed by a backward
transform (or vice versa) results in the original array scaled by `n'.
Note: What FFTW Really Computes, for the definition of DFT.

   A program using FFTW should be linked with `-lfftw -lm' on Unix
systems, or with the FFTW and standard math libraries in general.


automatically generated by info2www version 1.2.2.9