cannam@95: cannam@95: cannam@95: Complex One-Dimensional DFTs - FFTW 3.3.3 cannam@95: cannam@95: cannam@95: cannam@95: cannam@95: cannam@95: cannam@95: cannam@95: cannam@95: cannam@95: cannam@95: cannam@95: cannam@95: cannam@95:
cannam@95: cannam@95: cannam@95:

cannam@95: Next: , cannam@95: Previous: Tutorial, cannam@95: Up: Tutorial cannam@95:


cannam@95:
cannam@95: cannam@95:

2.1 Complex One-Dimensional DFTs

cannam@95: cannam@95:
cannam@95: Plan: To bother about the best method of accomplishing an accidental result. cannam@95: [Ambrose Bierce, The Enlarged Devil's Dictionary.] cannam@95:
cannam@95: cannam@95:

The basic usage of FFTW to compute a one-dimensional DFT of size cannam@95: N is simple, and it typically looks something like this code: cannam@95: cannam@95:

     #include <fftw3.h>
cannam@95:      ...
cannam@95:      {
cannam@95:          fftw_complex *in, *out;
cannam@95:          fftw_plan p;
cannam@95:          ...
cannam@95:          in = (fftw_complex*) fftw_malloc(sizeof(fftw_complex) * N);
cannam@95:          out = (fftw_complex*) fftw_malloc(sizeof(fftw_complex) * N);
cannam@95:          p = fftw_plan_dft_1d(N, in, out, FFTW_FORWARD, FFTW_ESTIMATE);
cannam@95:          ...
cannam@95:          fftw_execute(p); /* repeat as needed */
cannam@95:          ...
cannam@95:          fftw_destroy_plan(p);
cannam@95:          fftw_free(in); fftw_free(out);
cannam@95:      }
cannam@95: 
cannam@95:

You must link this code with the fftw3 library. On Unix systems, cannam@95: link with -lfftw3 -lm. cannam@95: cannam@95:

The example code first allocates the input and output arrays. You can cannam@95: allocate them in any way that you like, but we recommend using cannam@95: fftw_malloc, which behaves like cannam@95: malloc except that it properly aligns the array when SIMD cannam@95: instructions (such as SSE and Altivec) are available (see SIMD alignment and fftw_malloc). [Alternatively, we provide a convenient wrapper function fftw_alloc_complex(N) which has the same effect.] cannam@95: cannam@95: cannam@95:

The data is an array of type fftw_complex, which is by default a cannam@95: double[2] composed of the real (in[i][0]) and imaginary cannam@95: (in[i][1]) parts of a complex number. cannam@95: cannam@95: The next step is to create a plan, which is an object cannam@95: that contains all the data that FFTW needs to compute the FFT. cannam@95: This function creates the plan: cannam@95: cannam@95:

     fftw_plan fftw_plan_dft_1d(int n, fftw_complex *in, fftw_complex *out,
cannam@95:                                 int sign, unsigned flags);
cannam@95: 
cannam@95:

cannam@95: The first argument, n, is the size of the transform you are cannam@95: trying to compute. The size n can be any positive integer, but cannam@95: sizes that are products of small factors are transformed most cannam@95: efficiently (although prime sizes still use an O(n log n) algorithm). cannam@95: cannam@95:

The next two arguments are pointers to the input and output arrays of cannam@95: the transform. These pointers can be equal, indicating an cannam@95: in-place transform. cannam@95: cannam@95: cannam@95:

The fourth argument, sign, can be either FFTW_FORWARD cannam@95: (-1) or FFTW_BACKWARD (+1), cannam@95: and indicates the direction of the transform you are interested in; cannam@95: technically, it is the sign of the exponent in the transform. cannam@95: cannam@95:

The flags argument is usually either FFTW_MEASURE or cannam@95: FFTW_ESTIMATE. FFTW_MEASURE instructs FFTW to run cannam@95: and measure the execution time of several FFTs in order to find the cannam@95: best way to compute the transform of size n. This process takes cannam@95: some time (usually a few seconds), depending on your machine and on cannam@95: the size of the transform. FFTW_ESTIMATE, on the contrary, cannam@95: does not run any computation and just builds a cannam@95: reasonable plan that is probably sub-optimal. In short, if your cannam@95: program performs many transforms of the same size and initialization cannam@95: time is not important, use FFTW_MEASURE; otherwise use the cannam@95: estimate. cannam@95: cannam@95:

You must create the plan before initializing the input, because cannam@95: FFTW_MEASURE overwrites the in/out arrays. cannam@95: (Technically, FFTW_ESTIMATE does not touch your arrays, but you cannam@95: should always create plans first just to be sure.) cannam@95: cannam@95:

Once the plan has been created, you can use it as many times as you cannam@95: like for transforms on the specified in/out arrays, cannam@95: computing the actual transforms via fftw_execute(plan): cannam@95:

     void fftw_execute(const fftw_plan plan);
cannam@95: 
cannam@95:

cannam@95: The DFT results are stored in-order in the array out, with the cannam@95: zero-frequency (DC) component in out[0]. cannam@95: If in != out, the transform is out-of-place and the input cannam@95: array in is not modified. Otherwise, the input array is cannam@95: overwritten with the transform. cannam@95: cannam@95:

If you want to transform a different array of the same size, you cannam@95: can create a new plan with fftw_plan_dft_1d and FFTW cannam@95: automatically reuses the information from the previous plan, if cannam@95: possible. Alternatively, with the “guru” interface you can apply a cannam@95: given plan to a different array, if you are careful. cannam@95: See FFTW Reference. cannam@95: cannam@95:

When you are done with the plan, you deallocate it by calling cannam@95: fftw_destroy_plan(plan): cannam@95:

     void fftw_destroy_plan(fftw_plan plan);
cannam@95: 
cannam@95:

If you allocate an array with fftw_malloc() you must deallocate cannam@95: it with fftw_free(). Do not use free() or, heaven cannam@95: forbid, delete. cannam@95: cannam@95: FFTW computes an unnormalized DFT. Thus, computing a forward cannam@95: followed by a backward transform (or vice versa) results in the original cannam@95: array scaled by n. For the definition of the DFT, see What FFTW Really Computes. cannam@95: cannam@95: cannam@95:

If you have a C compiler, such as gcc, that supports the cannam@95: C99 standard, and you #include <complex.h> before cannam@95: <fftw3.h>, then fftw_complex is the native cannam@95: double-precision complex type and you can manipulate it with ordinary cannam@95: arithmetic. Otherwise, FFTW defines its own complex type, which is cannam@95: bit-compatible with the C99 complex type. See Complex numbers. cannam@95: (The C++ <complex> template class may also be usable via a cannam@95: typecast.) cannam@95: cannam@95: To use single or long-double precision versions of FFTW, replace the cannam@95: fftw_ prefix by fftwf_ or fftwl_ and link with cannam@95: -lfftw3f or -lfftw3l, but use the same cannam@95: <fftw3.h> header file. cannam@95: cannam@95: cannam@95:

Many more flags exist besides FFTW_MEASURE and cannam@95: FFTW_ESTIMATE. For example, use FFTW_PATIENT if you're cannam@95: willing to wait even longer for a possibly even faster plan (see FFTW Reference). cannam@95: You can also save plans for future use, as described by Words of Wisdom-Saving Plans. cannam@95: cannam@95: cannam@95: cannam@95: