d@0: d@0:
d@0:d@0: d@0: d@0: Next: Complex Multi-Dimensional DFTs, d@0: Previous: Tutorial, d@0: Up: Tutorial d@0:
d@0: Plan: To bother about the best method of accomplishing an accidental result. d@0: [Ambrose Bierce, The Enlarged Devil's Dictionary.] d@0:d@0: d@0:
The basic usage of FFTW to compute a one-dimensional DFT of size
d@0: N
is simple, and it typically looks something like this code:
d@0:
d@0:
#include <fftw3.h>
d@0: ...
d@0: {
d@0: fftw_complex *in, *out;
d@0: fftw_plan p;
d@0: ...
d@0: in = (fftw_complex*) fftw_malloc(sizeof(fftw_complex) * N);
d@0: out = (fftw_complex*) fftw_malloc(sizeof(fftw_complex) * N);
d@0: p = fftw_plan_dft_1d(N, in, out, FFTW_FORWARD, FFTW_ESTIMATE);
d@0: ...
d@0: fftw_execute(p); /* repeat as needed */
d@0: ...
d@0: fftw_destroy_plan(p);
d@0: fftw_free(in); fftw_free(out);
d@0: }
d@0:
d@0: (When you compile, you must also link with the fftw3
library,
d@0: e.g. -lfftw3 -lm
on Unix systems.)
d@0:
d@0:
First you allocate the input and output arrays. You can allocate them
d@0: in any way that you like, but we recommend using fftw_malloc
,
d@0: which behaves like
d@0: malloc
except that it properly aligns the array when SIMD
d@0: instructions (such as SSE and Altivec) are available (see SIMD alignment and fftw_malloc).
d@0:
d@0: The data is an array of type fftw_complex
, which is by default a
d@0: double[2]
composed of the real (in[i][0]
) and imaginary
d@0: (in[i][1]
) parts of a complex number.
d@0:
d@0: The next step is to create a plan, which is an object
d@0: that contains all the data that FFTW needs to compute the FFT.
d@0: This function creates the plan:
d@0:
d@0:
fftw_plan fftw_plan_dft_1d(int n, fftw_complex *in, fftw_complex *out, d@0: int sign, unsigned flags); d@0:d@0:
d@0: The first argument, n
, is the size of the transform you are
d@0: trying to compute. The size n
can be any positive integer, but
d@0: sizes that are products of small factors are transformed most
d@0: efficiently (although prime sizes still use an O(n log n) algorithm).
d@0:
d@0:
The next two arguments are pointers to the input and output arrays of
d@0: the transform. These pointers can be equal, indicating an
d@0: in-place transform.
d@0:
d@0: The fourth argument, sign
, can be either FFTW_FORWARD
d@0: (-1
) or FFTW_BACKWARD
(+1
),
d@0: and indicates the direction of the transform you are interested in;
d@0: technically, it is the sign of the exponent in the transform.
d@0:
d@0:
The flags
argument is usually either FFTW_MEASURE
or
d@0: FFTW_ESTIMATE
. FFTW_MEASURE
instructs FFTW to run
d@0: and measure the execution time of several FFTs in order to find the
d@0: best way to compute the transform of size n
. This process takes
d@0: some time (usually a few seconds), depending on your machine and on
d@0: the size of the transform. FFTW_ESTIMATE
, on the contrary,
d@0: does not run any computation and just builds a
d@0: reasonable plan that is probably sub-optimal. In short, if your
d@0: program performs many transforms of the same size and initialization
d@0: time is not important, use FFTW_MEASURE
; otherwise use the
d@0: estimate. The data in the in
/out
arrays is
d@0: overwritten during FFTW_MEASURE
planning, so such
d@0: planning should be done before the input is initialized by the
d@0: user.
d@0:
d@0:
Once the plan has been created, you can use it as many times as you
d@0: like for transforms on the specified in
/out
arrays,
d@0: computing the actual transforms via fftw_execute(plan)
:
d@0:
void fftw_execute(const fftw_plan plan); d@0:d@0:
d@0: If you want to transform a different array of the same size, you
d@0: can create a new plan with fftw_plan_dft_1d
and FFTW
d@0: automatically reuses the information from the previous plan, if
d@0: possible. (Alternatively, with the “guru” interface you can apply a
d@0: given plan to a different array, if you are careful.
d@0: See FFTW Reference.)
d@0:
d@0:
When you are done with the plan, you deallocate it by calling
d@0: fftw_destroy_plan(plan)
:
d@0:
void fftw_destroy_plan(fftw_plan plan); d@0:d@0:
Arrays allocated with fftw_malloc
should be deallocated by
d@0: fftw_free
rather than the ordinary free
(or, heaven
d@0: forbid, delete
).
d@0:
d@0: The DFT results are stored in-order in the array out
, with the
d@0: zero-frequency (DC) component in out[0]
.
d@0: If in != out
, the transform is out-of-place and the input
d@0: array in
is not modified. Otherwise, the input array is
d@0: overwritten with the transform.
d@0:
d@0:
Users should note that FFTW computes an unnormalized DFT.
d@0: Thus, computing a forward followed by a backward transform (or vice
d@0: versa) results in the original array scaled by n
. For the
d@0: definition of the DFT, see What FFTW Really Computes.
d@0:
d@0: If you have a C compiler, such as gcc
, that supports the
d@0: recent C99 standard, and you #include <complex.h>
before
d@0: <fftw3.h>
, then fftw_complex
is the native
d@0: double-precision complex type and you can manipulate it with ordinary
d@0: arithmetic. Otherwise, FFTW defines its own complex type, which is
d@0: bit-compatible with the C99 complex type. See Complex numbers.
d@0: (The C++ <complex>
template class may also be usable via a
d@0: typecast.)
d@0:
d@0: Single and long-double precision versions of FFTW may be installed; to
d@0: use them, replace the fftw_
prefix by fftwf_
or
d@0: fftwl_
and link with -lfftw3f
or -lfftw3l
, but
d@0: use the same <fftw3.h>
header file.
d@0:
d@0: Many more flags exist besides FFTW_MEASURE
and
d@0: FFTW_ESTIMATE
. For example, use FFTW_PATIENT
if you're
d@0: willing to wait even longer for a possibly even faster plan (see FFTW Reference).
d@0: You can also save plans for future use, as described by Words of Wisdom-Saving Plans.
d@0:
d@0:
d@0:
d@0: