annotate src/fftw-3.3.3/doc/tutorial.texi @ 169:223a55898ab9 tip default

Add null config files
author Chris Cannam <cannam@all-day-breakfast.com>
date Mon, 02 Mar 2020 14:03:47 +0000
parents 89f5e221ed7b
children
rev   line source
cannam@95 1 @node Tutorial, Other Important Topics, Introduction, Top
cannam@95 2 @chapter Tutorial
cannam@95 3 @menu
cannam@95 4 * Complex One-Dimensional DFTs::
cannam@95 5 * Complex Multi-Dimensional DFTs::
cannam@95 6 * One-Dimensional DFTs of Real Data::
cannam@95 7 * Multi-Dimensional DFTs of Real Data::
cannam@95 8 * More DFTs of Real Data::
cannam@95 9 @end menu
cannam@95 10
cannam@95 11 This chapter describes the basic usage of FFTW, i.e., how to compute
cannam@95 12 @cindex basic interface
cannam@95 13 the Fourier transform of a single array. This chapter tells the
cannam@95 14 truth, but not the @emph{whole} truth. Specifically, FFTW implements
cannam@95 15 additional routines and flags that are not documented here, although
cannam@95 16 in many cases we try to indicate where added capabilities exist. For
cannam@95 17 more complete information, see @ref{FFTW Reference}. (Note that you
cannam@95 18 need to compile and install FFTW before you can use it in a program.
cannam@95 19 For the details of the installation, see @ref{Installation and
cannam@95 20 Customization}.)
cannam@95 21
cannam@95 22 We recommend that you read this tutorial in order.@footnote{You can
cannam@95 23 read the tutorial in bit-reversed order after computing your first
cannam@95 24 transform.} At the least, read the first section (@pxref{Complex
cannam@95 25 One-Dimensional DFTs}) before reading any of the others, even if your
cannam@95 26 main interest lies in one of the other transform types.
cannam@95 27
cannam@95 28 Users of FFTW version 2 and earlier may also want to read @ref{Upgrading
cannam@95 29 from FFTW version 2}.
cannam@95 30
cannam@95 31 @c ------------------------------------------------------------
cannam@95 32 @node Complex One-Dimensional DFTs, Complex Multi-Dimensional DFTs, Tutorial, Tutorial
cannam@95 33 @section Complex One-Dimensional DFTs
cannam@95 34
cannam@95 35 @quotation
cannam@95 36 Plan: To bother about the best method of accomplishing an accidental result.
cannam@95 37 [Ambrose Bierce, @cite{The Enlarged Devil's Dictionary}.]
cannam@95 38 @cindex Devil
cannam@95 39 @end quotation
cannam@95 40
cannam@95 41 @iftex
cannam@95 42 @medskip
cannam@95 43 @end iftex
cannam@95 44
cannam@95 45 The basic usage of FFTW to compute a one-dimensional DFT of size
cannam@95 46 @code{N} is simple, and it typically looks something like this code:
cannam@95 47
cannam@95 48 @example
cannam@95 49 #include <fftw3.h>
cannam@95 50 ...
cannam@95 51 @{
cannam@95 52 fftw_complex *in, *out;
cannam@95 53 fftw_plan p;
cannam@95 54 ...
cannam@95 55 in = (fftw_complex*) fftw_malloc(sizeof(fftw_complex) * N);
cannam@95 56 out = (fftw_complex*) fftw_malloc(sizeof(fftw_complex) * N);
cannam@95 57 p = fftw_plan_dft_1d(N, in, out, FFTW_FORWARD, FFTW_ESTIMATE);
cannam@95 58 ...
cannam@95 59 fftw_execute(p); /* @r{repeat as needed} */
cannam@95 60 ...
cannam@95 61 fftw_destroy_plan(p);
cannam@95 62 fftw_free(in); fftw_free(out);
cannam@95 63 @}
cannam@95 64 @end example
cannam@95 65
cannam@95 66 You must link this code with the @code{fftw3} library. On Unix systems,
cannam@95 67 link with @code{-lfftw3 -lm}.
cannam@95 68
cannam@95 69 The example code first allocates the input and output arrays. You can
cannam@95 70 allocate them in any way that you like, but we recommend using
cannam@95 71 @code{fftw_malloc}, which behaves like
cannam@95 72 @findex fftw_malloc
cannam@95 73 @code{malloc} except that it properly aligns the array when SIMD
cannam@95 74 instructions (such as SSE and Altivec) are available (@pxref{SIMD
cannam@95 75 alignment and fftw_malloc}). [Alternatively, we provide a convenient wrapper function @code{fftw_alloc_complex(N)} which has the same effect.]
cannam@95 76 @findex fftw_alloc_complex
cannam@95 77 @cindex SIMD
cannam@95 78
cannam@95 79
cannam@95 80 The data is an array of type @code{fftw_complex}, which is by default a
cannam@95 81 @code{double[2]} composed of the real (@code{in[i][0]}) and imaginary
cannam@95 82 (@code{in[i][1]}) parts of a complex number.
cannam@95 83 @tindex fftw_complex
cannam@95 84
cannam@95 85 The next step is to create a @dfn{plan}, which is an object
cannam@95 86 @cindex plan
cannam@95 87 that contains all the data that FFTW needs to compute the FFT.
cannam@95 88 This function creates the plan:
cannam@95 89
cannam@95 90 @example
cannam@95 91 fftw_plan fftw_plan_dft_1d(int n, fftw_complex *in, fftw_complex *out,
cannam@95 92 int sign, unsigned flags);
cannam@95 93 @end example
cannam@95 94 @findex fftw_plan_dft_1d
cannam@95 95 @tindex fftw_plan
cannam@95 96
cannam@95 97 The first argument, @code{n}, is the size of the transform you are
cannam@95 98 trying to compute. The size @code{n} can be any positive integer, but
cannam@95 99 sizes that are products of small factors are transformed most
cannam@95 100 efficiently (although prime sizes still use an @Onlogn{} algorithm).
cannam@95 101
cannam@95 102 The next two arguments are pointers to the input and output arrays of
cannam@95 103 the transform. These pointers can be equal, indicating an
cannam@95 104 @dfn{in-place} transform.
cannam@95 105 @cindex in-place
cannam@95 106
cannam@95 107
cannam@95 108 The fourth argument, @code{sign}, can be either @code{FFTW_FORWARD}
cannam@95 109 (@code{-1}) or @code{FFTW_BACKWARD} (@code{+1}),
cannam@95 110 @ctindex FFTW_FORWARD
cannam@95 111 @ctindex FFTW_BACKWARD
cannam@95 112 and indicates the direction of the transform you are interested in;
cannam@95 113 technically, it is the sign of the exponent in the transform.
cannam@95 114
cannam@95 115 The @code{flags} argument is usually either @code{FFTW_MEASURE} or
cannam@95 116 @cindex flags
cannam@95 117 @code{FFTW_ESTIMATE}. @code{FFTW_MEASURE} instructs FFTW to run
cannam@95 118 @ctindex FFTW_MEASURE
cannam@95 119 and measure the execution time of several FFTs in order to find the
cannam@95 120 best way to compute the transform of size @code{n}. This process takes
cannam@95 121 some time (usually a few seconds), depending on your machine and on
cannam@95 122 the size of the transform. @code{FFTW_ESTIMATE}, on the contrary,
cannam@95 123 does not run any computation and just builds a
cannam@95 124 @ctindex FFTW_ESTIMATE
cannam@95 125 reasonable plan that is probably sub-optimal. In short, if your
cannam@95 126 program performs many transforms of the same size and initialization
cannam@95 127 time is not important, use @code{FFTW_MEASURE}; otherwise use the
cannam@95 128 estimate.
cannam@95 129
cannam@95 130 @emph{You must create the plan before initializing the input}, because
cannam@95 131 @code{FFTW_MEASURE} overwrites the @code{in}/@code{out} arrays.
cannam@95 132 (Technically, @code{FFTW_ESTIMATE} does not touch your arrays, but you
cannam@95 133 should always create plans first just to be sure.)
cannam@95 134
cannam@95 135 Once the plan has been created, you can use it as many times as you
cannam@95 136 like for transforms on the specified @code{in}/@code{out} arrays,
cannam@95 137 computing the actual transforms via @code{fftw_execute(plan)}:
cannam@95 138 @example
cannam@95 139 void fftw_execute(const fftw_plan plan);
cannam@95 140 @end example
cannam@95 141 @findex fftw_execute
cannam@95 142
cannam@95 143 The DFT results are stored in-order in the array @code{out}, with the
cannam@95 144 zero-frequency (DC) component in @code{out[0]}.
cannam@95 145 @cindex frequency
cannam@95 146 If @code{in != out}, the transform is @dfn{out-of-place} and the input
cannam@95 147 array @code{in} is not modified. Otherwise, the input array is
cannam@95 148 overwritten with the transform.
cannam@95 149
cannam@95 150 @cindex execute
cannam@95 151 If you want to transform a @emph{different} array of the same size, you
cannam@95 152 can create a new plan with @code{fftw_plan_dft_1d} and FFTW
cannam@95 153 automatically reuses the information from the previous plan, if
cannam@95 154 possible. Alternatively, with the ``guru'' interface you can apply a
cannam@95 155 given plan to a different array, if you are careful.
cannam@95 156 @xref{FFTW Reference}.
cannam@95 157
cannam@95 158 When you are done with the plan, you deallocate it by calling
cannam@95 159 @code{fftw_destroy_plan(plan)}:
cannam@95 160 @example
cannam@95 161 void fftw_destroy_plan(fftw_plan plan);
cannam@95 162 @end example
cannam@95 163 @findex fftw_destroy_plan
cannam@95 164 If you allocate an array with @code{fftw_malloc()} you must deallocate
cannam@95 165 it with @code{fftw_free()}. Do not use @code{free()} or, heaven
cannam@95 166 forbid, @code{delete}.
cannam@95 167 @findex fftw_free
cannam@95 168
cannam@95 169 FFTW computes an @emph{unnormalized} DFT. Thus, computing a forward
cannam@95 170 followed by a backward transform (or vice versa) results in the original
cannam@95 171 array scaled by @code{n}. For the definition of the DFT, see @ref{What
cannam@95 172 FFTW Really Computes}.
cannam@95 173 @cindex DFT
cannam@95 174 @cindex normalization
cannam@95 175
cannam@95 176
cannam@95 177 If you have a C compiler, such as @code{gcc}, that supports the
cannam@95 178 C99 standard, and you @code{#include <complex.h>} @emph{before}
cannam@95 179 @code{<fftw3.h>}, then @code{fftw_complex} is the native
cannam@95 180 double-precision complex type and you can manipulate it with ordinary
cannam@95 181 arithmetic. Otherwise, FFTW defines its own complex type, which is
cannam@95 182 bit-compatible with the C99 complex type. @xref{Complex numbers}.
cannam@95 183 (The C++ @code{<complex>} template class may also be usable via a
cannam@95 184 typecast.)
cannam@95 185 @cindex C++
cannam@95 186
cannam@95 187 To use single or long-double precision versions of FFTW, replace the
cannam@95 188 @code{fftw_} prefix by @code{fftwf_} or @code{fftwl_} and link with
cannam@95 189 @code{-lfftw3f} or @code{-lfftw3l}, but use the @emph{same}
cannam@95 190 @code{<fftw3.h>} header file.
cannam@95 191 @cindex precision
cannam@95 192
cannam@95 193
cannam@95 194 Many more flags exist besides @code{FFTW_MEASURE} and
cannam@95 195 @code{FFTW_ESTIMATE}. For example, use @code{FFTW_PATIENT} if you're
cannam@95 196 willing to wait even longer for a possibly even faster plan (@pxref{FFTW
cannam@95 197 Reference}).
cannam@95 198 @ctindex FFTW_PATIENT
cannam@95 199 You can also save plans for future use, as described by @ref{Words of
cannam@95 200 Wisdom-Saving Plans}.
cannam@95 201
cannam@95 202 @c ------------------------------------------------------------
cannam@95 203 @node Complex Multi-Dimensional DFTs, One-Dimensional DFTs of Real Data, Complex One-Dimensional DFTs, Tutorial
cannam@95 204 @section Complex Multi-Dimensional DFTs
cannam@95 205
cannam@95 206 Multi-dimensional transforms work much the same way as one-dimensional
cannam@95 207 transforms: you allocate arrays of @code{fftw_complex} (preferably
cannam@95 208 using @code{fftw_malloc}), create an @code{fftw_plan}, execute it as
cannam@95 209 many times as you want with @code{fftw_execute(plan)}, and clean up
cannam@95 210 with @code{fftw_destroy_plan(plan)} (and @code{fftw_free}).
cannam@95 211
cannam@95 212 FFTW provides two routines for creating plans for 2d and 3d transforms,
cannam@95 213 and one routine for creating plans of arbitrary dimensionality.
cannam@95 214 The 2d and 3d routines have the following signature:
cannam@95 215 @example
cannam@95 216 fftw_plan fftw_plan_dft_2d(int n0, int n1,
cannam@95 217 fftw_complex *in, fftw_complex *out,
cannam@95 218 int sign, unsigned flags);
cannam@95 219 fftw_plan fftw_plan_dft_3d(int n0, int n1, int n2,
cannam@95 220 fftw_complex *in, fftw_complex *out,
cannam@95 221 int sign, unsigned flags);
cannam@95 222 @end example
cannam@95 223 @findex fftw_plan_dft_2d
cannam@95 224 @findex fftw_plan_dft_3d
cannam@95 225
cannam@95 226 These routines create plans for @code{n0} by @code{n1} two-dimensional
cannam@95 227 (2d) transforms and @code{n0} by @code{n1} by @code{n2} 3d transforms,
cannam@95 228 respectively. All of these transforms operate on contiguous arrays in
cannam@95 229 the C-standard @dfn{row-major} order, so that the last dimension has the
cannam@95 230 fastest-varying index in the array. This layout is described further in
cannam@95 231 @ref{Multi-dimensional Array Format}.
cannam@95 232
cannam@95 233 FFTW can also compute transforms of higher dimensionality. In order to
cannam@95 234 avoid confusion between the various meanings of the the word
cannam@95 235 ``dimension'', we use the term @emph{rank}
cannam@95 236 @cindex rank
cannam@95 237 to denote the number of independent indices in an array.@footnote{The
cannam@95 238 term ``rank'' is commonly used in the APL, FORTRAN, and Common Lisp
cannam@95 239 traditions, although it is not so common in the C@tie{}world.} For
cannam@95 240 example, we say that a 2d transform has rank@tie{}2, a 3d transform has
cannam@95 241 rank@tie{}3, and so on. You can plan transforms of arbitrary rank by
cannam@95 242 means of the following function:
cannam@95 243
cannam@95 244 @example
cannam@95 245 fftw_plan fftw_plan_dft(int rank, const int *n,
cannam@95 246 fftw_complex *in, fftw_complex *out,
cannam@95 247 int sign, unsigned flags);
cannam@95 248 @end example
cannam@95 249 @findex fftw_plan_dft
cannam@95 250
cannam@95 251 Here, @code{n} is a pointer to an array @code{n[rank]} denoting an
cannam@95 252 @code{n[0]} by @code{n[1]} by @dots{} by @code{n[rank-1]} transform.
cannam@95 253 Thus, for example, the call
cannam@95 254 @example
cannam@95 255 fftw_plan_dft_2d(n0, n1, in, out, sign, flags);
cannam@95 256 @end example
cannam@95 257 is equivalent to the following code fragment:
cannam@95 258 @example
cannam@95 259 int n[2];
cannam@95 260 n[0] = n0;
cannam@95 261 n[1] = n1;
cannam@95 262 fftw_plan_dft(2, n, in, out, sign, flags);
cannam@95 263 @end example
cannam@95 264 @code{fftw_plan_dft} is not restricted to 2d and 3d transforms,
cannam@95 265 however, but it can plan transforms of arbitrary rank.
cannam@95 266
cannam@95 267 You may have noticed that all the planner routines described so far
cannam@95 268 have overlapping functionality. For example, you can plan a 1d or 2d
cannam@95 269 transform by using @code{fftw_plan_dft} with a @code{rank} of @code{1}
cannam@95 270 or @code{2}, or even by calling @code{fftw_plan_dft_3d} with @code{n0}
cannam@95 271 and/or @code{n1} equal to @code{1} (with no loss in efficiency). This
cannam@95 272 pattern continues, and FFTW's planning routines in general form a
cannam@95 273 ``partial order,'' sequences of
cannam@95 274 @cindex partial order
cannam@95 275 interfaces with strictly increasing generality but correspondingly
cannam@95 276 greater complexity.
cannam@95 277
cannam@95 278 @code{fftw_plan_dft} is the most general complex-DFT routine that we
cannam@95 279 describe in this tutorial, but there are also the advanced and guru interfaces,
cannam@95 280 @cindex advanced interface
cannam@95 281 @cindex guru interface
cannam@95 282 which allow one to efficiently combine multiple/strided transforms
cannam@95 283 into a single FFTW plan, transform a subset of a larger
cannam@95 284 multi-dimensional array, and/or to handle more general complex-number
cannam@95 285 formats. For more information, see @ref{FFTW Reference}.
cannam@95 286
cannam@95 287 @c ------------------------------------------------------------
cannam@95 288 @node One-Dimensional DFTs of Real Data, Multi-Dimensional DFTs of Real Data, Complex Multi-Dimensional DFTs, Tutorial
cannam@95 289 @section One-Dimensional DFTs of Real Data
cannam@95 290
cannam@95 291 In many practical applications, the input data @code{in[i]} are purely
cannam@95 292 real numbers, in which case the DFT output satisfies the ``Hermitian''
cannam@95 293 @cindex Hermitian
cannam@95 294 redundancy: @code{out[i]} is the conjugate of @code{out[n-i]}. It is
cannam@95 295 possible to take advantage of these circumstances in order to achieve
cannam@95 296 roughly a factor of two improvement in both speed and memory usage.
cannam@95 297
cannam@95 298 In exchange for these speed and space advantages, the user sacrifices
cannam@95 299 some of the simplicity of FFTW's complex transforms. First of all, the
cannam@95 300 input and output arrays are of @emph{different sizes and types}: the
cannam@95 301 input is @code{n} real numbers, while the output is @code{n/2+1}
cannam@95 302 complex numbers (the non-redundant outputs); this also requires slight
cannam@95 303 ``padding'' of the input array for
cannam@95 304 @cindex padding
cannam@95 305 in-place transforms. Second, the inverse transform (complex to real)
cannam@95 306 has the side-effect of @emph{overwriting its input array}, by default.
cannam@95 307 Neither of these inconveniences should pose a serious problem for
cannam@95 308 users, but it is important to be aware of them.
cannam@95 309
cannam@95 310 The routines to perform real-data transforms are almost the same as
cannam@95 311 those for complex transforms: you allocate arrays of @code{double}
cannam@95 312 and/or @code{fftw_complex} (preferably using @code{fftw_malloc} or
cannam@95 313 @code{fftw_alloc_complex}), create an @code{fftw_plan}, execute it as
cannam@95 314 many times as you want with @code{fftw_execute(plan)}, and clean up
cannam@95 315 with @code{fftw_destroy_plan(plan)} (and @code{fftw_free}). The only
cannam@95 316 differences are that the input (or output) is of type @code{double}
cannam@95 317 and there are new routines to create the plan. In one dimension:
cannam@95 318
cannam@95 319 @example
cannam@95 320 fftw_plan fftw_plan_dft_r2c_1d(int n, double *in, fftw_complex *out,
cannam@95 321 unsigned flags);
cannam@95 322 fftw_plan fftw_plan_dft_c2r_1d(int n, fftw_complex *in, double *out,
cannam@95 323 unsigned flags);
cannam@95 324 @end example
cannam@95 325 @findex fftw_plan_dft_r2c_1d
cannam@95 326 @findex fftw_plan_dft_c2r_1d
cannam@95 327
cannam@95 328 for the real input to complex-Hermitian output (@dfn{r2c}) and
cannam@95 329 complex-Hermitian input to real output (@dfn{c2r}) transforms.
cannam@95 330 @cindex r2c
cannam@95 331 @cindex c2r
cannam@95 332 Unlike the complex DFT planner, there is no @code{sign} argument.
cannam@95 333 Instead, r2c DFTs are always @code{FFTW_FORWARD} and c2r DFTs are
cannam@95 334 always @code{FFTW_BACKWARD}.
cannam@95 335 @ctindex FFTW_FORWARD
cannam@95 336 @ctindex FFTW_BACKWARD
cannam@95 337 (For single/long-double precision
cannam@95 338 @code{fftwf} and @code{fftwl}, @code{double} should be replaced by
cannam@95 339 @code{float} and @code{long double}, respectively.)
cannam@95 340 @cindex precision
cannam@95 341
cannam@95 342
cannam@95 343 Here, @code{n} is the ``logical'' size of the DFT, not necessarily the
cannam@95 344 physical size of the array. In particular, the real (@code{double})
cannam@95 345 array has @code{n} elements, while the complex (@code{fftw_complex})
cannam@95 346 array has @code{n/2+1} elements (where the division is rounded down).
cannam@95 347 For an in-place transform,
cannam@95 348 @cindex in-place
cannam@95 349 @code{in} and @code{out} are aliased to the same array, which must be
cannam@95 350 big enough to hold both; so, the real array would actually have
cannam@95 351 @code{2*(n/2+1)} elements, where the elements beyond the first
cannam@95 352 @code{n} are unused padding. (Note that this is very different from
cannam@95 353 the concept of ``zero-padding'' a transform to a larger length, which
cannam@95 354 changes the logical size of the DFT by actually adding new input
cannam@95 355 data.) The @math{k}th element of the complex array is exactly the
cannam@95 356 same as the @math{k}th element of the corresponding complex DFT. All
cannam@95 357 positive @code{n} are supported; products of small factors are most
cannam@95 358 efficient, but an @Onlogn algorithm is used even for prime sizes.
cannam@95 359
cannam@95 360 As noted above, the c2r transform destroys its input array even for
cannam@95 361 out-of-place transforms. This can be prevented, if necessary, by
cannam@95 362 including @code{FFTW_PRESERVE_INPUT} in the @code{flags}, with
cannam@95 363 unfortunately some sacrifice in performance.
cannam@95 364 @cindex flags
cannam@95 365 @ctindex FFTW_PRESERVE_INPUT
cannam@95 366 This flag is also not currently supported for multi-dimensional real
cannam@95 367 DFTs (next section).
cannam@95 368
cannam@95 369 Readers familiar with DFTs of real data will recall that the 0th (the
cannam@95 370 ``DC'') and @code{n/2}-th (the ``Nyquist'' frequency, when @code{n} is
cannam@95 371 even) elements of the complex output are purely real. Some
cannam@95 372 implementations therefore store the Nyquist element where the DC
cannam@95 373 imaginary part would go, in order to make the input and output arrays
cannam@95 374 the same size. Such packing, however, does not generalize well to
cannam@95 375 multi-dimensional transforms, and the space savings are miniscule in
cannam@95 376 any case; FFTW does not support it.
cannam@95 377
cannam@95 378 An alternative interface for one-dimensional r2c and c2r DFTs can be
cannam@95 379 found in the @samp{r2r} interface (@pxref{The Halfcomplex-format
cannam@95 380 DFT}), with ``halfcomplex''-format output that @emph{is} the same size
cannam@95 381 (and type) as the input array.
cannam@95 382 @cindex halfcomplex format
cannam@95 383 That interface, although it is not very useful for multi-dimensional
cannam@95 384 transforms, may sometimes yield better performance.
cannam@95 385
cannam@95 386 @c ------------------------------------------------------------
cannam@95 387 @node Multi-Dimensional DFTs of Real Data, More DFTs of Real Data, One-Dimensional DFTs of Real Data, Tutorial
cannam@95 388 @section Multi-Dimensional DFTs of Real Data
cannam@95 389
cannam@95 390 Multi-dimensional DFTs of real data use the following planner routines:
cannam@95 391
cannam@95 392 @example
cannam@95 393 fftw_plan fftw_plan_dft_r2c_2d(int n0, int n1,
cannam@95 394 double *in, fftw_complex *out,
cannam@95 395 unsigned flags);
cannam@95 396 fftw_plan fftw_plan_dft_r2c_3d(int n0, int n1, int n2,
cannam@95 397 double *in, fftw_complex *out,
cannam@95 398 unsigned flags);
cannam@95 399 fftw_plan fftw_plan_dft_r2c(int rank, const int *n,
cannam@95 400 double *in, fftw_complex *out,
cannam@95 401 unsigned flags);
cannam@95 402 @end example
cannam@95 403 @findex fftw_plan_dft_r2c_2d
cannam@95 404 @findex fftw_plan_dft_r2c_3d
cannam@95 405 @findex fftw_plan_dft_r2c
cannam@95 406
cannam@95 407 as well as the corresponding @code{c2r} routines with the input/output
cannam@95 408 types swapped. These routines work similarly to their complex
cannam@95 409 analogues, except for the fact that here the complex output array is cut
cannam@95 410 roughly in half and the real array requires padding for in-place
cannam@95 411 transforms (as in 1d, above).
cannam@95 412
cannam@95 413 As before, @code{n} is the logical size of the array, and the
cannam@95 414 consequences of this on the the format of the complex arrays deserve
cannam@95 415 careful attention.
cannam@95 416 @cindex r2c/c2r multi-dimensional array format
cannam@95 417 Suppose that the real data has dimensions @ndims (in row-major order).
cannam@95 418 Then, after an r2c transform, the output is an @ndimshalf array of
cannam@95 419 @code{fftw_complex} values in row-major order, corresponding to slightly
cannam@95 420 over half of the output of the corresponding complex DFT. (The division
cannam@95 421 is rounded down.) The ordering of the data is otherwise exactly the
cannam@95 422 same as in the complex-DFT case.
cannam@95 423
cannam@95 424 For out-of-place transforms, this is the end of the story: the real
cannam@95 425 data is stored as a row-major array of size @ndims and the complex
cannam@95 426 data is stored as a row-major array of size @ndimshalf{}.
cannam@95 427
cannam@95 428 For in-place transforms, however, extra padding of the real-data array
cannam@95 429 is necessary because the complex array is larger than the real array,
cannam@95 430 and the two arrays share the same memory locations. Thus, for
cannam@95 431 in-place transforms, the final dimension of the real-data array must
cannam@95 432 be padded with extra values to accommodate the size of the complex
cannam@95 433 data---two values if the last dimension is even and one if it is odd.
cannam@95 434 @cindex padding
cannam@95 435 That is, the last dimension of the real data must physically contain
cannam@95 436 @tex
cannam@95 437 $2 (n_{d-1}/2+1)$
cannam@95 438 @end tex
cannam@95 439 @ifinfo
cannam@95 440 2 * (n[d-1]/2+1)
cannam@95 441 @end ifinfo
cannam@95 442 @html
cannam@95 443 2 * (n<sub>d-1</sub>/2+1)
cannam@95 444 @end html
cannam@95 445 @code{double} values (exactly enough to hold the complex data).
cannam@95 446 This physical array size does not, however, change the @emph{logical}
cannam@95 447 array size---only
cannam@95 448 @tex
cannam@95 449 $n_{d-1}$
cannam@95 450 @end tex
cannam@95 451 @ifinfo
cannam@95 452 n[d-1]
cannam@95 453 @end ifinfo
cannam@95 454 @html
cannam@95 455 n<sub>d-1</sub>
cannam@95 456 @end html
cannam@95 457 values are actually stored in the last dimension, and
cannam@95 458 @tex
cannam@95 459 $n_{d-1}$
cannam@95 460 @end tex
cannam@95 461 @ifinfo
cannam@95 462 n[d-1]
cannam@95 463 @end ifinfo
cannam@95 464 @html
cannam@95 465 n<sub>d-1</sub>
cannam@95 466 @end html
cannam@95 467 is the last dimension passed to the plan-creation routine.
cannam@95 468
cannam@95 469 For example, consider the transform of a two-dimensional real array of
cannam@95 470 size @code{n0} by @code{n1}. The output of the r2c transform is a
cannam@95 471 two-dimensional complex array of size @code{n0} by @code{n1/2+1}, where
cannam@95 472 the @code{y} dimension has been cut nearly in half because of
cannam@95 473 redundancies in the output. Because @code{fftw_complex} is twice the
cannam@95 474 size of @code{double}, the output array is slightly bigger than the
cannam@95 475 input array. Thus, if we want to compute the transform in place, we
cannam@95 476 must @emph{pad} the input array so that it is of size @code{n0} by
cannam@95 477 @code{2*(n1/2+1)}. If @code{n1} is even, then there are two padding
cannam@95 478 elements at the end of each row (which need not be initialized, as they
cannam@95 479 are only used for output).
cannam@95 480
cannam@95 481 @ifhtml
cannam@95 482 The following illustration depicts the input and output arrays just
cannam@95 483 described, for both the out-of-place and in-place transforms (with the
cannam@95 484 arrows indicating consecutive memory locations):
cannam@95 485 @image{rfftwnd-for-html}
cannam@95 486 @end ifhtml
cannam@95 487 @ifnotinfo
cannam@95 488 @ifnothtml
cannam@95 489 @float Figure,fig:rfftwnd
cannam@95 490 @center @image{rfftwnd}
cannam@95 491 @caption{Illustration of the data layout for a 2d @code{nx} by @code{ny}
cannam@95 492 real-to-complex transform.}
cannam@95 493 @end float
cannam@95 494 @ref{fig:rfftwnd} depicts the input and output arrays just
cannam@95 495 described, for both the out-of-place and in-place transforms (with the
cannam@95 496 arrows indicating consecutive memory locations):
cannam@95 497 @end ifnothtml
cannam@95 498 @end ifnotinfo
cannam@95 499
cannam@95 500 These transforms are unnormalized, so an r2c followed by a c2r
cannam@95 501 transform (or vice versa) will result in the original data scaled by
cannam@95 502 the number of real data elements---that is, the product of the
cannam@95 503 (logical) dimensions of the real data.
cannam@95 504 @cindex normalization
cannam@95 505
cannam@95 506
cannam@95 507 (Because the last dimension is treated specially, if it is equal to
cannam@95 508 @code{1} the transform is @emph{not} equivalent to a lower-dimensional
cannam@95 509 r2c/c2r transform. In that case, the last complex dimension also has
cannam@95 510 size @code{1} (@code{=1/2+1}), and no advantage is gained over the
cannam@95 511 complex transforms.)
cannam@95 512
cannam@95 513 @c ------------------------------------------------------------
cannam@95 514 @node More DFTs of Real Data, , Multi-Dimensional DFTs of Real Data, Tutorial
cannam@95 515 @section More DFTs of Real Data
cannam@95 516 @menu
cannam@95 517 * The Halfcomplex-format DFT::
cannam@95 518 * Real even/odd DFTs (cosine/sine transforms)::
cannam@95 519 * The Discrete Hartley Transform::
cannam@95 520 @end menu
cannam@95 521
cannam@95 522 FFTW supports several other transform types via a unified @dfn{r2r}
cannam@95 523 (real-to-real) interface,
cannam@95 524 @cindex r2r
cannam@95 525 so called because it takes a real (@code{double}) array and outputs a
cannam@95 526 real array of the same size. These r2r transforms currently fall into
cannam@95 527 three categories: DFTs of real input and complex-Hermitian output in
cannam@95 528 halfcomplex format, DFTs of real input with even/odd symmetry
cannam@95 529 (a.k.a. discrete cosine/sine transforms, DCTs/DSTs), and discrete
cannam@95 530 Hartley transforms (DHTs), all described in more detail by the
cannam@95 531 following sections.
cannam@95 532
cannam@95 533 The r2r transforms follow the by now familiar interface of creating an
cannam@95 534 @code{fftw_plan}, executing it with @code{fftw_execute(plan)}, and
cannam@95 535 destroying it with @code{fftw_destroy_plan(plan)}. Furthermore, all
cannam@95 536 r2r transforms share the same planner interface:
cannam@95 537
cannam@95 538 @example
cannam@95 539 fftw_plan fftw_plan_r2r_1d(int n, double *in, double *out,
cannam@95 540 fftw_r2r_kind kind, unsigned flags);
cannam@95 541 fftw_plan fftw_plan_r2r_2d(int n0, int n1, double *in, double *out,
cannam@95 542 fftw_r2r_kind kind0, fftw_r2r_kind kind1,
cannam@95 543 unsigned flags);
cannam@95 544 fftw_plan fftw_plan_r2r_3d(int n0, int n1, int n2,
cannam@95 545 double *in, double *out,
cannam@95 546 fftw_r2r_kind kind0,
cannam@95 547 fftw_r2r_kind kind1,
cannam@95 548 fftw_r2r_kind kind2,
cannam@95 549 unsigned flags);
cannam@95 550 fftw_plan fftw_plan_r2r(int rank, const int *n, double *in, double *out,
cannam@95 551 const fftw_r2r_kind *kind, unsigned flags);
cannam@95 552 @end example
cannam@95 553 @findex fftw_plan_r2r_1d
cannam@95 554 @findex fftw_plan_r2r_2d
cannam@95 555 @findex fftw_plan_r2r_3d
cannam@95 556 @findex fftw_plan_r2r
cannam@95 557
cannam@95 558 Just as for the complex DFT, these plan 1d/2d/3d/multi-dimensional
cannam@95 559 transforms for contiguous arrays in row-major order, transforming (real)
cannam@95 560 input to output of the same size, where @code{n} specifies the
cannam@95 561 @emph{physical} dimensions of the arrays. All positive @code{n} are
cannam@95 562 supported (with the exception of @code{n=1} for the @code{FFTW_REDFT00}
cannam@95 563 kind, noted in the real-even subsection below); products of small
cannam@95 564 factors are most efficient (factorizing @code{n-1} and @code{n+1} for
cannam@95 565 @code{FFTW_REDFT00} and @code{FFTW_RODFT00} kinds, described below), but
cannam@95 566 an @Onlogn algorithm is used even for prime sizes.
cannam@95 567
cannam@95 568 Each dimension has a @dfn{kind} parameter, of type
cannam@95 569 @code{fftw_r2r_kind}, specifying the kind of r2r transform to be used
cannam@95 570 for that dimension.
cannam@95 571 @cindex kind (r2r)
cannam@95 572 @tindex fftw_r2r_kind
cannam@95 573 (In the case of @code{fftw_plan_r2r}, this is an array @code{kind[rank]}
cannam@95 574 where @code{kind[i]} is the transform kind for the dimension
cannam@95 575 @code{n[i]}.) The kind can be one of a set of predefined constants,
cannam@95 576 defined in the following subsections.
cannam@95 577
cannam@95 578 In other words, FFTW computes the separable product of the specified
cannam@95 579 r2r transforms over each dimension, which can be used e.g. for partial
cannam@95 580 differential equations with mixed boundary conditions. (For some r2r
cannam@95 581 kinds, notably the halfcomplex DFT and the DHT, such a separable
cannam@95 582 product is somewhat problematic in more than one dimension, however,
cannam@95 583 as is described below.)
cannam@95 584
cannam@95 585 In the current version of FFTW, all r2r transforms except for the
cannam@95 586 halfcomplex type are computed via pre- or post-processing of
cannam@95 587 halfcomplex transforms, and they are therefore not as fast as they
cannam@95 588 could be. Since most other general DCT/DST codes employ a similar
cannam@95 589 algorithm, however, FFTW's implementation should provide at least
cannam@95 590 competitive performance.
cannam@95 591
cannam@95 592 @c =========>
cannam@95 593 @node The Halfcomplex-format DFT, Real even/odd DFTs (cosine/sine transforms), More DFTs of Real Data, More DFTs of Real Data
cannam@95 594 @subsection The Halfcomplex-format DFT
cannam@95 595
cannam@95 596 An r2r kind of @code{FFTW_R2HC} (@dfn{r2hc}) corresponds to an r2c DFT
cannam@95 597 @ctindex FFTW_R2HC
cannam@95 598 @cindex r2c
cannam@95 599 @cindex r2hc
cannam@95 600 (@pxref{One-Dimensional DFTs of Real Data}) but with ``halfcomplex''
cannam@95 601 format output, and may sometimes be faster and/or more convenient than
cannam@95 602 the latter.
cannam@95 603 @cindex halfcomplex format
cannam@95 604 The inverse @dfn{hc2r} transform is of kind @code{FFTW_HC2R}.
cannam@95 605 @ctindex FFTW_HC2R
cannam@95 606 @cindex hc2r
cannam@95 607 This consists of the non-redundant half of the complex output for a 1d
cannam@95 608 real-input DFT of size @code{n}, stored as a sequence of @code{n} real
cannam@95 609 numbers (@code{double}) in the format:
cannam@95 610
cannam@95 611 @tex
cannam@95 612 $$
cannam@95 613 r_0, r_1, r_2, \ldots, r_{n/2}, i_{(n+1)/2-1}, \ldots, i_2, i_1
cannam@95 614 $$
cannam@95 615 @end tex
cannam@95 616 @ifinfo
cannam@95 617 r0, r1, r2, r(n/2), i((n+1)/2-1), ..., i2, i1
cannam@95 618 @end ifinfo
cannam@95 619 @html
cannam@95 620 <p align=center>
cannam@95 621 r<sub>0</sub>, r<sub>1</sub>, r<sub>2</sub>, ..., r<sub>n/2</sub>, i<sub>(n+1)/2-1</sub>, ..., i<sub>2</sub>, i<sub>1</sub>
cannam@95 622 </p>
cannam@95 623 @end html
cannam@95 624
cannam@95 625 Here,
cannam@95 626 @ifinfo
cannam@95 627 rk
cannam@95 628 @end ifinfo
cannam@95 629 @tex
cannam@95 630 $r_k$
cannam@95 631 @end tex
cannam@95 632 @html
cannam@95 633 r<sub>k</sub>
cannam@95 634 @end html
cannam@95 635 is the real part of the @math{k}th output, and
cannam@95 636 @ifinfo
cannam@95 637 ik
cannam@95 638 @end ifinfo
cannam@95 639 @tex
cannam@95 640 $i_k$
cannam@95 641 @end tex
cannam@95 642 @html
cannam@95 643 i<sub>k</sub>
cannam@95 644 @end html
cannam@95 645 is the imaginary part. (Division by 2 is rounded down.) For a
cannam@95 646 halfcomplex array @code{hc[n]}, the @math{k}th component thus has its
cannam@95 647 real part in @code{hc[k]} and its imaginary part in @code{hc[n-k]}, with
cannam@95 648 the exception of @code{k} @code{==} @code{0} or @code{n/2} (the latter
cannam@95 649 only if @code{n} is even)---in these two cases, the imaginary part is
cannam@95 650 zero due to symmetries of the real-input DFT, and is not stored.
cannam@95 651 Thus, the r2hc transform of @code{n} real values is a halfcomplex array of
cannam@95 652 length @code{n}, and vice versa for hc2r.
cannam@95 653 @cindex normalization
cannam@95 654
cannam@95 655
cannam@95 656 Aside from the differing format, the output of
cannam@95 657 @code{FFTW_R2HC}/@code{FFTW_HC2R} is otherwise exactly the same as for
cannam@95 658 the corresponding 1d r2c/c2r transform
cannam@95 659 (i.e. @code{FFTW_FORWARD}/@code{FFTW_BACKWARD} transforms, respectively).
cannam@95 660 Recall that these transforms are unnormalized, so r2hc followed by hc2r
cannam@95 661 will result in the original data multiplied by @code{n}. Furthermore,
cannam@95 662 like the c2r transform, an out-of-place hc2r transform will
cannam@95 663 @emph{destroy its input} array.
cannam@95 664
cannam@95 665 Although these halfcomplex transforms can be used with the
cannam@95 666 multi-dimensional r2r interface, the interpretation of such a separable
cannam@95 667 product of transforms along each dimension is problematic. For example,
cannam@95 668 consider a two-dimensional @code{n0} by @code{n1}, r2hc by r2hc
cannam@95 669 transform planned by @code{fftw_plan_r2r_2d(n0, n1, in, out, FFTW_R2HC,
cannam@95 670 FFTW_R2HC, FFTW_MEASURE)}. Conceptually, FFTW first transforms the rows
cannam@95 671 (of size @code{n1}) to produce halfcomplex rows, and then transforms the
cannam@95 672 columns (of size @code{n0}). Half of these column transforms, however,
cannam@95 673 are of imaginary parts, and should therefore be multiplied by @math{i}
cannam@95 674 and combined with the r2hc transforms of the real columns to produce the
cannam@95 675 2d DFT amplitudes; FFTW's r2r transform does @emph{not} perform this
cannam@95 676 combination for you. Thus, if a multi-dimensional real-input/output DFT
cannam@95 677 is required, we recommend using the ordinary r2c/c2r
cannam@95 678 interface (@pxref{Multi-Dimensional DFTs of Real Data}).
cannam@95 679
cannam@95 680 @c =========>
cannam@95 681 @node Real even/odd DFTs (cosine/sine transforms), The Discrete Hartley Transform, The Halfcomplex-format DFT, More DFTs of Real Data
cannam@95 682 @subsection Real even/odd DFTs (cosine/sine transforms)
cannam@95 683
cannam@95 684 The Fourier transform of a real-even function @math{f(-x) = f(x)} is
cannam@95 685 real-even, and @math{i} times the Fourier transform of a real-odd
cannam@95 686 function @math{f(-x) = -f(x)} is real-odd. Similar results hold for a
cannam@95 687 discrete Fourier transform, and thus for these symmetries the need for
cannam@95 688 complex inputs/outputs is entirely eliminated. Moreover, one gains a
cannam@95 689 factor of two in speed/space from the fact that the data are real, and
cannam@95 690 an additional factor of two from the even/odd symmetry: only the
cannam@95 691 non-redundant (first) half of the array need be stored. The result is
cannam@95 692 the real-even DFT (@dfn{REDFT}) and the real-odd DFT (@dfn{RODFT}), also
cannam@95 693 known as the discrete cosine and sine transforms (@dfn{DCT} and
cannam@95 694 @dfn{DST}), respectively.
cannam@95 695 @cindex real-even DFT
cannam@95 696 @cindex REDFT
cannam@95 697 @cindex real-odd DFT
cannam@95 698 @cindex RODFT
cannam@95 699 @cindex discrete cosine transform
cannam@95 700 @cindex DCT
cannam@95 701 @cindex discrete sine transform
cannam@95 702 @cindex DST
cannam@95 703
cannam@95 704
cannam@95 705 (In this section, we describe the 1d transforms; multi-dimensional
cannam@95 706 transforms are just a separable product of these transforms operating
cannam@95 707 along each dimension.)
cannam@95 708
cannam@95 709 Because of the discrete sampling, one has an additional choice: is the
cannam@95 710 data even/odd around a sampling point, or around the point halfway
cannam@95 711 between two samples? The latter corresponds to @emph{shifting} the
cannam@95 712 samples by @emph{half} an interval, and gives rise to several transform
cannam@95 713 variants denoted by REDFT@math{ab} and RODFT@math{ab}: @math{a} and
cannam@95 714 @math{b} are @math{0} or @math{1}, and indicate whether the input
cannam@95 715 (@math{a}) and/or output (@math{b}) are shifted by half a sample
cannam@95 716 (@math{1} means it is shifted). These are also known as types I-IV of
cannam@95 717 the DCT and DST, and all four types are supported by FFTW's r2r
cannam@95 718 interface.@footnote{There are also type V-VIII transforms, which
cannam@95 719 correspond to a logical DFT of @emph{odd} size @math{N}, independent of
cannam@95 720 whether the physical size @code{n} is odd, but we do not support these
cannam@95 721 variants.}
cannam@95 722
cannam@95 723 The r2r kinds for the various REDFT and RODFT types supported by FFTW,
cannam@95 724 along with the boundary conditions at both ends of the @emph{input}
cannam@95 725 array (@code{n} real numbers @code{in[j=0..n-1]}), are:
cannam@95 726
cannam@95 727 @itemize @bullet
cannam@95 728
cannam@95 729 @item
cannam@95 730 @code{FFTW_REDFT00} (DCT-I): even around @math{j=0} and even around @math{j=n-1}.
cannam@95 731 @ctindex FFTW_REDFT00
cannam@95 732
cannam@95 733 @item
cannam@95 734 @code{FFTW_REDFT10} (DCT-II, ``the'' DCT): even around @math{j=-0.5} and even around @math{j=n-0.5}.
cannam@95 735 @ctindex FFTW_REDFT10
cannam@95 736
cannam@95 737 @item
cannam@95 738 @code{FFTW_REDFT01} (DCT-III, ``the'' IDCT): even around @math{j=0} and odd around @math{j=n}.
cannam@95 739 @ctindex FFTW_REDFT01
cannam@95 740 @cindex IDCT
cannam@95 741
cannam@95 742 @item
cannam@95 743 @code{FFTW_REDFT11} (DCT-IV): even around @math{j=-0.5} and odd around @math{j=n-0.5}.
cannam@95 744 @ctindex FFTW_REDFT11
cannam@95 745
cannam@95 746 @item
cannam@95 747 @code{FFTW_RODFT00} (DST-I): odd around @math{j=-1} and odd around @math{j=n}.
cannam@95 748 @ctindex FFTW_RODFT00
cannam@95 749
cannam@95 750 @item
cannam@95 751 @code{FFTW_RODFT10} (DST-II): odd around @math{j=-0.5} and odd around @math{j=n-0.5}.
cannam@95 752 @ctindex FFTW_RODFT10
cannam@95 753
cannam@95 754 @item
cannam@95 755 @code{FFTW_RODFT01} (DST-III): odd around @math{j=-1} and even around @math{j=n-1}.
cannam@95 756 @ctindex FFTW_RODFT01
cannam@95 757
cannam@95 758 @item
cannam@95 759 @code{FFTW_RODFT11} (DST-IV): odd around @math{j=-0.5} and even around @math{j=n-0.5}.
cannam@95 760 @ctindex FFTW_RODFT11
cannam@95 761
cannam@95 762 @end itemize
cannam@95 763
cannam@95 764 Note that these symmetries apply to the ``logical'' array being
cannam@95 765 transformed; @strong{there are no constraints on your physical input
cannam@95 766 data}. So, for example, if you specify a size-5 REDFT00 (DCT-I) of the
cannam@95 767 data @math{abcde}, it corresponds to the DFT of the logical even array
cannam@95 768 @math{abcdedcb} of size 8. A size-4 REDFT10 (DCT-II) of the data
cannam@95 769 @math{abcd} corresponds to the size-8 logical DFT of the even array
cannam@95 770 @math{abcddcba}, shifted by half a sample.
cannam@95 771
cannam@95 772 All of these transforms are invertible. The inverse of R*DFT00 is
cannam@95 773 R*DFT00; of R*DFT10 is R*DFT01 and vice versa (these are often called
cannam@95 774 simply ``the'' DCT and IDCT, respectively); and of R*DFT11 is R*DFT11.
cannam@95 775 However, the transforms computed by FFTW are unnormalized, exactly
cannam@95 776 like the corresponding real and complex DFTs, so computing a transform
cannam@95 777 followed by its inverse yields the original array scaled by @math{N},
cannam@95 778 where @math{N} is the @emph{logical} DFT size. For REDFT00,
cannam@95 779 @math{N=2(n-1)}; for RODFT00, @math{N=2(n+1)}; otherwise, @math{N=2n}.
cannam@95 780 @cindex normalization
cannam@95 781 @cindex IDCT
cannam@95 782
cannam@95 783
cannam@95 784 Note that the boundary conditions of the transform output array are
cannam@95 785 given by the input boundary conditions of the inverse transform.
cannam@95 786 Thus, the above transforms are all inequivalent in terms of
cannam@95 787 input/output boundary conditions, even neglecting the 0.5 shift
cannam@95 788 difference.
cannam@95 789
cannam@95 790 FFTW is most efficient when @math{N} is a product of small factors; note
cannam@95 791 that this @emph{differs} from the factorization of the physical size
cannam@95 792 @code{n} for REDFT00 and RODFT00! There is another oddity: @code{n=1}
cannam@95 793 REDFT00 transforms correspond to @math{N=0}, and so are @emph{not
cannam@95 794 defined} (the planner will return @code{NULL}). Otherwise, any positive
cannam@95 795 @code{n} is supported.
cannam@95 796
cannam@95 797 For the precise mathematical definitions of these transforms as used by
cannam@95 798 FFTW, see @ref{What FFTW Really Computes}. (For people accustomed to
cannam@95 799 the DCT/DST, FFTW's definitions have a coefficient of @math{2} in front
cannam@95 800 of the cos/sin functions so that they correspond precisely to an
cannam@95 801 even/odd DFT of size @math{N}. Some authors also include additional
cannam@95 802 multiplicative factors of
cannam@95 803 @ifinfo
cannam@95 804 sqrt(2)
cannam@95 805 @end ifinfo
cannam@95 806 @html
cannam@95 807 &radic;2
cannam@95 808 @end html
cannam@95 809 @tex
cannam@95 810 $\sqrt{2}$
cannam@95 811 @end tex
cannam@95 812 for selected inputs and outputs; this makes
cannam@95 813 the transform orthogonal, but sacrifices the direct equivalence to a
cannam@95 814 symmetric DFT.)
cannam@95 815
cannam@95 816 @subsubheading Which type do you need?
cannam@95 817
cannam@95 818 Since the required flavor of even/odd DFT depends upon your problem,
cannam@95 819 you are the best judge of this choice, but we can make a few comments
cannam@95 820 on relative efficiency to help you in your selection. In particular,
cannam@95 821 R*DFT01 and R*DFT10 tend to be slightly faster than R*DFT11
cannam@95 822 (especially for odd sizes), while the R*DFT00 transforms are sometimes
cannam@95 823 significantly slower (especially for even sizes).@footnote{R*DFT00 is
cannam@95 824 sometimes slower in FFTW because we discovered that the standard
cannam@95 825 algorithm for computing this by a pre/post-processed real DFT---the
cannam@95 826 algorithm used in FFTPACK, Numerical Recipes, and other sources for
cannam@95 827 decades now---has serious numerical problems: it already loses several
cannam@95 828 decimal places of accuracy for 16k sizes. There seem to be only two
cannam@95 829 alternatives in the literature that do not suffer similarly: a
cannam@95 830 recursive decomposition into smaller DCTs, which would require a large
cannam@95 831 set of codelets for efficiency and generality, or sacrificing a factor of
cannam@95 832 @tex
cannam@95 833 $\sim 2$
cannam@95 834 @end tex
cannam@95 835 @ifnottex
cannam@95 836 2
cannam@95 837 @end ifnottex
cannam@95 838 in speed to use a real DFT of twice the size. We currently
cannam@95 839 employ the latter technique for general @math{n}, as well as a limited
cannam@95 840 form of the former method: a split-radix decomposition when @math{n}
cannam@95 841 is odd (@math{N} a multiple of 4). For @math{N} containing many
cannam@95 842 factors of 2, the split-radix method seems to recover most of the
cannam@95 843 speed of the standard algorithm without the accuracy tradeoff.}
cannam@95 844
cannam@95 845 Thus, if only the boundary conditions on the transform inputs are
cannam@95 846 specified, we generally recommend R*DFT10 over R*DFT00 and R*DFT01 over
cannam@95 847 R*DFT11 (unless the half-sample shift or the self-inverse property is
cannam@95 848 significant for your problem).
cannam@95 849
cannam@95 850 If performance is important to you and you are using only small sizes
cannam@95 851 (say @math{n<200}), e.g. for multi-dimensional transforms, then you
cannam@95 852 might consider generating hard-coded transforms of those sizes and types
cannam@95 853 that you are interested in (@pxref{Generating your own code}).
cannam@95 854
cannam@95 855 We are interested in hearing what types of symmetric transforms you find
cannam@95 856 most useful.
cannam@95 857
cannam@95 858 @c =========>
cannam@95 859 @node The Discrete Hartley Transform, , Real even/odd DFTs (cosine/sine transforms), More DFTs of Real Data
cannam@95 860 @subsection The Discrete Hartley Transform
cannam@95 861
cannam@95 862 If you are planning to use the DHT because you've heard that it is
cannam@95 863 ``faster'' than the DFT (FFT), @strong{stop here}. The DHT is not
cannam@95 864 faster than the DFT. That story is an old but enduring misconception
cannam@95 865 that was debunked in 1987.
cannam@95 866
cannam@95 867 The discrete Hartley transform (DHT) is an invertible linear transform
cannam@95 868 closely related to the DFT. In the DFT, one multiplies each input by
cannam@95 869 @math{cos - i * sin} (a complex exponential), whereas in the DHT each
cannam@95 870 input is multiplied by simply @math{cos + sin}. Thus, the DHT
cannam@95 871 transforms @code{n} real numbers to @code{n} real numbers, and has the
cannam@95 872 convenient property of being its own inverse. In FFTW, a DHT (of any
cannam@95 873 positive @code{n}) can be specified by an r2r kind of @code{FFTW_DHT}.
cannam@95 874 @ctindex FFTW_DHT
cannam@95 875 @cindex discrete Hartley transform
cannam@95 876 @cindex DHT
cannam@95 877
cannam@95 878 Like the DFT, in FFTW the DHT is unnormalized, so computing a DHT of
cannam@95 879 size @code{n} followed by another DHT of the same size will result in
cannam@95 880 the original array multiplied by @code{n}.
cannam@95 881 @cindex normalization
cannam@95 882
cannam@95 883 The DHT was originally proposed as a more efficient alternative to the
cannam@95 884 DFT for real data, but it was subsequently shown that a specialized DFT
cannam@95 885 (such as FFTW's r2hc or r2c transforms) could be just as fast. In FFTW,
cannam@95 886 the DHT is actually computed by post-processing an r2hc transform, so
cannam@95 887 there is ordinarily no reason to prefer it from a performance
cannam@95 888 perspective.@footnote{We provide the DHT mainly as a byproduct of some
cannam@95 889 internal algorithms. FFTW computes a real input/output DFT of
cannam@95 890 @emph{prime} size by re-expressing it as a DHT plus post/pre-processing
cannam@95 891 and then using Rader's prime-DFT algorithm adapted to the DHT.}
cannam@95 892 However, we have heard rumors that the DHT might be the most appropriate
cannam@95 893 transform in its own right for certain applications, and we would be
cannam@95 894 very interested to hear from anyone who finds it useful.
cannam@95 895
cannam@95 896 If @code{FFTW_DHT} is specified for multiple dimensions of a
cannam@95 897 multi-dimensional transform, FFTW computes the separable product of 1d
cannam@95 898 DHTs along each dimension. Unfortunately, this is not quite the same
cannam@95 899 thing as a true multi-dimensional DHT; you can compute the latter, if
cannam@95 900 necessary, with at most @code{rank-1} post-processing passes
cannam@95 901 [see e.g. H. Hao and R. N. Bracewell, @i{Proc. IEEE} @b{75}, 264--266 (1987)].
cannam@95 902
cannam@95 903 For the precise mathematical definition of the DHT as used by FFTW, see
cannam@95 904 @ref{What FFTW Really Computes}.
cannam@95 905