Chris@10: Chris@10: Chris@10: One-Dimensional DFTs of Real Data - FFTW 3.3.3 Chris@10: Chris@10: Chris@10: Chris@10: Chris@10: Chris@10: Chris@10: Chris@10: Chris@10: Chris@10: Chris@10: Chris@10: Chris@10: Chris@10:
Chris@10: Chris@10: Chris@10:

Chris@10: Next: , Chris@10: Previous: Complex Multi-Dimensional DFTs, Chris@10: Up: Tutorial Chris@10:


Chris@10:
Chris@10: Chris@10:

2.3 One-Dimensional DFTs of Real Data

Chris@10: Chris@10:

In many practical applications, the input data in[i] are purely Chris@10: real numbers, in which case the DFT output satisfies the “Hermitian” Chris@10: redundancy: out[i] is the conjugate of out[n-i]. It is Chris@10: possible to take advantage of these circumstances in order to achieve Chris@10: roughly a factor of two improvement in both speed and memory usage. Chris@10: Chris@10:

In exchange for these speed and space advantages, the user sacrifices Chris@10: some of the simplicity of FFTW's complex transforms. First of all, the Chris@10: input and output arrays are of different sizes and types: the Chris@10: input is n real numbers, while the output is n/2+1 Chris@10: complex numbers (the non-redundant outputs); this also requires slight Chris@10: “padding” of the input array for Chris@10: in-place transforms. Second, the inverse transform (complex to real) Chris@10: has the side-effect of overwriting its input array, by default. Chris@10: Neither of these inconveniences should pose a serious problem for Chris@10: users, but it is important to be aware of them. Chris@10: Chris@10:

The routines to perform real-data transforms are almost the same as Chris@10: those for complex transforms: you allocate arrays of double Chris@10: and/or fftw_complex (preferably using fftw_malloc or Chris@10: fftw_alloc_complex), create an fftw_plan, execute it as Chris@10: many times as you want with fftw_execute(plan), and clean up Chris@10: with fftw_destroy_plan(plan) (and fftw_free). The only Chris@10: differences are that the input (or output) is of type double Chris@10: and there are new routines to create the plan. In one dimension: Chris@10: Chris@10:

     fftw_plan fftw_plan_dft_r2c_1d(int n, double *in, fftw_complex *out,
Chris@10:                                     unsigned flags);
Chris@10:      fftw_plan fftw_plan_dft_c2r_1d(int n, fftw_complex *in, double *out,
Chris@10:                                     unsigned flags);
Chris@10: 
Chris@10:

Chris@10: for the real input to complex-Hermitian output (r2c) and Chris@10: complex-Hermitian input to real output (c2r) transforms. Chris@10: Unlike the complex DFT planner, there is no sign argument. Chris@10: Instead, r2c DFTs are always FFTW_FORWARD and c2r DFTs are Chris@10: always FFTW_BACKWARD. Chris@10: (For single/long-double precision Chris@10: fftwf and fftwl, double should be replaced by Chris@10: float and long double, respectively.) Chris@10: Chris@10: Chris@10:

Here, n is the “logical” size of the DFT, not necessarily the Chris@10: physical size of the array. In particular, the real (double) Chris@10: array has n elements, while the complex (fftw_complex) Chris@10: array has n/2+1 elements (where the division is rounded down). Chris@10: For an in-place transform, Chris@10: in and out are aliased to the same array, which must be Chris@10: big enough to hold both; so, the real array would actually have Chris@10: 2*(n/2+1) elements, where the elements beyond the first Chris@10: n are unused padding. (Note that this is very different from Chris@10: the concept of “zero-padding” a transform to a larger length, which Chris@10: changes the logical size of the DFT by actually adding new input Chris@10: data.) The kth element of the complex array is exactly the Chris@10: same as the kth element of the corresponding complex DFT. All Chris@10: positive n are supported; products of small factors are most Chris@10: efficient, but an O(n log n) algorithm is used even for prime sizes. Chris@10: Chris@10:

As noted above, the c2r transform destroys its input array even for Chris@10: out-of-place transforms. This can be prevented, if necessary, by Chris@10: including FFTW_PRESERVE_INPUT in the flags, with Chris@10: unfortunately some sacrifice in performance. Chris@10: This flag is also not currently supported for multi-dimensional real Chris@10: DFTs (next section). Chris@10: Chris@10:

Readers familiar with DFTs of real data will recall that the 0th (the Chris@10: “DC”) and n/2-th (the “Nyquist” frequency, when n is Chris@10: even) elements of the complex output are purely real. Some Chris@10: implementations therefore store the Nyquist element where the DC Chris@10: imaginary part would go, in order to make the input and output arrays Chris@10: the same size. Such packing, however, does not generalize well to Chris@10: multi-dimensional transforms, and the space savings are miniscule in Chris@10: any case; FFTW does not support it. Chris@10: Chris@10:

An alternative interface for one-dimensional r2c and c2r DFTs can be Chris@10: found in the ‘r2r’ interface (see The Halfcomplex-format DFT), with “halfcomplex”-format output that is the same size Chris@10: (and type) as the input array. Chris@10: That interface, although it is not very useful for multi-dimensional Chris@10: transforms, may sometimes yield better performance. Chris@10: Chris@10: Chris@10: Chris@10: