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