Chris@42: Chris@42: Chris@42: Chris@42: Chris@42:
Chris@42:Chris@42: Next: Multi-Dimensional DFTs of Real Data, Previous: Complex Multi-Dimensional DFTs, Up: Tutorial [Contents][Index]
Chris@42:In many practical applications, the input data in[i] are purely
Chris@42: real numbers, in which case the DFT output satisfies the “Hermitian”
Chris@42: 
Chris@42: redundancy: out[i] is the conjugate of out[n-i].  It is
Chris@42: possible to take advantage of these circumstances in order to achieve
Chris@42: roughly a factor of two improvement in both speed and memory usage.
Chris@42: 
In exchange for these speed and space advantages, the user sacrifices
Chris@42: some of the simplicity of FFTW’s complex transforms. First of all, the
Chris@42: input and output arrays are of different sizes and types: the
Chris@42: input is n real numbers, while the output is n/2+1
Chris@42: complex numbers (the non-redundant outputs); this also requires slight
Chris@42: “padding” of the input array for
Chris@42: 
Chris@42: in-place transforms.  Second, the inverse transform (complex to real)
Chris@42: has the side-effect of overwriting its input array, by default.
Chris@42: Neither of these inconveniences should pose a serious problem for
Chris@42: users, but it is important to be aware of them.
Chris@42: 
The routines to perform real-data transforms are almost the same as
Chris@42: those for complex transforms: you allocate arrays of double
Chris@42: and/or fftw_complex (preferably using fftw_malloc or
Chris@42: fftw_alloc_complex), create an fftw_plan, execute it as
Chris@42: many times as you want with fftw_execute(plan), and clean up
Chris@42: with fftw_destroy_plan(plan) (and fftw_free).  The only
Chris@42: differences are that the input (or output) is of type double
Chris@42: and there are new routines to create the plan.  In one dimension:
Chris@42: 
fftw_plan fftw_plan_dft_r2c_1d(int n, double *in, fftw_complex *out, Chris@42: unsigned flags); Chris@42: fftw_plan fftw_plan_dft_c2r_1d(int n, fftw_complex *in, double *out, Chris@42: unsigned flags); Chris@42:
for the real input to complex-Hermitian output (r2c) and
Chris@42: complex-Hermitian input to real output (c2r) transforms.
Chris@42: 
Chris@42: 
Chris@42: Unlike the complex DFT planner, there is no sign argument.
Chris@42: Instead, r2c DFTs are always FFTW_FORWARD and c2r DFTs are
Chris@42: always FFTW_BACKWARD.
Chris@42: 
Chris@42: 
Chris@42: (For single/long-double precision
Chris@42: fftwf and fftwl, double should be replaced by
Chris@42: float and long double, respectively.)
Chris@42: 
Chris@42: 
Here, n is the “logical” size of the DFT, not necessarily the
Chris@42: physical size of the array.  In particular, the real (double)
Chris@42: array has n elements, while the complex (fftw_complex)
Chris@42: array has n/2+1 elements (where the division is rounded down).
Chris@42: For an in-place transform,
Chris@42: 
Chris@42: in and out are aliased to the same array, which must be
Chris@42: big enough to hold both; so, the real array would actually have
Chris@42: 2*(n/2+1) elements, where the elements beyond the first
Chris@42: n are unused padding.  (Note that this is very different from
Chris@42: the concept of “zero-padding” a transform to a larger length, which
Chris@42: changes the logical size of the DFT by actually adding new input
Chris@42: data.)  The kth element of the complex array is exactly the
Chris@42: same as the kth element of the corresponding complex DFT.  All
Chris@42: positive n are supported; products of small factors are most
Chris@42: efficient, but an O(n log n) algorithm is used even for prime sizes.
Chris@42: 
As noted above, the c2r transform destroys its input array even for
Chris@42: out-of-place transforms.  This can be prevented, if necessary, by
Chris@42: including FFTW_PRESERVE_INPUT in the flags, with
Chris@42: unfortunately some sacrifice in performance.
Chris@42: 
Chris@42: 
Chris@42: This flag is also not currently supported for multi-dimensional real
Chris@42: DFTs (next section).
Chris@42: 
Readers familiar with DFTs of real data will recall that the 0th (the
Chris@42: “DC”) and n/2-th (the “Nyquist” frequency, when n is
Chris@42: even) elements of the complex output are purely real.  Some
Chris@42: implementations therefore store the Nyquist element where the DC
Chris@42: imaginary part would go, in order to make the input and output arrays
Chris@42: the same size.  Such packing, however, does not generalize well to
Chris@42: multi-dimensional transforms, and the space savings are miniscule in
Chris@42: any case; FFTW does not support it.
Chris@42: 
An alternative interface for one-dimensional r2c and c2r DFTs can be Chris@42: found in the ‘r2r’ interface (see The Halfcomplex-format DFT), with “halfcomplex”-format output that is the same size Chris@42: (and type) as the input array. Chris@42: Chris@42: That interface, although it is not very useful for multi-dimensional Chris@42: transforms, may sometimes yield better performance. Chris@42:
Chris@42:Chris@42: Next: Multi-Dimensional DFTs of Real Data, Previous: Complex Multi-Dimensional DFTs, Up: Tutorial [Contents][Index]
Chris@42: