Previous: Multi-Dimensional DFTs of Real Data, Up: Tutorial [Contents][Index]
• The Halfcomplex-format DFT: | ||
• Real even/odd DFTs (cosine/sine transforms): | ||
• The Discrete Hartley Transform: |
FFTW supports several other transform types via a unified r2r
(real-to-real) interface,
so called because it takes a real (double
) array and outputs a
real array of the same size. These r2r transforms currently fall into
three categories: DFTs of real input and complex-Hermitian output in
halfcomplex format, DFTs of real input with even/odd symmetry
(a.k.a. discrete cosine/sine transforms, DCTs/DSTs), and discrete
Hartley transforms (DHTs), all described in more detail by the
following sections.
The r2r transforms follow the by now familiar interface of creating an
fftw_plan
, executing it with fftw_execute(plan)
, and
destroying it with fftw_destroy_plan(plan)
. Furthermore, all
r2r transforms share the same planner interface:
fftw_plan fftw_plan_r2r_1d(int n, double *in, double *out, fftw_r2r_kind kind, unsigned flags); fftw_plan fftw_plan_r2r_2d(int n0, int n1, double *in, double *out, fftw_r2r_kind kind0, fftw_r2r_kind kind1, unsigned flags); fftw_plan fftw_plan_r2r_3d(int n0, int n1, int n2, double *in, double *out, fftw_r2r_kind kind0, fftw_r2r_kind kind1, fftw_r2r_kind kind2, unsigned flags); fftw_plan fftw_plan_r2r(int rank, const int *n, double *in, double *out, const fftw_r2r_kind *kind, unsigned flags);
Just as for the complex DFT, these plan 1d/2d/3d/multi-dimensional
transforms for contiguous arrays in row-major order, transforming (real)
input to output of the same size, where n
specifies the
physical dimensions of the arrays. All positive n
are
supported (with the exception of n=1
for the FFTW_REDFT00
kind, noted in the real-even subsection below); products of small
factors are most efficient (factorizing n-1
and n+1
for
FFTW_REDFT00
and FFTW_RODFT00
kinds, described below), but
an O(n log n) algorithm is used even for prime sizes.
Each dimension has a kind parameter, of type
fftw_r2r_kind
, specifying the kind of r2r transform to be used
for that dimension.
(In the case of fftw_plan_r2r
, this is an array kind[rank]
where kind[i]
is the transform kind for the dimension
n[i]
.) The kind can be one of a set of predefined constants,
defined in the following subsections.
In other words, FFTW computes the separable product of the specified r2r transforms over each dimension, which can be used e.g. for partial differential equations with mixed boundary conditions. (For some r2r kinds, notably the halfcomplex DFT and the DHT, such a separable product is somewhat problematic in more than one dimension, however, as is described below.)
In the current version of FFTW, all r2r transforms except for the halfcomplex type are computed via pre- or post-processing of halfcomplex transforms, and they are therefore not as fast as they could be. Since most other general DCT/DST codes employ a similar algorithm, however, FFTW’s implementation should provide at least competitive performance.
Previous: Multi-Dimensional DFTs of Real Data, Up: Tutorial [Contents][Index]