This chapter describes the basic usage of FFTW for computing complex one-dimensional Fourier Transforms. This chapter tells the truth, but not the whole truth. See section FFTW Reference, for more information.
Using FFTW is very simple. First of all, include <fftw.h>
at the top of your program:
#include <fftw.h>
Then, in order to use FFTW, you must first create a plan. A plan
is an object that contains all the data FFTW needs to compute the
Fourier Transform. The plan is created by fftw_create_plan
:
fftw_plan fftw_create_plan(int n, fftw_direction dir, int flags);
The function fftw_create_plan
accepts three arguments. The first
argument, n
, is the size of the transform you are trying to
compute. It can be any integer greater than 0. The second argument,
dir
, can be either FFTW_FORWARD
or FFTW_BACKWARD
,
and indicates which transform you are interested in. Alternatively, you
can use the sign of the exponent in the transform, or
, which correspond to FFTW_FORWARD
and
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 your clock. FFTW_ESTIMATE
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.
Once the plan has been created, you can use it as many times as you like
for arrays of the same size. The fftw
function accepts several
arguments:
void fftw(fftw_plan plan, int howmany, FFTW_COMPLEX *in, int istride, int idist, FFTW_COMPLEX *out, int ostride, int odist);
This means: compute howmany
FFTs(1) according to the plan plan
(the size of
the transform is specified by the plan). The input array is pointed to
by in
, and consists of howmany
subarrays, placed at a
distance idist
from each other, and whose elements are
istride
apart. In other words, the i
-th element of the
j
-th subarray will be in position in[i * istride + j *
idist]
, where i
and j
are both zero-based.. The result
will be stored in the out
array, which consists of howmany
subarrays, placed at a distance odist
from each other, and whose
elements are ostride
apart. The input and output arrays
must be distinct. For example, in order to compute just one transform
of an array, use
fftw(plan, 1, in, 1, 0, out, 1, 0);
Finally, don't forget to destroy the plan once you are done using it:
void fftw_destroy_plan(fftw_plan plan)
This frees all the memory and data structures associated with the plan.
FFTW uses the FFTPACK/Matlab convention that . In other words, applying the forward and then the backward transform will multiply the input by .
More precisely, the forward transform of an array of size computes an array , where
The backward transform computes
The following simple program illustrates the intended use of FFTW (for one-dimensional transforms).
#include <fftw.h> #include <stdlib.h> #include <stdio.h> void example(void) { int n = 16; int i; FFTW_COMPLEX *in = malloc(n * sizeof(FFTW_COMPLEX)); FFTW_COMPLEX *out = malloc(n * sizeof(FFTW_COMPLEX)); fftw_plan plan; plan = fftw_create_plan(n, FFTW_FORWARD, FFTW_MEASURE); if (plan == NULL) { fprintf(stderr, "Error creating plan\n"); exit(-1); } /* initialize to arbitrary values */ for (i = 0; i < n; ++i) { c_re(in[i]) = i; c_im(in[i]) = -i; } fftw(plan, 1, in, 1, 0, out, 1, 0); /* now use the output array */ /* ... */ }
The program should be linked with -lfftw -lm
.