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