annotate src/fftw-3.3.5/doc/reference.texi @ 127:7867fa7e1b6b

Current fftw source
author Chris Cannam <cannam@all-day-breakfast.com>
date Tue, 18 Oct 2016 13:40:26 +0100
parents
children
rev   line source
cannam@127 1 @node FFTW Reference, Multi-threaded FFTW, Other Important Topics, Top
cannam@127 2 @chapter FFTW Reference
cannam@127 3
cannam@127 4 This chapter provides a complete reference for all sequential (i.e.,
cannam@127 5 one-processor) FFTW functions. Parallel transforms are described in
cannam@127 6 later chapters.
cannam@127 7
cannam@127 8 @menu
cannam@127 9 * Data Types and Files::
cannam@127 10 * Using Plans::
cannam@127 11 * Basic Interface::
cannam@127 12 * Advanced Interface::
cannam@127 13 * Guru Interface::
cannam@127 14 * New-array Execute Functions::
cannam@127 15 * Wisdom::
cannam@127 16 * What FFTW Really Computes::
cannam@127 17 @end menu
cannam@127 18
cannam@127 19 @c ------------------------------------------------------------
cannam@127 20 @node Data Types and Files, Using Plans, FFTW Reference, FFTW Reference
cannam@127 21 @section Data Types and Files
cannam@127 22
cannam@127 23 All programs using FFTW should include its header file:
cannam@127 24
cannam@127 25 @example
cannam@127 26 #include <fftw3.h>
cannam@127 27 @end example
cannam@127 28
cannam@127 29 You must also link to the FFTW library. On Unix, this
cannam@127 30 means adding @code{-lfftw3 -lm} at the @emph{end} of the link command.
cannam@127 31
cannam@127 32 @menu
cannam@127 33 * Complex numbers::
cannam@127 34 * Precision::
cannam@127 35 * Memory Allocation::
cannam@127 36 @end menu
cannam@127 37
cannam@127 38 @c =========>
cannam@127 39 @node Complex numbers, Precision, Data Types and Files, Data Types and Files
cannam@127 40 @subsection Complex numbers
cannam@127 41
cannam@127 42 The default FFTW interface uses @code{double} precision for all
cannam@127 43 floating-point numbers, and defines a @code{fftw_complex} type to hold
cannam@127 44 complex numbers as:
cannam@127 45
cannam@127 46 @example
cannam@127 47 typedef double fftw_complex[2];
cannam@127 48 @end example
cannam@127 49 @tindex fftw_complex
cannam@127 50
cannam@127 51 Here, the @code{[0]} element holds the real part and the @code{[1]}
cannam@127 52 element holds the imaginary part.
cannam@127 53
cannam@127 54 Alternatively, if you have a C compiler (such as @code{gcc}) that
cannam@127 55 supports the C99 revision of the ANSI C standard, you can use C's new
cannam@127 56 native complex type (which is binary-compatible with the typedef above).
cannam@127 57 In particular, if you @code{#include <complex.h>} @emph{before}
cannam@127 58 @code{<fftw3.h>}, then @code{fftw_complex} is defined to be the native
cannam@127 59 complex type and you can manipulate it with ordinary arithmetic
cannam@127 60 (e.g. @code{x = y * (3+4*I)}, where @code{x} and @code{y} are
cannam@127 61 @code{fftw_complex} and @code{I} is the standard symbol for the
cannam@127 62 imaginary unit);
cannam@127 63 @cindex C99
cannam@127 64
cannam@127 65
cannam@127 66 C++ has its own @code{complex<T>} template class, defined in the
cannam@127 67 standard @code{<complex>} header file. Reportedly, the C++ standards
cannam@127 68 committee has recently agreed to mandate that the storage format used
cannam@127 69 for this type be binary-compatible with the C99 type, i.e. an array
cannam@127 70 @code{T[2]} with consecutive real @code{[0]} and imaginary @code{[1]}
cannam@127 71 parts. (See report
cannam@127 72 @uref{http://www.open-std.org/jtc1/sc22/WG21/docs/papers/2002/n1388.pdf
cannam@127 73 WG21/N1388}.) Although not part of the official standard as of this
cannam@127 74 writing, the proposal stated that: ``This solution has been tested with
cannam@127 75 all current major implementations of the standard library and shown to
cannam@127 76 be working.'' To the extent that this is true, if you have a variable
cannam@127 77 @code{complex<double> *x}, you can pass it directly to FFTW via
cannam@127 78 @code{reinterpret_cast<fftw_complex*>(x)}.
cannam@127 79 @cindex C++
cannam@127 80 @cindex portability
cannam@127 81
cannam@127 82 @c =========>
cannam@127 83 @node Precision, Memory Allocation, Complex numbers, Data Types and Files
cannam@127 84 @subsection Precision
cannam@127 85 @cindex precision
cannam@127 86
cannam@127 87 You can install single and long-double precision versions of FFTW,
cannam@127 88 which replace @code{double} with @code{float} and @code{long double},
cannam@127 89 respectively (@pxref{Installation and Customization}). To use these
cannam@127 90 interfaces, you:
cannam@127 91
cannam@127 92 @itemize @bullet
cannam@127 93
cannam@127 94 @item
cannam@127 95 Link to the single/long-double libraries; on Unix, @code{-lfftw3f} or
cannam@127 96 @code{-lfftw3l} instead of (or in addition to) @code{-lfftw3}. (You
cannam@127 97 can link to the different-precision libraries simultaneously.)
cannam@127 98
cannam@127 99 @item
cannam@127 100 Include the @emph{same} @code{<fftw3.h>} header file.
cannam@127 101
cannam@127 102 @item
cannam@127 103 Replace all lowercase instances of @samp{fftw_} with @samp{fftwf_} or
cannam@127 104 @samp{fftwl_} for single or long-double precision, respectively.
cannam@127 105 (@code{fftw_complex} becomes @code{fftwf_complex}, @code{fftw_execute}
cannam@127 106 becomes @code{fftwf_execute}, etcetera.)
cannam@127 107
cannam@127 108 @item
cannam@127 109 Uppercase names, i.e. names beginning with @samp{FFTW_}, remain the
cannam@127 110 same.
cannam@127 111
cannam@127 112 @item
cannam@127 113 Replace @code{double} with @code{float} or @code{long double} for
cannam@127 114 subroutine parameters.
cannam@127 115
cannam@127 116 @end itemize
cannam@127 117
cannam@127 118 Depending upon your compiler and/or hardware, @code{long double} may not
cannam@127 119 be any more precise than @code{double} (or may not be supported at all,
cannam@127 120 although it is standard in C99).
cannam@127 121 @cindex C99
cannam@127 122
cannam@127 123
cannam@127 124 We also support using the nonstandard @code{__float128}
cannam@127 125 quadruple-precision type provided by recent versions of @code{gcc} on
cannam@127 126 32- and 64-bit x86 hardware (@pxref{Installation and Customization}).
cannam@127 127 To use this type, link with @code{-lfftw3q -lquadmath -lm} (the
cannam@127 128 @code{libquadmath} library provided by @code{gcc} is needed for
cannam@127 129 quadruple-precision trigonometric functions) and use @samp{fftwq_}
cannam@127 130 identifiers.
cannam@127 131
cannam@127 132 @c =========>
cannam@127 133 @node Memory Allocation, , Precision, Data Types and Files
cannam@127 134 @subsection Memory Allocation
cannam@127 135
cannam@127 136 @example
cannam@127 137 void *fftw_malloc(size_t n);
cannam@127 138 void fftw_free(void *p);
cannam@127 139 @end example
cannam@127 140 @findex fftw_malloc
cannam@127 141 @findex fftw_free
cannam@127 142
cannam@127 143 These are functions that behave identically to @code{malloc} and
cannam@127 144 @code{free}, except that they guarantee that the returned pointer obeys
cannam@127 145 any special alignment restrictions imposed by any algorithm in FFTW
cannam@127 146 (e.g. for SIMD acceleration). @xref{SIMD alignment and fftw_malloc}.
cannam@127 147 @cindex alignment
cannam@127 148
cannam@127 149
cannam@127 150 Data allocated by @code{fftw_malloc} @emph{must} be deallocated by
cannam@127 151 @code{fftw_free} and not by the ordinary @code{free}.
cannam@127 152
cannam@127 153 These routines simply call through to your operating system's
cannam@127 154 @code{malloc} or, if necessary, its aligned equivalent
cannam@127 155 (e.g. @code{memalign}), so you normally need not worry about any
cannam@127 156 significant time or space overhead. You are @emph{not required} to use
cannam@127 157 them to allocate your data, but we strongly recommend it.
cannam@127 158
cannam@127 159 Note: in C++, just as with ordinary @code{malloc}, you must typecast
cannam@127 160 the output of @code{fftw_malloc} to whatever pointer type you are
cannam@127 161 allocating.
cannam@127 162 @cindex C++
cannam@127 163
cannam@127 164
cannam@127 165 We also provide the following two convenience functions to allocate
cannam@127 166 real and complex arrays with @code{n} elements, which are equivalent
cannam@127 167 to @code{(double *) fftw_malloc(sizeof(double) * n)} and
cannam@127 168 @code{(fftw_complex *) fftw_malloc(sizeof(fftw_complex) * n)},
cannam@127 169 respectively:
cannam@127 170
cannam@127 171 @example
cannam@127 172 double *fftw_alloc_real(size_t n);
cannam@127 173 fftw_complex *fftw_alloc_complex(size_t n);
cannam@127 174 @end example
cannam@127 175 @findex fftw_alloc_real
cannam@127 176 @findex fftw_alloc_complex
cannam@127 177
cannam@127 178 The equivalent functions in other precisions allocate arrays of @code{n}
cannam@127 179 elements in that precision. e.g. @code{fftwf_alloc_real(n)} is
cannam@127 180 equivalent to @code{(float *) fftwf_malloc(sizeof(float) * n)}.
cannam@127 181 @cindex precision
cannam@127 182
cannam@127 183 @c ------------------------------------------------------------
cannam@127 184 @node Using Plans, Basic Interface, Data Types and Files, FFTW Reference
cannam@127 185 @section Using Plans
cannam@127 186
cannam@127 187 Plans for all transform types in FFTW are stored as type
cannam@127 188 @code{fftw_plan} (an opaque pointer type), and are created by one of the
cannam@127 189 various planning routines described in the following sections.
cannam@127 190 @tindex fftw_plan
cannam@127 191 An @code{fftw_plan} contains all information necessary to compute the
cannam@127 192 transform, including the pointers to the input and output arrays.
cannam@127 193
cannam@127 194 @example
cannam@127 195 void fftw_execute(const fftw_plan plan);
cannam@127 196 @end example
cannam@127 197 @findex fftw_execute
cannam@127 198
cannam@127 199 This executes the @code{plan}, to compute the corresponding transform on
cannam@127 200 the arrays for which it was planned (which must still exist). The plan
cannam@127 201 is not modified, and @code{fftw_execute} can be called as many times as
cannam@127 202 desired.
cannam@127 203
cannam@127 204 To apply a given plan to a different array, you can use the new-array execute
cannam@127 205 interface. @xref{New-array Execute Functions}.
cannam@127 206
cannam@127 207 @code{fftw_execute} (and equivalents) is the only function in FFTW
cannam@127 208 guaranteed to be thread-safe; see @ref{Thread safety}.
cannam@127 209
cannam@127 210 This function:
cannam@127 211 @example
cannam@127 212 void fftw_destroy_plan(fftw_plan plan);
cannam@127 213 @end example
cannam@127 214 @findex fftw_destroy_plan
cannam@127 215 deallocates the @code{plan} and all its associated data.
cannam@127 216
cannam@127 217 FFTW's planner saves some other persistent data, such as the
cannam@127 218 accumulated wisdom and a list of algorithms available in the current
cannam@127 219 configuration. If you want to deallocate all of that and reset FFTW
cannam@127 220 to the pristine state it was in when you started your program, you can
cannam@127 221 call:
cannam@127 222
cannam@127 223 @example
cannam@127 224 void fftw_cleanup(void);
cannam@127 225 @end example
cannam@127 226 @findex fftw_cleanup
cannam@127 227
cannam@127 228 After calling @code{fftw_cleanup}, all existing plans become undefined,
cannam@127 229 and you should not attempt to execute them nor to destroy them. You can
cannam@127 230 however create and execute/destroy new plans, in which case FFTW starts
cannam@127 231 accumulating wisdom information again.
cannam@127 232
cannam@127 233 @code{fftw_cleanup} does not deallocate your plans, however. To prevent
cannam@127 234 memory leaks, you must still call @code{fftw_destroy_plan} before
cannam@127 235 executing @code{fftw_cleanup}.
cannam@127 236
cannam@127 237 Occasionally, it may useful to know FFTW's internal ``cost'' metric
cannam@127 238 that it uses to compare plans to one another; this cost is
cannam@127 239 proportional to an execution time of the plan, in undocumented units,
cannam@127 240 if the plan was created with the @code{FFTW_MEASURE} or other
cannam@127 241 timing-based options, or alternatively is a heuristic cost function
cannam@127 242 for @code{FFTW_ESTIMATE} plans. (The cost values of measured and
cannam@127 243 estimated plans are not comparable, being in different units. Also,
cannam@127 244 costs from different FFTW versions or the same version compiled
cannam@127 245 differently may not be in the same units. Plans created from wisdom
cannam@127 246 have a cost of 0 since no timing measurement is performed for them.
cannam@127 247 Finally, certain problems for which only one top-level algorithm was
cannam@127 248 possible may have required no measurements of the cost of the whole
cannam@127 249 plan, in which case @code{fftw_cost} will also return 0.) The cost
cannam@127 250 metric for a given plan is returned by:
cannam@127 251
cannam@127 252 @example
cannam@127 253 double fftw_cost(const fftw_plan plan);
cannam@127 254 @end example
cannam@127 255 @findex fftw_cost
cannam@127 256
cannam@127 257 The following two routines are provided purely for academic purposes
cannam@127 258 (that is, for entertainment).
cannam@127 259
cannam@127 260 @example
cannam@127 261 void fftw_flops(const fftw_plan plan,
cannam@127 262 double *add, double *mul, double *fma);
cannam@127 263 @end example
cannam@127 264 @findex fftw_flops
cannam@127 265
cannam@127 266 Given a @code{plan}, set @code{add}, @code{mul}, and @code{fma} to an
cannam@127 267 exact count of the number of floating-point additions, multiplications,
cannam@127 268 and fused multiply-add operations involved in the plan's execution. The
cannam@127 269 total number of floating-point operations (flops) is @code{add + mul +
cannam@127 270 2*fma}, or @code{add + mul + fma} if the hardware supports fused
cannam@127 271 multiply-add instructions (although the number of FMA operations is only
cannam@127 272 approximate because of compiler voodoo). (The number of operations
cannam@127 273 should be an integer, but we use @code{double} to avoid overflowing
cannam@127 274 @code{int} for large transforms; the arguments are of type @code{double}
cannam@127 275 even for single and long-double precision versions of FFTW.)
cannam@127 276
cannam@127 277 @example
cannam@127 278 void fftw_fprint_plan(const fftw_plan plan, FILE *output_file);
cannam@127 279 void fftw_print_plan(const fftw_plan plan);
cannam@127 280 char *fftw_sprint_plan(const fftw_plan plan);
cannam@127 281 @end example
cannam@127 282 @findex fftw_fprint_plan
cannam@127 283 @findex fftw_print_plan
cannam@127 284
cannam@127 285 This outputs a ``nerd-readable'' representation of the @code{plan} to
cannam@127 286 the given file, to @code{stdout}, or two a newly allocated
cannam@127 287 NUL-terminated string (which the caller is responsible for deallocating
cannam@127 288 with @code{free}), respectively.
cannam@127 289
cannam@127 290 @c ------------------------------------------------------------
cannam@127 291 @node Basic Interface, Advanced Interface, Using Plans, FFTW Reference
cannam@127 292 @section Basic Interface
cannam@127 293 @cindex basic interface
cannam@127 294
cannam@127 295 Recall that the FFTW API is divided into three parts@footnote{@i{Gallia est
cannam@127 296 omnis divisa in partes tres} (Julius Caesar).}: the @dfn{basic interface}
cannam@127 297 computes a single transform of contiguous data, the @dfn{advanced
cannam@127 298 interface} computes transforms of multiple or strided arrays, and the
cannam@127 299 @dfn{guru interface} supports the most general data layouts,
cannam@127 300 multiplicities, and strides. This section describes the the basic
cannam@127 301 interface, which we expect to satisfy the needs of most users.
cannam@127 302
cannam@127 303 @menu
cannam@127 304 * Complex DFTs::
cannam@127 305 * Planner Flags::
cannam@127 306 * Real-data DFTs::
cannam@127 307 * Real-data DFT Array Format::
cannam@127 308 * Real-to-Real Transforms::
cannam@127 309 * Real-to-Real Transform Kinds::
cannam@127 310 @end menu
cannam@127 311
cannam@127 312 @c =========>
cannam@127 313 @node Complex DFTs, Planner Flags, Basic Interface, Basic Interface
cannam@127 314 @subsection Complex DFTs
cannam@127 315
cannam@127 316 @example
cannam@127 317 fftw_plan fftw_plan_dft_1d(int n0,
cannam@127 318 fftw_complex *in, fftw_complex *out,
cannam@127 319 int sign, unsigned flags);
cannam@127 320 fftw_plan fftw_plan_dft_2d(int n0, int n1,
cannam@127 321 fftw_complex *in, fftw_complex *out,
cannam@127 322 int sign, unsigned flags);
cannam@127 323 fftw_plan fftw_plan_dft_3d(int n0, int n1, int n2,
cannam@127 324 fftw_complex *in, fftw_complex *out,
cannam@127 325 int sign, unsigned flags);
cannam@127 326 fftw_plan fftw_plan_dft(int rank, const int *n,
cannam@127 327 fftw_complex *in, fftw_complex *out,
cannam@127 328 int sign, unsigned flags);
cannam@127 329 @end example
cannam@127 330 @findex fftw_plan_dft_1d
cannam@127 331 @findex fftw_plan_dft_2d
cannam@127 332 @findex fftw_plan_dft_3d
cannam@127 333 @findex fftw_plan_dft
cannam@127 334
cannam@127 335 Plan a complex input/output discrete Fourier transform (DFT) in zero or
cannam@127 336 more dimensions, returning an @code{fftw_plan} (@pxref{Using Plans}).
cannam@127 337
cannam@127 338 Once you have created a plan for a certain transform type and
cannam@127 339 parameters, then creating another plan of the same type and parameters,
cannam@127 340 but for different arrays, is fast and shares constant data with the
cannam@127 341 first plan (if it still exists).
cannam@127 342
cannam@127 343 The planner returns @code{NULL} if the plan cannot be created. In the
cannam@127 344 standard FFTW distribution, the basic interface is guaranteed to return
cannam@127 345 a non-@code{NULL} plan. A plan may be @code{NULL}, however, if you are
cannam@127 346 using a customized FFTW configuration supporting a restricted set of
cannam@127 347 transforms.
cannam@127 348
cannam@127 349 @subsubheading Arguments
cannam@127 350 @itemize @bullet
cannam@127 351
cannam@127 352 @item
cannam@127 353 @code{rank} is the rank of the transform (it should be the size of the
cannam@127 354 array @code{*n}), and can be any non-negative integer. (@xref{Complex
cannam@127 355 Multi-Dimensional DFTs}, for the definition of ``rank''.) The
cannam@127 356 @samp{_1d}, @samp{_2d}, and @samp{_3d} planners correspond to a
cannam@127 357 @code{rank} of @code{1}, @code{2}, and @code{3}, respectively. The rank
cannam@127 358 may be zero, which is equivalent to a rank-1 transform of size 1, i.e. a
cannam@127 359 copy of one number from input to output.
cannam@127 360
cannam@127 361 @item
cannam@127 362 @code{n0}, @code{n1}, @code{n2}, or @code{n[0..rank-1]} (as appropriate
cannam@127 363 for each routine) specify the size of the transform dimensions. They
cannam@127 364 can be any positive integer.
cannam@127 365
cannam@127 366 @itemize @minus
cannam@127 367 @item
cannam@127 368 @cindex row-major
cannam@127 369 Multi-dimensional arrays are stored in row-major order with dimensions:
cannam@127 370 @code{n0} x @code{n1}; or @code{n0} x @code{n1} x @code{n2}; or
cannam@127 371 @code{n[0]} x @code{n[1]} x ... x @code{n[rank-1]}.
cannam@127 372 @xref{Multi-dimensional Array Format}.
cannam@127 373 @item
cannam@127 374 FFTW is best at handling sizes of the form
cannam@127 375 @ifinfo
cannam@127 376 @math{2^a 3^b 5^c 7^d 11^e 13^f},
cannam@127 377 @end ifinfo
cannam@127 378 @tex
cannam@127 379 $2^a 3^b 5^c 7^d 11^e 13^f$,
cannam@127 380 @end tex
cannam@127 381 @html
cannam@127 382 2<sup>a</sup> 3<sup>b</sup> 5<sup>c</sup> 7<sup>d</sup>
cannam@127 383 11<sup>e</sup> 13<sup>f</sup>,
cannam@127 384 @end html
cannam@127 385 where @math{e+f} is either @math{0} or @math{1}, and the other exponents
cannam@127 386 are arbitrary. Other sizes are computed by means of a slow,
cannam@127 387 general-purpose algorithm (which nevertheless retains @Onlogn{} performance even for prime sizes). It is possible to customize FFTW
cannam@127 388 for different array sizes; see @ref{Installation and Customization}.
cannam@127 389 Transforms whose sizes are powers of @math{2} are especially fast.
cannam@127 390 @end itemize
cannam@127 391
cannam@127 392 @item
cannam@127 393 @code{in} and @code{out} point to the input and output arrays of the
cannam@127 394 transform, which may be the same (yielding an in-place transform).
cannam@127 395 @cindex in-place
cannam@127 396 These arrays are overwritten during planning, unless
cannam@127 397 @code{FFTW_ESTIMATE} is used in the flags. (The arrays need not be
cannam@127 398 initialized, but they must be allocated.)
cannam@127 399
cannam@127 400 If @code{in == out}, the transform is @dfn{in-place} and the input
cannam@127 401 array is overwritten. If @code{in != out}, the two arrays must
cannam@127 402 not overlap (but FFTW does not check for this condition).
cannam@127 403
cannam@127 404 @item
cannam@127 405 @ctindex FFTW_FORWARD
cannam@127 406 @ctindex FFTW_BACKWARD
cannam@127 407 @code{sign} is the sign of the exponent in the formula that defines the
cannam@127 408 Fourier transform. It can be @math{-1} (= @code{FFTW_FORWARD}) or
cannam@127 409 @math{+1} (= @code{FFTW_BACKWARD}).
cannam@127 410
cannam@127 411 @item
cannam@127 412 @cindex flags
cannam@127 413 @code{flags} is a bitwise OR (@samp{|}) of zero or more planner flags,
cannam@127 414 as defined in @ref{Planner Flags}.
cannam@127 415
cannam@127 416 @end itemize
cannam@127 417
cannam@127 418 FFTW computes an unnormalized transform: computing a forward followed by
cannam@127 419 a backward transform (or vice versa) will result in the original data
cannam@127 420 multiplied by the size of the transform (the product of the dimensions).
cannam@127 421 @cindex normalization
cannam@127 422 For more information, see @ref{What FFTW Really Computes}.
cannam@127 423
cannam@127 424 @c =========>
cannam@127 425 @node Planner Flags, Real-data DFTs, Complex DFTs, Basic Interface
cannam@127 426 @subsection Planner Flags
cannam@127 427
cannam@127 428 All of the planner routines in FFTW accept an integer @code{flags}
cannam@127 429 argument, which is a bitwise OR (@samp{|}) of zero or more of the flag
cannam@127 430 constants defined below. These flags control the rigor (and time) of
cannam@127 431 the planning process, and can also impose (or lift) restrictions on the
cannam@127 432 type of transform algorithm that is employed.
cannam@127 433
cannam@127 434 @emph{Important:} the planner overwrites the input array during
cannam@127 435 planning unless a saved plan (@pxref{Wisdom}) is available for that
cannam@127 436 problem, so you should initialize your input data after creating the
cannam@127 437 plan. The only exceptions to this are the @code{FFTW_ESTIMATE} and
cannam@127 438 @code{FFTW_WISDOM_ONLY} flags, as mentioned below.
cannam@127 439
cannam@127 440 In all cases, if wisdom is available for the given problem that was
cannam@127 441 created with equal-or-greater planning rigor, then the more rigorous
cannam@127 442 wisdom is used. For example, in @code{FFTW_ESTIMATE} mode any available
cannam@127 443 wisdom is used, whereas in @code{FFTW_PATIENT} mode only wisdom created
cannam@127 444 in patient or exhaustive mode can be used. @xref{Words of Wisdom-Saving
cannam@127 445 Plans}.
cannam@127 446
cannam@127 447 @subsubheading Planning-rigor flags
cannam@127 448 @itemize @bullet
cannam@127 449
cannam@127 450 @item
cannam@127 451 @ctindex FFTW_ESTIMATE
cannam@127 452 @code{FFTW_ESTIMATE} specifies that, instead of actual measurements of
cannam@127 453 different algorithms, a simple heuristic is used to pick a (probably
cannam@127 454 sub-optimal) plan quickly. With this flag, the input/output arrays are
cannam@127 455 not overwritten during planning.
cannam@127 456
cannam@127 457 @item
cannam@127 458 @ctindex FFTW_MEASURE
cannam@127 459 @code{FFTW_MEASURE} tells FFTW to find an optimized plan by actually
cannam@127 460 @emph{computing} several FFTs and measuring their execution time.
cannam@127 461 Depending on your machine, this can take some time (often a few
cannam@127 462 seconds). @code{FFTW_MEASURE} is the default planning option.
cannam@127 463
cannam@127 464 @item
cannam@127 465 @ctindex FFTW_PATIENT
cannam@127 466 @code{FFTW_PATIENT} is like @code{FFTW_MEASURE}, but considers a wider
cannam@127 467 range of algorithms and often produces a ``more optimal'' plan
cannam@127 468 (especially for large transforms), but at the expense of several times
cannam@127 469 longer planning time (especially for large transforms).
cannam@127 470
cannam@127 471 @item
cannam@127 472 @ctindex FFTW_EXHAUSTIVE
cannam@127 473 @code{FFTW_EXHAUSTIVE} is like @code{FFTW_PATIENT}, but considers an
cannam@127 474 even wider range of algorithms, including many that we think are
cannam@127 475 unlikely to be fast, to produce the most optimal plan but with a
cannam@127 476 substantially increased planning time.
cannam@127 477
cannam@127 478 @item
cannam@127 479 @ctindex FFTW_WISDOM_ONLY
cannam@127 480 @code{FFTW_WISDOM_ONLY} is a special planning mode in which the plan
cannam@127 481 is only created if wisdom is available for the given problem, and
cannam@127 482 otherwise a @code{NULL} plan is returned. This can be combined with
cannam@127 483 other flags, e.g. @samp{FFTW_WISDOM_ONLY | FFTW_PATIENT} creates a
cannam@127 484 plan only if wisdom is available that was created in
cannam@127 485 @code{FFTW_PATIENT} or @code{FFTW_EXHAUSTIVE} mode. The
cannam@127 486 @code{FFTW_WISDOM_ONLY} flag is intended for users who need to detect
cannam@127 487 whether wisdom is available; for example, if wisdom is not available
cannam@127 488 one may wish to allocate new arrays for planning so that user data is
cannam@127 489 not overwritten.
cannam@127 490
cannam@127 491 @end itemize
cannam@127 492
cannam@127 493 @subsubheading Algorithm-restriction flags
cannam@127 494 @itemize @bullet
cannam@127 495
cannam@127 496 @item
cannam@127 497 @ctindex FFTW_DESTROY_INPUT
cannam@127 498 @code{FFTW_DESTROY_INPUT} specifies that an out-of-place transform is
cannam@127 499 allowed to @emph{overwrite its input} array with arbitrary data; this
cannam@127 500 can sometimes allow more efficient algorithms to be employed.
cannam@127 501 @cindex out-of-place
cannam@127 502
cannam@127 503 @item
cannam@127 504 @ctindex FFTW_PRESERVE_INPUT
cannam@127 505 @code{FFTW_PRESERVE_INPUT} specifies that an out-of-place transform must
cannam@127 506 @emph{not change its input} array. This is ordinarily the
cannam@127 507 @emph{default}, except for c2r and hc2r (i.e. complex-to-real)
cannam@127 508 transforms for which @code{FFTW_DESTROY_INPUT} is the default. In the
cannam@127 509 latter cases, passing @code{FFTW_PRESERVE_INPUT} will attempt to use
cannam@127 510 algorithms that do not destroy the input, at the expense of worse
cannam@127 511 performance; for multi-dimensional c2r transforms, however, no
cannam@127 512 input-preserving algorithms are implemented and the planner will return
cannam@127 513 @code{NULL} if one is requested.
cannam@127 514 @cindex c2r
cannam@127 515 @cindex hc2r
cannam@127 516
cannam@127 517 @item
cannam@127 518 @ctindex FFTW_UNALIGNED
cannam@127 519 @cindex alignment
cannam@127 520 @findex fftw_malloc
cannam@127 521 @findex fftw_alignment_of
cannam@127 522 @code{FFTW_UNALIGNED} specifies that the algorithm may not impose any
cannam@127 523 unusual alignment requirements on the input/output arrays (i.e. no
cannam@127 524 SIMD may be used). This flag is normally @emph{not necessary}, since
cannam@127 525 the planner automatically detects misaligned arrays. The only use for
cannam@127 526 this flag is if you want to use the new-array execute interface to
cannam@127 527 execute a given plan on a different array that may not be aligned like
cannam@127 528 the original. (Using @code{fftw_malloc} makes this flag unnecessary
cannam@127 529 even then. You can also use @code{fftw_alignment_of} to detect
cannam@127 530 whether two arrays are equivalently aligned.)
cannam@127 531
cannam@127 532 @end itemize
cannam@127 533
cannam@127 534 @subsubheading Limiting planning time
cannam@127 535
cannam@127 536 @example
cannam@127 537 extern void fftw_set_timelimit(double seconds);
cannam@127 538 @end example
cannam@127 539 @findex fftw_set_timelimit
cannam@127 540
cannam@127 541 This function instructs FFTW to spend at most @code{seconds} seconds
cannam@127 542 (approximately) in the planner. If @code{seconds ==
cannam@127 543 FFTW_NO_TIMELIMIT} (the default value, which is negative), then
cannam@127 544 planning time is unbounded. Otherwise, FFTW plans with a
cannam@127 545 progressively wider range of algorithms until the the given time limit
cannam@127 546 is reached or the given range of algorithms is explored, returning the
cannam@127 547 best available plan.
cannam@127 548 @ctindex FFTW_NO_TIMELIMIT
cannam@127 549
cannam@127 550
cannam@127 551 For example, specifying @code{FFTW_PATIENT} first plans in
cannam@127 552 @code{FFTW_ESTIMATE} mode, then in @code{FFTW_MEASURE} mode, then
cannam@127 553 finally (time permitting) in @code{FFTW_PATIENT}. If
cannam@127 554 @code{FFTW_EXHAUSTIVE} is specified instead, the planner will further
cannam@127 555 progress to @code{FFTW_EXHAUSTIVE} mode.
cannam@127 556
cannam@127 557 Note that the @code{seconds} argument specifies only a rough limit; in
cannam@127 558 practice, the planner may use somewhat more time if the time limit is
cannam@127 559 reached when the planner is in the middle of an operation that cannot
cannam@127 560 be interrupted. At the very least, the planner will complete planning
cannam@127 561 in @code{FFTW_ESTIMATE} mode (which is thus equivalent to a time limit
cannam@127 562 of 0).
cannam@127 563
cannam@127 564
cannam@127 565 @c =========>
cannam@127 566 @node Real-data DFTs, Real-data DFT Array Format, Planner Flags, Basic Interface
cannam@127 567 @subsection Real-data DFTs
cannam@127 568
cannam@127 569 @example
cannam@127 570 fftw_plan fftw_plan_dft_r2c_1d(int n0,
cannam@127 571 double *in, fftw_complex *out,
cannam@127 572 unsigned flags);
cannam@127 573 fftw_plan fftw_plan_dft_r2c_2d(int n0, int n1,
cannam@127 574 double *in, fftw_complex *out,
cannam@127 575 unsigned flags);
cannam@127 576 fftw_plan fftw_plan_dft_r2c_3d(int n0, int n1, int n2,
cannam@127 577 double *in, fftw_complex *out,
cannam@127 578 unsigned flags);
cannam@127 579 fftw_plan fftw_plan_dft_r2c(int rank, const int *n,
cannam@127 580 double *in, fftw_complex *out,
cannam@127 581 unsigned flags);
cannam@127 582 @end example
cannam@127 583 @findex fftw_plan_dft_r2c_1d
cannam@127 584 @findex fftw_plan_dft_r2c_2d
cannam@127 585 @findex fftw_plan_dft_r2c_3d
cannam@127 586 @findex fftw_plan_dft_r2c
cannam@127 587 @cindex r2c
cannam@127 588
cannam@127 589 Plan a real-input/complex-output discrete Fourier transform (DFT) in
cannam@127 590 zero or more dimensions, returning an @code{fftw_plan} (@pxref{Using
cannam@127 591 Plans}).
cannam@127 592
cannam@127 593 Once you have created a plan for a certain transform type and
cannam@127 594 parameters, then creating another plan of the same type and parameters,
cannam@127 595 but for different arrays, is fast and shares constant data with the
cannam@127 596 first plan (if it still exists).
cannam@127 597
cannam@127 598 The planner returns @code{NULL} if the plan cannot be created. A
cannam@127 599 non-@code{NULL} plan is always returned by the basic interface unless
cannam@127 600 you are using a customized FFTW configuration supporting a restricted
cannam@127 601 set of transforms, or if you use the @code{FFTW_PRESERVE_INPUT} flag
cannam@127 602 with a multi-dimensional out-of-place c2r transform (see below).
cannam@127 603
cannam@127 604 @subsubheading Arguments
cannam@127 605 @itemize @bullet
cannam@127 606
cannam@127 607 @item
cannam@127 608 @code{rank} is the rank of the transform (it should be the size of the
cannam@127 609 array @code{*n}), and can be any non-negative integer. (@xref{Complex
cannam@127 610 Multi-Dimensional DFTs}, for the definition of ``rank''.) The
cannam@127 611 @samp{_1d}, @samp{_2d}, and @samp{_3d} planners correspond to a
cannam@127 612 @code{rank} of @code{1}, @code{2}, and @code{3}, respectively. The rank
cannam@127 613 may be zero, which is equivalent to a rank-1 transform of size 1, i.e. a
cannam@127 614 copy of one real number (with zero imaginary part) from input to output.
cannam@127 615
cannam@127 616 @item
cannam@127 617 @code{n0}, @code{n1}, @code{n2}, or @code{n[0..rank-1]}, (as appropriate
cannam@127 618 for each routine) specify the size of the transform dimensions. They
cannam@127 619 can be any positive integer. This is different in general from the
cannam@127 620 @emph{physical} array dimensions, which are described in @ref{Real-data
cannam@127 621 DFT Array Format}.
cannam@127 622
cannam@127 623 @itemize @minus
cannam@127 624 @item
cannam@127 625 FFTW is best at handling sizes of the form
cannam@127 626 @ifinfo
cannam@127 627 @math{2^a 3^b 5^c 7^d 11^e 13^f},
cannam@127 628 @end ifinfo
cannam@127 629 @tex
cannam@127 630 $2^a 3^b 5^c 7^d 11^e 13^f$,
cannam@127 631 @end tex
cannam@127 632 @html
cannam@127 633 2<sup>a</sup> 3<sup>b</sup> 5<sup>c</sup> 7<sup>d</sup>
cannam@127 634 11<sup>e</sup> 13<sup>f</sup>,
cannam@127 635 @end html
cannam@127 636 where @math{e+f} is either @math{0} or @math{1}, and the other exponents
cannam@127 637 are arbitrary. Other sizes are computed by means of a slow,
cannam@127 638 general-purpose algorithm (which nevertheless retains @Onlogn{} performance even for prime sizes). (It is possible to customize FFTW
cannam@127 639 for different array sizes; see @ref{Installation and Customization}.)
cannam@127 640 Transforms whose sizes are powers of @math{2} are especially fast, and
cannam@127 641 it is generally beneficial for the @emph{last} dimension of an r2c/c2r
cannam@127 642 transform to be @emph{even}.
cannam@127 643 @end itemize
cannam@127 644
cannam@127 645 @item
cannam@127 646 @code{in} and @code{out} point to the input and output arrays of the
cannam@127 647 transform, which may be the same (yielding an in-place transform).
cannam@127 648 @cindex in-place
cannam@127 649 These arrays are overwritten during planning, unless
cannam@127 650 @code{FFTW_ESTIMATE} is used in the flags. (The arrays need not be
cannam@127 651 initialized, but they must be allocated.) For an in-place transform, it
cannam@127 652 is important to remember that the real array will require padding,
cannam@127 653 described in @ref{Real-data DFT Array Format}.
cannam@127 654 @cindex padding
cannam@127 655
cannam@127 656 @item
cannam@127 657 @cindex flags
cannam@127 658 @code{flags} is a bitwise OR (@samp{|}) of zero or more planner flags,
cannam@127 659 as defined in @ref{Planner Flags}.
cannam@127 660
cannam@127 661 @end itemize
cannam@127 662
cannam@127 663 The inverse transforms, taking complex input (storing the non-redundant
cannam@127 664 half of a logically Hermitian array) to real output, are given by:
cannam@127 665
cannam@127 666 @example
cannam@127 667 fftw_plan fftw_plan_dft_c2r_1d(int n0,
cannam@127 668 fftw_complex *in, double *out,
cannam@127 669 unsigned flags);
cannam@127 670 fftw_plan fftw_plan_dft_c2r_2d(int n0, int n1,
cannam@127 671 fftw_complex *in, double *out,
cannam@127 672 unsigned flags);
cannam@127 673 fftw_plan fftw_plan_dft_c2r_3d(int n0, int n1, int n2,
cannam@127 674 fftw_complex *in, double *out,
cannam@127 675 unsigned flags);
cannam@127 676 fftw_plan fftw_plan_dft_c2r(int rank, const int *n,
cannam@127 677 fftw_complex *in, double *out,
cannam@127 678 unsigned flags);
cannam@127 679 @end example
cannam@127 680 @findex fftw_plan_dft_c2r_1d
cannam@127 681 @findex fftw_plan_dft_c2r_2d
cannam@127 682 @findex fftw_plan_dft_c2r_3d
cannam@127 683 @findex fftw_plan_dft_c2r
cannam@127 684 @cindex c2r
cannam@127 685
cannam@127 686 The arguments are the same as for the r2c transforms, except that the
cannam@127 687 input and output data formats are reversed.
cannam@127 688
cannam@127 689 FFTW computes an unnormalized transform: computing an r2c followed by a
cannam@127 690 c2r transform (or vice versa) will result in the original data
cannam@127 691 multiplied by the size of the transform (the product of the logical
cannam@127 692 dimensions).
cannam@127 693 @cindex normalization
cannam@127 694 An r2c transform produces the same output as a @code{FFTW_FORWARD}
cannam@127 695 complex DFT of the same input, and a c2r transform is correspondingly
cannam@127 696 equivalent to @code{FFTW_BACKWARD}. For more information, see @ref{What
cannam@127 697 FFTW Really Computes}.
cannam@127 698
cannam@127 699 @c =========>
cannam@127 700 @node Real-data DFT Array Format, Real-to-Real Transforms, Real-data DFTs, Basic Interface
cannam@127 701 @subsection Real-data DFT Array Format
cannam@127 702 @cindex r2c/c2r multi-dimensional array format
cannam@127 703
cannam@127 704 The output of a DFT of real data (r2c) contains symmetries that, in
cannam@127 705 principle, make half of the outputs redundant (@pxref{What FFTW Really
cannam@127 706 Computes}). (Similarly for the input of an inverse c2r transform.) In
cannam@127 707 practice, it is not possible to entirely realize these savings in an
cannam@127 708 efficient and understandable format that generalizes to
cannam@127 709 multi-dimensional transforms. Instead, the output of the r2c
cannam@127 710 transforms is @emph{slightly} over half of the output of the
cannam@127 711 corresponding complex transform. We do not ``pack'' the data in any
cannam@127 712 way, but store it as an ordinary array of @code{fftw_complex} values.
cannam@127 713 In fact, this data is simply a subsection of what would be the array in
cannam@127 714 the corresponding complex transform.
cannam@127 715
cannam@127 716 Specifically, for a real transform of @math{d} (= @code{rank})
cannam@127 717 dimensions @ndims{}, the complex data is an @ndimshalf array of
cannam@127 718 @code{fftw_complex} values in row-major order (with the division rounded
cannam@127 719 down). That is, we only store the @emph{lower} half (non-negative
cannam@127 720 frequencies), plus one element, of the last dimension of the data from
cannam@127 721 the ordinary complex transform. (We could have instead taken half of
cannam@127 722 any other dimension, but implementation turns out to be simpler if the
cannam@127 723 last, contiguous, dimension is used.)
cannam@127 724
cannam@127 725 @cindex out-of-place
cannam@127 726 For an out-of-place transform, the real data is simply an array with
cannam@127 727 physical dimensions @ndims in row-major order.
cannam@127 728
cannam@127 729 @cindex in-place
cannam@127 730 @cindex padding
cannam@127 731 For an in-place transform, some complications arise since the complex data
cannam@127 732 is slightly larger than the real data. In this case, the final
cannam@127 733 dimension of the real data must be @emph{padded} with extra values to
cannam@127 734 accommodate the size of the complex data---two extra if the last
cannam@127 735 dimension is even and one if it is odd. That is, the last dimension of
cannam@127 736 the real data must physically contain
cannam@127 737 @tex
cannam@127 738 $2 (n_{d-1}/2+1)$
cannam@127 739 @end tex
cannam@127 740 @ifinfo
cannam@127 741 2 * (n[d-1]/2+1)
cannam@127 742 @end ifinfo
cannam@127 743 @html
cannam@127 744 2 * (n<sub>d-1</sub>/2+1)
cannam@127 745 @end html
cannam@127 746 @code{double} values (exactly enough to hold the complex data). This
cannam@127 747 physical array size does not, however, change the @emph{logical} array
cannam@127 748 size---only
cannam@127 749 @tex
cannam@127 750 $n_{d-1}$
cannam@127 751 @end tex
cannam@127 752 @ifinfo
cannam@127 753 n[d-1]
cannam@127 754 @end ifinfo
cannam@127 755 @html
cannam@127 756 n<sub>d-1</sub>
cannam@127 757 @end html
cannam@127 758 values are actually stored in the last dimension, and
cannam@127 759 @tex
cannam@127 760 $n_{d-1}$
cannam@127 761 @end tex
cannam@127 762 @ifinfo
cannam@127 763 n[d-1]
cannam@127 764 @end ifinfo
cannam@127 765 @html
cannam@127 766 n<sub>d-1</sub>
cannam@127 767 @end html
cannam@127 768 is the last dimension passed to the planner.
cannam@127 769
cannam@127 770 @c =========>
cannam@127 771 @node Real-to-Real Transforms, Real-to-Real Transform Kinds, Real-data DFT Array Format, Basic Interface
cannam@127 772 @subsection Real-to-Real Transforms
cannam@127 773 @cindex r2r
cannam@127 774
cannam@127 775 @example
cannam@127 776 fftw_plan fftw_plan_r2r_1d(int n, double *in, double *out,
cannam@127 777 fftw_r2r_kind kind, unsigned flags);
cannam@127 778 fftw_plan fftw_plan_r2r_2d(int n0, int n1, double *in, double *out,
cannam@127 779 fftw_r2r_kind kind0, fftw_r2r_kind kind1,
cannam@127 780 unsigned flags);
cannam@127 781 fftw_plan fftw_plan_r2r_3d(int n0, int n1, int n2,
cannam@127 782 double *in, double *out,
cannam@127 783 fftw_r2r_kind kind0,
cannam@127 784 fftw_r2r_kind kind1,
cannam@127 785 fftw_r2r_kind kind2,
cannam@127 786 unsigned flags);
cannam@127 787 fftw_plan fftw_plan_r2r(int rank, const int *n, double *in, double *out,
cannam@127 788 const fftw_r2r_kind *kind, unsigned flags);
cannam@127 789 @end example
cannam@127 790 @findex fftw_plan_r2r_1d
cannam@127 791 @findex fftw_plan_r2r_2d
cannam@127 792 @findex fftw_plan_r2r_3d
cannam@127 793 @findex fftw_plan_r2r
cannam@127 794
cannam@127 795 Plan a real input/output (r2r) transform of various kinds in zero or
cannam@127 796 more dimensions, returning an @code{fftw_plan} (@pxref{Using Plans}).
cannam@127 797
cannam@127 798 Once you have created a plan for a certain transform type and
cannam@127 799 parameters, then creating another plan of the same type and parameters,
cannam@127 800 but for different arrays, is fast and shares constant data with the
cannam@127 801 first plan (if it still exists).
cannam@127 802
cannam@127 803 The planner returns @code{NULL} if the plan cannot be created. A
cannam@127 804 non-@code{NULL} plan is always returned by the basic interface unless
cannam@127 805 you are using a customized FFTW configuration supporting a restricted
cannam@127 806 set of transforms, or for size-1 @code{FFTW_REDFT00} kinds (which are
cannam@127 807 not defined).
cannam@127 808 @ctindex FFTW_REDFT00
cannam@127 809
cannam@127 810 @subsubheading Arguments
cannam@127 811 @itemize @bullet
cannam@127 812
cannam@127 813 @item
cannam@127 814 @code{rank} is the dimensionality of the transform (it should be the
cannam@127 815 size of the arrays @code{*n} and @code{*kind}), and can be any
cannam@127 816 non-negative integer. The @samp{_1d}, @samp{_2d}, and @samp{_3d}
cannam@127 817 planners correspond to a @code{rank} of @code{1}, @code{2}, and
cannam@127 818 @code{3}, respectively. A @code{rank} of zero is equivalent to a copy
cannam@127 819 of one number from input to output.
cannam@127 820
cannam@127 821 @item
cannam@127 822 @code{n}, or @code{n0}/@code{n1}/@code{n2}, or @code{n[rank]},
cannam@127 823 respectively, gives the (physical) size of the transform dimensions.
cannam@127 824 They can be any positive integer.
cannam@127 825
cannam@127 826 @itemize @minus
cannam@127 827 @item
cannam@127 828 @cindex row-major
cannam@127 829 Multi-dimensional arrays are stored in row-major order with dimensions:
cannam@127 830 @code{n0} x @code{n1}; or @code{n0} x @code{n1} x @code{n2}; or
cannam@127 831 @code{n[0]} x @code{n[1]} x ... x @code{n[rank-1]}.
cannam@127 832 @xref{Multi-dimensional Array Format}.
cannam@127 833 @item
cannam@127 834 FFTW is generally best at handling sizes of the form
cannam@127 835 @ifinfo
cannam@127 836 @math{2^a 3^b 5^c 7^d 11^e 13^f},
cannam@127 837 @end ifinfo
cannam@127 838 @tex
cannam@127 839 $2^a 3^b 5^c 7^d 11^e 13^f$,
cannam@127 840 @end tex
cannam@127 841 @html
cannam@127 842 2<sup>a</sup> 3<sup>b</sup> 5<sup>c</sup> 7<sup>d</sup>
cannam@127 843 11<sup>e</sup> 13<sup>f</sup>,
cannam@127 844 @end html
cannam@127 845 where @math{e+f} is either @math{0} or @math{1}, and the other exponents
cannam@127 846 are arbitrary. Other sizes are computed by means of a slow,
cannam@127 847 general-purpose algorithm (which nevertheless retains @Onlogn{} performance even for prime sizes). (It is possible to customize FFTW
cannam@127 848 for different array sizes; see @ref{Installation and Customization}.)
cannam@127 849 Transforms whose sizes are powers of @math{2} are especially fast.
cannam@127 850 @item
cannam@127 851 For a @code{REDFT00} or @code{RODFT00} transform kind in a dimension of
cannam@127 852 size @math{n}, it is @math{n-1} or @math{n+1}, respectively, that
cannam@127 853 should be factorizable in the above form.
cannam@127 854 @end itemize
cannam@127 855
cannam@127 856 @item
cannam@127 857 @code{in} and @code{out} point to the input and output arrays of the
cannam@127 858 transform, which may be the same (yielding an in-place transform).
cannam@127 859 @cindex in-place
cannam@127 860 These arrays are overwritten during planning, unless
cannam@127 861 @code{FFTW_ESTIMATE} is used in the flags. (The arrays need not be
cannam@127 862 initialized, but they must be allocated.)
cannam@127 863
cannam@127 864 @item
cannam@127 865 @code{kind}, or @code{kind0}/@code{kind1}/@code{kind2}, or
cannam@127 866 @code{kind[rank]}, is the kind of r2r transform used for the
cannam@127 867 corresponding dimension. The valid kind constants are described in
cannam@127 868 @ref{Real-to-Real Transform Kinds}. In a multi-dimensional transform,
cannam@127 869 what is computed is the separable product formed by taking each
cannam@127 870 transform kind along the corresponding dimension, one dimension after
cannam@127 871 another.
cannam@127 872
cannam@127 873 @item
cannam@127 874 @cindex flags
cannam@127 875 @code{flags} is a bitwise OR (@samp{|}) of zero or more planner flags,
cannam@127 876 as defined in @ref{Planner Flags}.
cannam@127 877
cannam@127 878 @end itemize
cannam@127 879
cannam@127 880 @c =========>
cannam@127 881 @node Real-to-Real Transform Kinds, , Real-to-Real Transforms, Basic Interface
cannam@127 882 @subsection Real-to-Real Transform Kinds
cannam@127 883 @cindex kind (r2r)
cannam@127 884
cannam@127 885 FFTW currently supports 11 different r2r transform kinds, specified by
cannam@127 886 one of the constants below. For the precise definitions of these
cannam@127 887 transforms, see @ref{What FFTW Really Computes}. For a more colloquial
cannam@127 888 introduction to these transform kinds, see @ref{More DFTs of Real Data}.
cannam@127 889
cannam@127 890 For dimension of size @code{n}, there is a corresponding ``logical''
cannam@127 891 dimension @code{N} that determines the normalization (and the optimal
cannam@127 892 factorization); the formula for @code{N} is given for each kind below.
cannam@127 893 Also, with each transform kind is listed its corrsponding inverse
cannam@127 894 transform. FFTW computes unnormalized transforms: a transform followed
cannam@127 895 by its inverse will result in the original data multiplied by @code{N}
cannam@127 896 (or the product of the @code{N}'s for each dimension, in
cannam@127 897 multi-dimensions).
cannam@127 898 @cindex normalization
cannam@127 899
cannam@127 900 @itemize @bullet
cannam@127 901
cannam@127 902 @item
cannam@127 903 @ctindex FFTW_R2HC
cannam@127 904 @code{FFTW_R2HC} computes a real-input DFT with output in
cannam@127 905 ``halfcomplex'' format, i.e. real and imaginary parts for a transform of
cannam@127 906 size @code{n} stored as:
cannam@127 907 @tex
cannam@127 908 $$
cannam@127 909 r_0, r_1, r_2, \ldots, r_{n/2}, i_{(n+1)/2-1}, \ldots, i_2, i_1
cannam@127 910 $$
cannam@127 911 @end tex
cannam@127 912 @ifinfo
cannam@127 913 r0, r1, r2, r(n/2), i((n+1)/2-1), ..., i2, i1
cannam@127 914 @end ifinfo
cannam@127 915 @html
cannam@127 916 <p align=center>
cannam@127 917 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@127 918 </p>
cannam@127 919 @end html
cannam@127 920 (Logical @code{N=n}, inverse is @code{FFTW_HC2R}.)
cannam@127 921
cannam@127 922 @item
cannam@127 923 @ctindex FFTW_HC2R
cannam@127 924 @code{FFTW_HC2R} computes the reverse of @code{FFTW_R2HC}, above.
cannam@127 925 (Logical @code{N=n}, inverse is @code{FFTW_R2HC}.)
cannam@127 926
cannam@127 927 @item
cannam@127 928 @ctindex FFTW_DHT
cannam@127 929 @code{FFTW_DHT} computes a discrete Hartley transform.
cannam@127 930 (Logical @code{N=n}, inverse is @code{FFTW_DHT}.)
cannam@127 931 @cindex discrete Hartley transform
cannam@127 932
cannam@127 933 @item
cannam@127 934 @ctindex FFTW_REDFT00
cannam@127 935 @code{FFTW_REDFT00} computes an REDFT00 transform, i.e. a DCT-I.
cannam@127 936 (Logical @code{N=2*(n-1)}, inverse is @code{FFTW_REDFT00}.)
cannam@127 937 @cindex discrete cosine transform
cannam@127 938 @cindex DCT
cannam@127 939
cannam@127 940 @item
cannam@127 941 @ctindex FFTW_REDFT10
cannam@127 942 @code{FFTW_REDFT10} computes an REDFT10 transform, i.e. a DCT-II (sometimes called ``the'' DCT).
cannam@127 943 (Logical @code{N=2*n}, inverse is @code{FFTW_REDFT01}.)
cannam@127 944
cannam@127 945 @item
cannam@127 946 @ctindex FFTW_REDFT01
cannam@127 947 @code{FFTW_REDFT01} computes an REDFT01 transform, i.e. a DCT-III (sometimes called ``the'' IDCT, being the inverse of DCT-II).
cannam@127 948 (Logical @code{N=2*n}, inverse is @code{FFTW_REDFT=10}.)
cannam@127 949 @cindex IDCT
cannam@127 950
cannam@127 951 @item
cannam@127 952 @ctindex FFTW_REDFT11
cannam@127 953 @code{FFTW_REDFT11} computes an REDFT11 transform, i.e. a DCT-IV.
cannam@127 954 (Logical @code{N=2*n}, inverse is @code{FFTW_REDFT11}.)
cannam@127 955
cannam@127 956 @item
cannam@127 957 @ctindex FFTW_RODFT00
cannam@127 958 @code{FFTW_RODFT00} computes an RODFT00 transform, i.e. a DST-I.
cannam@127 959 (Logical @code{N=2*(n+1)}, inverse is @code{FFTW_RODFT00}.)
cannam@127 960 @cindex discrete sine transform
cannam@127 961 @cindex DST
cannam@127 962
cannam@127 963 @item
cannam@127 964 @ctindex FFTW_RODFT10
cannam@127 965 @code{FFTW_RODFT10} computes an RODFT10 transform, i.e. a DST-II.
cannam@127 966 (Logical @code{N=2*n}, inverse is @code{FFTW_RODFT01}.)
cannam@127 967
cannam@127 968 @item
cannam@127 969 @ctindex FFTW_RODFT01
cannam@127 970 @code{FFTW_RODFT01} computes an RODFT01 transform, i.e. a DST-III.
cannam@127 971 (Logical @code{N=2*n}, inverse is @code{FFTW_RODFT=10}.)
cannam@127 972
cannam@127 973 @item
cannam@127 974 @ctindex FFTW_RODFT11
cannam@127 975 @code{FFTW_RODFT11} computes an RODFT11 transform, i.e. a DST-IV.
cannam@127 976 (Logical @code{N=2*n}, inverse is @code{FFTW_RODFT11}.)
cannam@127 977
cannam@127 978 @end itemize
cannam@127 979
cannam@127 980 @c ------------------------------------------------------------
cannam@127 981 @node Advanced Interface, Guru Interface, Basic Interface, FFTW Reference
cannam@127 982 @section Advanced Interface
cannam@127 983 @cindex advanced interface
cannam@127 984
cannam@127 985 FFTW's ``advanced'' interface supplements the basic interface with four
cannam@127 986 new planner routines, providing a new level of flexibility: you can plan
cannam@127 987 a transform of multiple arrays simultaneously, operate on non-contiguous
cannam@127 988 (strided) data, and transform a subset of a larger multi-dimensional
cannam@127 989 array. Other than these additional features, the planner operates in
cannam@127 990 the same fashion as in the basic interface, and the resulting
cannam@127 991 @code{fftw_plan} is used in the same way (@pxref{Using Plans}).
cannam@127 992
cannam@127 993 @menu
cannam@127 994 * Advanced Complex DFTs::
cannam@127 995 * Advanced Real-data DFTs::
cannam@127 996 * Advanced Real-to-real Transforms::
cannam@127 997 @end menu
cannam@127 998
cannam@127 999 @c =========>
cannam@127 1000 @node Advanced Complex DFTs, Advanced Real-data DFTs, Advanced Interface, Advanced Interface
cannam@127 1001 @subsection Advanced Complex DFTs
cannam@127 1002
cannam@127 1003 @example
cannam@127 1004 fftw_plan fftw_plan_many_dft(int rank, const int *n, int howmany,
cannam@127 1005 fftw_complex *in, const int *inembed,
cannam@127 1006 int istride, int idist,
cannam@127 1007 fftw_complex *out, const int *onembed,
cannam@127 1008 int ostride, int odist,
cannam@127 1009 int sign, unsigned flags);
cannam@127 1010 @end example
cannam@127 1011 @findex fftw_plan_many_dft
cannam@127 1012
cannam@127 1013 This routine plans multiple multidimensional complex DFTs, and it
cannam@127 1014 extends the @code{fftw_plan_dft} routine (@pxref{Complex DFTs}) to
cannam@127 1015 compute @code{howmany} transforms, each having rank @code{rank} and size
cannam@127 1016 @code{n}. In addition, the transform data need not be contiguous, but
cannam@127 1017 it may be laid out in memory with an arbitrary stride. To account for
cannam@127 1018 these possibilities, @code{fftw_plan_many_dft} adds the new parameters
cannam@127 1019 @code{howmany}, @{@code{i},@code{o}@}@code{nembed},
cannam@127 1020 @{@code{i},@code{o}@}@code{stride}, and
cannam@127 1021 @{@code{i},@code{o}@}@code{dist}. The FFTW basic interface
cannam@127 1022 (@pxref{Complex DFTs}) provides routines specialized for ranks 1, 2,
cannam@127 1023 and@tie{}3, but the advanced interface handles only the general-rank
cannam@127 1024 case.
cannam@127 1025
cannam@127 1026 @code{howmany} is the number of transforms to compute. The resulting
cannam@127 1027 plan computes @code{howmany} transforms, where the input of the
cannam@127 1028 @code{k}-th transform is at location @code{in+k*idist} (in C pointer
cannam@127 1029 arithmetic), and its output is at location @code{out+k*odist}. Plans
cannam@127 1030 obtained in this way can often be faster than calling FFTW multiple
cannam@127 1031 times for the individual transforms. The basic @code{fftw_plan_dft}
cannam@127 1032 interface corresponds to @code{howmany=1} (in which case the @code{dist}
cannam@127 1033 parameters are ignored).
cannam@127 1034 @cindex howmany parameter
cannam@127 1035 @cindex dist
cannam@127 1036
cannam@127 1037
cannam@127 1038 Each of the @code{howmany} transforms has rank @code{rank} and size
cannam@127 1039 @code{n}, as in the basic interface. In addition, the advanced
cannam@127 1040 interface allows the input and output arrays of each transform to be
cannam@127 1041 row-major subarrays of larger rank-@code{rank} arrays, described by
cannam@127 1042 @code{inembed} and @code{onembed} parameters, respectively.
cannam@127 1043 @{@code{i},@code{o}@}@code{nembed} must be arrays of length @code{rank},
cannam@127 1044 and @code{n} should be elementwise less than or equal to
cannam@127 1045 @{@code{i},@code{o}@}@code{nembed}. Passing @code{NULL} for an
cannam@127 1046 @code{nembed} parameter is equivalent to passing @code{n} (i.e. same
cannam@127 1047 physical and logical dimensions, as in the basic interface.)
cannam@127 1048
cannam@127 1049 The @code{stride} parameters indicate that the @code{j}-th element of
cannam@127 1050 the input or output arrays is located at @code{j*istride} or
cannam@127 1051 @code{j*ostride}, respectively. (For a multi-dimensional array,
cannam@127 1052 @code{j} is the ordinary row-major index.) When combined with the
cannam@127 1053 @code{k}-th transform in a @code{howmany} loop, from above, this means
cannam@127 1054 that the (@code{j},@code{k})-th element is at @code{j*stride+k*dist}.
cannam@127 1055 (The basic @code{fftw_plan_dft} interface corresponds to a stride of 1.)
cannam@127 1056 @cindex stride
cannam@127 1057
cannam@127 1058
cannam@127 1059 For in-place transforms, the input and output @code{stride} and
cannam@127 1060 @code{dist} parameters should be the same; otherwise, the planner may
cannam@127 1061 return @code{NULL}.
cannam@127 1062
cannam@127 1063 Arrays @code{n}, @code{inembed}, and @code{onembed} are not used after
cannam@127 1064 this function returns. You can safely free or reuse them.
cannam@127 1065
cannam@127 1066 @strong{Examples}:
cannam@127 1067 One transform of one 5 by 6 array contiguous in memory:
cannam@127 1068 @example
cannam@127 1069 int rank = 2;
cannam@127 1070 int n[] = @{5, 6@};
cannam@127 1071 int howmany = 1;
cannam@127 1072 int idist = odist = 0; /* unused because howmany = 1 */
cannam@127 1073 int istride = ostride = 1; /* array is contiguous in memory */
cannam@127 1074 int *inembed = n, *onembed = n;
cannam@127 1075 @end example
cannam@127 1076
cannam@127 1077 Transform of three 5 by 6 arrays, each contiguous in memory,
cannam@127 1078 stored in memory one after another:
cannam@127 1079 @example
cannam@127 1080 int rank = 2;
cannam@127 1081 int n[] = @{5, 6@};
cannam@127 1082 int howmany = 3;
cannam@127 1083 int idist = odist = n[0]*n[1]; /* = 30, the distance in memory
cannam@127 1084 between the first element
cannam@127 1085 of the first array and the
cannam@127 1086 first element of the second array */
cannam@127 1087 int istride = ostride = 1; /* array is contiguous in memory */
cannam@127 1088 int *inembed = n, *onembed = n;
cannam@127 1089 @end example
cannam@127 1090
cannam@127 1091 Transform each column of a 2d array with 10 rows and 3 columns:
cannam@127 1092 @example
cannam@127 1093 int rank = 1; /* not 2: we are computing 1d transforms */
cannam@127 1094 int n[] = @{10@}; /* 1d transforms of length 10 */
cannam@127 1095 int howmany = 3;
cannam@127 1096 int idist = odist = 1;
cannam@127 1097 int istride = ostride = 3; /* distance between two elements in
cannam@127 1098 the same column */
cannam@127 1099 int *inembed = n, *onembed = n;
cannam@127 1100 @end example
cannam@127 1101
cannam@127 1102 @c =========>
cannam@127 1103 @node Advanced Real-data DFTs, Advanced Real-to-real Transforms, Advanced Complex DFTs, Advanced Interface
cannam@127 1104 @subsection Advanced Real-data DFTs
cannam@127 1105
cannam@127 1106 @example
cannam@127 1107 fftw_plan fftw_plan_many_dft_r2c(int rank, const int *n, int howmany,
cannam@127 1108 double *in, const int *inembed,
cannam@127 1109 int istride, int idist,
cannam@127 1110 fftw_complex *out, const int *onembed,
cannam@127 1111 int ostride, int odist,
cannam@127 1112 unsigned flags);
cannam@127 1113 fftw_plan fftw_plan_many_dft_c2r(int rank, const int *n, int howmany,
cannam@127 1114 fftw_complex *in, const int *inembed,
cannam@127 1115 int istride, int idist,
cannam@127 1116 double *out, const int *onembed,
cannam@127 1117 int ostride, int odist,
cannam@127 1118 unsigned flags);
cannam@127 1119 @end example
cannam@127 1120 @findex fftw_plan_many_dft_r2c
cannam@127 1121 @findex fftw_plan_many_dft_c2r
cannam@127 1122
cannam@127 1123 Like @code{fftw_plan_many_dft}, these two functions add @code{howmany},
cannam@127 1124 @code{nembed}, @code{stride}, and @code{dist} parameters to the
cannam@127 1125 @code{fftw_plan_dft_r2c} and @code{fftw_plan_dft_c2r} functions, but
cannam@127 1126 otherwise behave the same as the basic interface.
cannam@127 1127
cannam@127 1128 The interpretation of @code{howmany}, @code{stride}, and @code{dist} are
cannam@127 1129 the same as for @code{fftw_plan_many_dft}, above. Note that the
cannam@127 1130 @code{stride} and @code{dist} for the real array are in units of
cannam@127 1131 @code{double}, and for the complex array are in units of
cannam@127 1132 @code{fftw_complex}.
cannam@127 1133
cannam@127 1134 If an @code{nembed} parameter is @code{NULL}, it is interpreted as what
cannam@127 1135 it would be in the basic interface, as described in @ref{Real-data DFT
cannam@127 1136 Array Format}. That is, for the complex array the size is assumed to be
cannam@127 1137 the same as @code{n}, but with the last dimension cut roughly in half.
cannam@127 1138 For the real array, the size is assumed to be @code{n} if the transform
cannam@127 1139 is out-of-place, or @code{n} with the last dimension ``padded'' if the
cannam@127 1140 transform is in-place.
cannam@127 1141
cannam@127 1142 If an @code{nembed} parameter is non-@code{NULL}, it is interpreted as
cannam@127 1143 the physical size of the corresponding array, in row-major order, just
cannam@127 1144 as for @code{fftw_plan_many_dft}. In this case, each dimension of
cannam@127 1145 @code{nembed} should be @code{>=} what it would be in the basic
cannam@127 1146 interface (e.g. the halved or padded @code{n}).
cannam@127 1147
cannam@127 1148 Arrays @code{n}, @code{inembed}, and @code{onembed} are not used after
cannam@127 1149 this function returns. You can safely free or reuse them.
cannam@127 1150
cannam@127 1151 @c =========>
cannam@127 1152 @node Advanced Real-to-real Transforms, , Advanced Real-data DFTs, Advanced Interface
cannam@127 1153 @subsection Advanced Real-to-real Transforms
cannam@127 1154
cannam@127 1155 @example
cannam@127 1156 fftw_plan fftw_plan_many_r2r(int rank, const int *n, int howmany,
cannam@127 1157 double *in, const int *inembed,
cannam@127 1158 int istride, int idist,
cannam@127 1159 double *out, const int *onembed,
cannam@127 1160 int ostride, int odist,
cannam@127 1161 const fftw_r2r_kind *kind, unsigned flags);
cannam@127 1162 @end example
cannam@127 1163 @findex fftw_plan_many_r2r
cannam@127 1164
cannam@127 1165 Like @code{fftw_plan_many_dft}, this functions adds @code{howmany},
cannam@127 1166 @code{nembed}, @code{stride}, and @code{dist} parameters to the
cannam@127 1167 @code{fftw_plan_r2r} function, but otherwise behave the same as the
cannam@127 1168 basic interface. The interpretation of those additional parameters are
cannam@127 1169 the same as for @code{fftw_plan_many_dft}. (Of course, the
cannam@127 1170 @code{stride} and @code{dist} parameters are now in units of
cannam@127 1171 @code{double}, not @code{fftw_complex}.)
cannam@127 1172
cannam@127 1173 Arrays @code{n}, @code{inembed}, @code{onembed}, and @code{kind} are not
cannam@127 1174 used after this function returns. You can safely free or reuse them.
cannam@127 1175
cannam@127 1176 @c ------------------------------------------------------------
cannam@127 1177 @node Guru Interface, New-array Execute Functions, Advanced Interface, FFTW Reference
cannam@127 1178 @section Guru Interface
cannam@127 1179 @cindex guru interface
cannam@127 1180
cannam@127 1181 The ``guru'' interface to FFTW is intended to expose as much as possible
cannam@127 1182 of the flexibility in the underlying FFTW architecture. It allows one
cannam@127 1183 to compute multi-dimensional ``vectors'' (loops) of multi-dimensional
cannam@127 1184 transforms, where each vector/transform dimension has an independent
cannam@127 1185 size and stride.
cannam@127 1186 @cindex vector
cannam@127 1187 One can also use more general complex-number formats, e.g. separate real
cannam@127 1188 and imaginary arrays.
cannam@127 1189
cannam@127 1190 For those users who require the flexibility of the guru interface, it is
cannam@127 1191 important that they pay special attention to the documentation lest they
cannam@127 1192 shoot themselves in the foot.
cannam@127 1193
cannam@127 1194 @menu
cannam@127 1195 * Interleaved and split arrays::
cannam@127 1196 * Guru vector and transform sizes::
cannam@127 1197 * Guru Complex DFTs::
cannam@127 1198 * Guru Real-data DFTs::
cannam@127 1199 * Guru Real-to-real Transforms::
cannam@127 1200 * 64-bit Guru Interface::
cannam@127 1201 @end menu
cannam@127 1202
cannam@127 1203 @c =========>
cannam@127 1204 @node Interleaved and split arrays, Guru vector and transform sizes, Guru Interface, Guru Interface
cannam@127 1205 @subsection Interleaved and split arrays
cannam@127 1206
cannam@127 1207 The guru interface supports two representations of complex numbers,
cannam@127 1208 which we call the interleaved and the split format.
cannam@127 1209
cannam@127 1210 The @dfn{interleaved} format is the same one used by the basic and
cannam@127 1211 advanced interfaces, and it is documented in @ref{Complex numbers}.
cannam@127 1212 In the interleaved format, you provide pointers to the real part of a
cannam@127 1213 complex number, and the imaginary part understood to be stored in the
cannam@127 1214 next memory location.
cannam@127 1215 @cindex interleaved format
cannam@127 1216
cannam@127 1217
cannam@127 1218 The @dfn{split} format allows separate pointers to the real and
cannam@127 1219 imaginary parts of a complex array.
cannam@127 1220 @cindex split format
cannam@127 1221
cannam@127 1222
cannam@127 1223 Technically, the interleaved format is redundant, because you can
cannam@127 1224 always express an interleaved array in terms of a split array with
cannam@127 1225 appropriate pointers and strides. On the other hand, the interleaved
cannam@127 1226 format is simpler to use, and it is common in practice. Hence, FFTW
cannam@127 1227 supports it as a special case.
cannam@127 1228
cannam@127 1229 @c =========>
cannam@127 1230 @node Guru vector and transform sizes, Guru Complex DFTs, Interleaved and split arrays, Guru Interface
cannam@127 1231 @subsection Guru vector and transform sizes
cannam@127 1232
cannam@127 1233 The guru interface introduces one basic new data structure,
cannam@127 1234 @code{fftw_iodim}, that is used to specify sizes and strides for
cannam@127 1235 multi-dimensional transforms and vectors:
cannam@127 1236
cannam@127 1237 @example
cannam@127 1238 typedef struct @{
cannam@127 1239 int n;
cannam@127 1240 int is;
cannam@127 1241 int os;
cannam@127 1242 @} fftw_iodim;
cannam@127 1243 @end example
cannam@127 1244 @tindex fftw_iodim
cannam@127 1245
cannam@127 1246 Here, @code{n} is the size of the dimension, and @code{is} and @code{os}
cannam@127 1247 are the strides of that dimension for the input and output arrays. (The
cannam@127 1248 stride is the separation of consecutive elements along this dimension.)
cannam@127 1249
cannam@127 1250 The meaning of the stride parameter depends on the type of the array
cannam@127 1251 that the stride refers to. @emph{If the array is interleaved complex,
cannam@127 1252 strides are expressed in units of complex numbers
cannam@127 1253 (@code{fftw_complex}). If the array is split complex or real, strides
cannam@127 1254 are expressed in units of real numbers (@code{double}).} This
cannam@127 1255 convention is consistent with the usual pointer arithmetic in the C
cannam@127 1256 language. An interleaved array is denoted by a pointer @code{p} to
cannam@127 1257 @code{fftw_complex}, so that @code{p+1} points to the next complex
cannam@127 1258 number. Split arrays are denoted by pointers to @code{double}, in
cannam@127 1259 which case pointer arithmetic operates in units of
cannam@127 1260 @code{sizeof(double)}.
cannam@127 1261 @cindex stride
cannam@127 1262
cannam@127 1263
cannam@127 1264 The guru planner interfaces all take a (@code{rank}, @code{dims[rank]})
cannam@127 1265 pair describing the transform size, and a (@code{howmany_rank},
cannam@127 1266 @code{howmany_dims[howmany_rank]}) pair describing the ``vector'' size (a
cannam@127 1267 multi-dimensional loop of transforms to perform), where @code{dims} and
cannam@127 1268 @code{howmany_dims} are arrays of @code{fftw_iodim}.
cannam@127 1269
cannam@127 1270 For example, the @code{howmany} parameter in the advanced complex-DFT
cannam@127 1271 interface corresponds to @code{howmany_rank} = 1,
cannam@127 1272 @code{howmany_dims[0].n} = @code{howmany}, @code{howmany_dims[0].is} =
cannam@127 1273 @code{idist}, and @code{howmany_dims[0].os} = @code{odist}.
cannam@127 1274 @cindex howmany loop
cannam@127 1275 @cindex dist
cannam@127 1276 (To compute a single transform, you can just use @code{howmany_rank} = 0.)
cannam@127 1277
cannam@127 1278
cannam@127 1279 A row-major multidimensional array with dimensions @code{n[rank]}
cannam@127 1280 (@pxref{Row-major Format}) corresponds to @code{dims[i].n} =
cannam@127 1281 @code{n[i]} and the recurrence @code{dims[i].is} = @code{n[i+1] *
cannam@127 1282 dims[i+1].is} (similarly for @code{os}). The stride of the last
cannam@127 1283 (@code{i=rank-1}) dimension is the overall stride of the array.
cannam@127 1284 e.g. to be equivalent to the advanced complex-DFT interface, you would
cannam@127 1285 have @code{dims[rank-1].is} = @code{istride} and
cannam@127 1286 @code{dims[rank-1].os} = @code{ostride}.
cannam@127 1287 @cindex row-major
cannam@127 1288
cannam@127 1289
cannam@127 1290 In general, we only guarantee FFTW to return a non-@code{NULL} plan if
cannam@127 1291 the vector and transform dimensions correspond to a set of distinct
cannam@127 1292 indices, and for in-place transforms the input/output strides should
cannam@127 1293 be the same.
cannam@127 1294
cannam@127 1295 @c =========>
cannam@127 1296 @node Guru Complex DFTs, Guru Real-data DFTs, Guru vector and transform sizes, Guru Interface
cannam@127 1297 @subsection Guru Complex DFTs
cannam@127 1298
cannam@127 1299 @example
cannam@127 1300 fftw_plan fftw_plan_guru_dft(
cannam@127 1301 int rank, const fftw_iodim *dims,
cannam@127 1302 int howmany_rank, const fftw_iodim *howmany_dims,
cannam@127 1303 fftw_complex *in, fftw_complex *out,
cannam@127 1304 int sign, unsigned flags);
cannam@127 1305
cannam@127 1306 fftw_plan fftw_plan_guru_split_dft(
cannam@127 1307 int rank, const fftw_iodim *dims,
cannam@127 1308 int howmany_rank, const fftw_iodim *howmany_dims,
cannam@127 1309 double *ri, double *ii, double *ro, double *io,
cannam@127 1310 unsigned flags);
cannam@127 1311 @end example
cannam@127 1312 @findex fftw_plan_guru_dft
cannam@127 1313 @findex fftw_plan_guru_split_dft
cannam@127 1314
cannam@127 1315 These two functions plan a complex-data, multi-dimensional DFT
cannam@127 1316 for the interleaved and split format, respectively.
cannam@127 1317 Transform dimensions are given by (@code{rank}, @code{dims}) over a
cannam@127 1318 multi-dimensional vector (loop) of dimensions (@code{howmany_rank},
cannam@127 1319 @code{howmany_dims}). @code{dims} and @code{howmany_dims} should point
cannam@127 1320 to @code{fftw_iodim} arrays of length @code{rank} and
cannam@127 1321 @code{howmany_rank}, respectively.
cannam@127 1322
cannam@127 1323 @cindex flags
cannam@127 1324 @code{flags} is a bitwise OR (@samp{|}) of zero or more planner flags,
cannam@127 1325 as defined in @ref{Planner Flags}.
cannam@127 1326
cannam@127 1327 In the @code{fftw_plan_guru_dft} function, the pointers @code{in} and
cannam@127 1328 @code{out} point to the interleaved input and output arrays,
cannam@127 1329 respectively. The sign can be either @math{-1} (=
cannam@127 1330 @code{FFTW_FORWARD}) or @math{+1} (= @code{FFTW_BACKWARD}). If the
cannam@127 1331 pointers are equal, the transform is in-place.
cannam@127 1332
cannam@127 1333 In the @code{fftw_plan_guru_split_dft} function,
cannam@127 1334 @code{ri} and @code{ii} point to the real and imaginary input arrays,
cannam@127 1335 and @code{ro} and @code{io} point to the real and imaginary output
cannam@127 1336 arrays. The input and output pointers may be the same, indicating an
cannam@127 1337 in-place transform. For example, for @code{fftw_complex} pointers
cannam@127 1338 @code{in} and @code{out}, the corresponding parameters are:
cannam@127 1339
cannam@127 1340 @example
cannam@127 1341 ri = (double *) in;
cannam@127 1342 ii = (double *) in + 1;
cannam@127 1343 ro = (double *) out;
cannam@127 1344 io = (double *) out + 1;
cannam@127 1345 @end example
cannam@127 1346
cannam@127 1347 Because @code{fftw_plan_guru_split_dft} accepts split arrays, strides
cannam@127 1348 are expressed in units of @code{double}. For a contiguous
cannam@127 1349 @code{fftw_complex} array, the overall stride of the transform should
cannam@127 1350 be 2, the distance between consecutive real parts or between
cannam@127 1351 consecutive imaginary parts; see @ref{Guru vector and transform
cannam@127 1352 sizes}. Note that the dimension strides are applied equally to the
cannam@127 1353 real and imaginary parts; real and imaginary arrays with different
cannam@127 1354 strides are not supported.
cannam@127 1355
cannam@127 1356 There is no @code{sign} parameter in @code{fftw_plan_guru_split_dft}.
cannam@127 1357 This function always plans for an @code{FFTW_FORWARD} transform. To
cannam@127 1358 plan for an @code{FFTW_BACKWARD} transform, you can exploit the
cannam@127 1359 identity that the backwards DFT is equal to the forwards DFT with the
cannam@127 1360 real and imaginary parts swapped. For example, in the case of the
cannam@127 1361 @code{fftw_complex} arrays above, the @code{FFTW_BACKWARD} transform
cannam@127 1362 is computed by the parameters:
cannam@127 1363
cannam@127 1364 @example
cannam@127 1365 ri = (double *) in + 1;
cannam@127 1366 ii = (double *) in;
cannam@127 1367 ro = (double *) out + 1;
cannam@127 1368 io = (double *) out;
cannam@127 1369 @end example
cannam@127 1370
cannam@127 1371 @c =========>
cannam@127 1372 @node Guru Real-data DFTs, Guru Real-to-real Transforms, Guru Complex DFTs, Guru Interface
cannam@127 1373 @subsection Guru Real-data DFTs
cannam@127 1374
cannam@127 1375 @example
cannam@127 1376 fftw_plan fftw_plan_guru_dft_r2c(
cannam@127 1377 int rank, const fftw_iodim *dims,
cannam@127 1378 int howmany_rank, const fftw_iodim *howmany_dims,
cannam@127 1379 double *in, fftw_complex *out,
cannam@127 1380 unsigned flags);
cannam@127 1381
cannam@127 1382 fftw_plan fftw_plan_guru_split_dft_r2c(
cannam@127 1383 int rank, const fftw_iodim *dims,
cannam@127 1384 int howmany_rank, const fftw_iodim *howmany_dims,
cannam@127 1385 double *in, double *ro, double *io,
cannam@127 1386 unsigned flags);
cannam@127 1387
cannam@127 1388 fftw_plan fftw_plan_guru_dft_c2r(
cannam@127 1389 int rank, const fftw_iodim *dims,
cannam@127 1390 int howmany_rank, const fftw_iodim *howmany_dims,
cannam@127 1391 fftw_complex *in, double *out,
cannam@127 1392 unsigned flags);
cannam@127 1393
cannam@127 1394 fftw_plan fftw_plan_guru_split_dft_c2r(
cannam@127 1395 int rank, const fftw_iodim *dims,
cannam@127 1396 int howmany_rank, const fftw_iodim *howmany_dims,
cannam@127 1397 double *ri, double *ii, double *out,
cannam@127 1398 unsigned flags);
cannam@127 1399 @end example
cannam@127 1400 @findex fftw_plan_guru_dft_r2c
cannam@127 1401 @findex fftw_plan_guru_split_dft_r2c
cannam@127 1402 @findex fftw_plan_guru_dft_c2r
cannam@127 1403 @findex fftw_plan_guru_split_dft_c2r
cannam@127 1404
cannam@127 1405 Plan a real-input (r2c) or real-output (c2r), multi-dimensional DFT with
cannam@127 1406 transform dimensions given by (@code{rank}, @code{dims}) over a
cannam@127 1407 multi-dimensional vector (loop) of dimensions (@code{howmany_rank},
cannam@127 1408 @code{howmany_dims}). @code{dims} and @code{howmany_dims} should point
cannam@127 1409 to @code{fftw_iodim} arrays of length @code{rank} and
cannam@127 1410 @code{howmany_rank}, respectively. As for the basic and advanced
cannam@127 1411 interfaces, an r2c transform is @code{FFTW_FORWARD} and a c2r transform
cannam@127 1412 is @code{FFTW_BACKWARD}.
cannam@127 1413
cannam@127 1414 The @emph{last} dimension of @code{dims} is interpreted specially:
cannam@127 1415 that dimension of the real array has size @code{dims[rank-1].n}, but
cannam@127 1416 that dimension of the complex array has size @code{dims[rank-1].n/2+1}
cannam@127 1417 (division rounded down). The strides, on the other hand, are taken to
cannam@127 1418 be exactly as specified. It is up to the user to specify the strides
cannam@127 1419 appropriately for the peculiar dimensions of the data, and we do not
cannam@127 1420 guarantee that the planner will succeed (return non-@code{NULL}) for
cannam@127 1421 any dimensions other than those described in @ref{Real-data DFT Array
cannam@127 1422 Format} and generalized in @ref{Advanced Real-data DFTs}. (That is,
cannam@127 1423 for an in-place transform, each individual dimension should be able to
cannam@127 1424 operate in place.)
cannam@127 1425 @cindex in-place
cannam@127 1426
cannam@127 1427
cannam@127 1428 @code{in} and @code{out} point to the input and output arrays for r2c
cannam@127 1429 and c2r transforms, respectively. For split arrays, @code{ri} and
cannam@127 1430 @code{ii} point to the real and imaginary input arrays for a c2r
cannam@127 1431 transform, and @code{ro} and @code{io} point to the real and imaginary
cannam@127 1432 output arrays for an r2c transform. @code{in} and @code{ro} or
cannam@127 1433 @code{ri} and @code{out} may be the same, indicating an in-place
cannam@127 1434 transform. (In-place transforms where @code{in} and @code{io} or
cannam@127 1435 @code{ii} and @code{out} are the same are not currently supported.)
cannam@127 1436
cannam@127 1437 @cindex flags
cannam@127 1438 @code{flags} is a bitwise OR (@samp{|}) of zero or more planner flags,
cannam@127 1439 as defined in @ref{Planner Flags}.
cannam@127 1440
cannam@127 1441 In-place transforms of rank greater than 1 are currently only
cannam@127 1442 supported for interleaved arrays. For split arrays, the planner will
cannam@127 1443 return @code{NULL}.
cannam@127 1444 @cindex in-place
cannam@127 1445
cannam@127 1446 @c =========>
cannam@127 1447 @node Guru Real-to-real Transforms, 64-bit Guru Interface, Guru Real-data DFTs, Guru Interface
cannam@127 1448 @subsection Guru Real-to-real Transforms
cannam@127 1449
cannam@127 1450 @example
cannam@127 1451 fftw_plan fftw_plan_guru_r2r(int rank, const fftw_iodim *dims,
cannam@127 1452 int howmany_rank,
cannam@127 1453 const fftw_iodim *howmany_dims,
cannam@127 1454 double *in, double *out,
cannam@127 1455 const fftw_r2r_kind *kind,
cannam@127 1456 unsigned flags);
cannam@127 1457 @end example
cannam@127 1458 @findex fftw_plan_guru_r2r
cannam@127 1459
cannam@127 1460 Plan a real-to-real (r2r) multi-dimensional @code{FFTW_FORWARD}
cannam@127 1461 transform with transform dimensions given by (@code{rank}, @code{dims})
cannam@127 1462 over a multi-dimensional vector (loop) of dimensions
cannam@127 1463 (@code{howmany_rank}, @code{howmany_dims}). @code{dims} and
cannam@127 1464 @code{howmany_dims} should point to @code{fftw_iodim} arrays of length
cannam@127 1465 @code{rank} and @code{howmany_rank}, respectively.
cannam@127 1466
cannam@127 1467 The transform kind of each dimension is given by the @code{kind}
cannam@127 1468 parameter, which should point to an array of length @code{rank}. Valid
cannam@127 1469 @code{fftw_r2r_kind} constants are given in @ref{Real-to-Real Transform
cannam@127 1470 Kinds}.
cannam@127 1471
cannam@127 1472 @code{in} and @code{out} point to the real input and output arrays; they
cannam@127 1473 may be the same, indicating an in-place transform.
cannam@127 1474
cannam@127 1475 @cindex flags
cannam@127 1476 @code{flags} is a bitwise OR (@samp{|}) of zero or more planner flags,
cannam@127 1477 as defined in @ref{Planner Flags}.
cannam@127 1478
cannam@127 1479 @c =========>
cannam@127 1480 @node 64-bit Guru Interface, , Guru Real-to-real Transforms, Guru Interface
cannam@127 1481 @subsection 64-bit Guru Interface
cannam@127 1482 @cindex 64-bit architecture
cannam@127 1483
cannam@127 1484 When compiled in 64-bit mode on a 64-bit architecture (where addresses
cannam@127 1485 are 64 bits wide), FFTW uses 64-bit quantities internally for all
cannam@127 1486 transform sizes, strides, and so on---you don't have to do anything
cannam@127 1487 special to exploit this. However, in the ordinary FFTW interfaces,
cannam@127 1488 you specify the transform size by an @code{int} quantity, which is
cannam@127 1489 normally only 32 bits wide. This means that, even though FFTW is
cannam@127 1490 using 64-bit sizes internally, you cannot specify a single transform
cannam@127 1491 dimension larger than
cannam@127 1492 @ifinfo
cannam@127 1493 2^31-1
cannam@127 1494 @end ifinfo
cannam@127 1495 @html
cannam@127 1496 2<sup><small>31</small></sup>&minus;1
cannam@127 1497 @end html
cannam@127 1498 @tex
cannam@127 1499 $2^{31}-1$
cannam@127 1500 @end tex
cannam@127 1501 numbers.
cannam@127 1502
cannam@127 1503 We expect that few users will require transforms larger than this, but,
cannam@127 1504 for those who do, we provide a 64-bit version of the guru interface in
cannam@127 1505 which all sizes are specified as integers of type @code{ptrdiff_t}
cannam@127 1506 instead of @code{int}. (@code{ptrdiff_t} is a signed integer type
cannam@127 1507 defined by the C standard to be wide enough to represent address
cannam@127 1508 differences, and thus must be at least 64 bits wide on a 64-bit
cannam@127 1509 machine.) We stress that there is @emph{no performance advantage} to
cannam@127 1510 using this interface---the same internal FFTW code is employed
cannam@127 1511 regardless---and it is only necessary if you want to specify very
cannam@127 1512 large transform sizes.
cannam@127 1513 @tindex ptrdiff_t
cannam@127 1514
cannam@127 1515
cannam@127 1516 In particular, the 64-bit guru interface is a set of planner routines
cannam@127 1517 that are exactly the same as the guru planner routines, except that
cannam@127 1518 they are named with @samp{guru64} instead of @samp{guru} and they take
cannam@127 1519 arguments of type @code{fftw_iodim64} instead of @code{fftw_iodim}.
cannam@127 1520 For example, instead of @code{fftw_plan_guru_dft}, we have
cannam@127 1521 @code{fftw_plan_guru64_dft}.
cannam@127 1522
cannam@127 1523 @example
cannam@127 1524 fftw_plan fftw_plan_guru64_dft(
cannam@127 1525 int rank, const fftw_iodim64 *dims,
cannam@127 1526 int howmany_rank, const fftw_iodim64 *howmany_dims,
cannam@127 1527 fftw_complex *in, fftw_complex *out,
cannam@127 1528 int sign, unsigned flags);
cannam@127 1529 @end example
cannam@127 1530 @findex fftw_plan_guru64_dft
cannam@127 1531
cannam@127 1532 The @code{fftw_iodim64} type is similar to @code{fftw_iodim}, with the
cannam@127 1533 same interpretation, except that it uses type @code{ptrdiff_t} instead
cannam@127 1534 of type @code{int}.
cannam@127 1535
cannam@127 1536 @example
cannam@127 1537 typedef struct @{
cannam@127 1538 ptrdiff_t n;
cannam@127 1539 ptrdiff_t is;
cannam@127 1540 ptrdiff_t os;
cannam@127 1541 @} fftw_iodim64;
cannam@127 1542 @end example
cannam@127 1543 @tindex fftw_iodim64
cannam@127 1544
cannam@127 1545 Every other @samp{fftw_plan_guru} function also has a
cannam@127 1546 @samp{fftw_plan_guru64} equivalent, but we do not repeat their
cannam@127 1547 documentation here since they are identical to the 32-bit versions
cannam@127 1548 except as noted above.
cannam@127 1549
cannam@127 1550 @c -----------------------------------------------------------
cannam@127 1551 @node New-array Execute Functions, Wisdom, Guru Interface, FFTW Reference
cannam@127 1552 @section New-array Execute Functions
cannam@127 1553 @cindex execute
cannam@127 1554 @cindex new-array execution
cannam@127 1555
cannam@127 1556 Normally, one executes a plan for the arrays with which the plan was
cannam@127 1557 created, by calling @code{fftw_execute(plan)} as described in @ref{Using
cannam@127 1558 Plans}.
cannam@127 1559 @findex fftw_execute
cannam@127 1560 However, it is possible for sophisticated users to apply a given plan
cannam@127 1561 to a @emph{different} array using the ``new-array execute'' functions
cannam@127 1562 detailed below, provided that the following conditions are met:
cannam@127 1563
cannam@127 1564 @itemize @bullet
cannam@127 1565
cannam@127 1566 @item
cannam@127 1567 The array size, strides, etcetera are the same (since those are set by
cannam@127 1568 the plan).
cannam@127 1569
cannam@127 1570 @item
cannam@127 1571 The input and output arrays are the same (in-place) or different
cannam@127 1572 (out-of-place) if the plan was originally created to be in-place or
cannam@127 1573 out-of-place, respectively.
cannam@127 1574
cannam@127 1575 @item
cannam@127 1576 For split arrays, the separations between the real and imaginary
cannam@127 1577 parts, @code{ii-ri} and @code{io-ro}, are the same as they were for
cannam@127 1578 the input and output arrays when the plan was created. (This
cannam@127 1579 condition is automatically satisfied for interleaved arrays.)
cannam@127 1580
cannam@127 1581 @item
cannam@127 1582 The @dfn{alignment} of the new input/output arrays is the same as that
cannam@127 1583 of the input/output arrays when the plan was created, unless the plan
cannam@127 1584 was created with the @code{FFTW_UNALIGNED} flag.
cannam@127 1585 @ctindex FFTW_UNALIGNED
cannam@127 1586 Here, the alignment is a platform-dependent quantity (for example, it is
cannam@127 1587 the address modulo 16 if SSE SIMD instructions are used, but the address
cannam@127 1588 modulo 4 for non-SIMD single-precision FFTW on the same machine). In
cannam@127 1589 general, only arrays allocated with @code{fftw_malloc} are guaranteed to
cannam@127 1590 be equally aligned (@pxref{SIMD alignment and fftw_malloc}).
cannam@127 1591
cannam@127 1592 @end itemize
cannam@127 1593
cannam@127 1594 @cindex alignment
cannam@127 1595 The alignment issue is especially critical, because if you don't use
cannam@127 1596 @code{fftw_malloc} then you may have little control over the alignment
cannam@127 1597 of arrays in memory. For example, neither the C++ @code{new} function
cannam@127 1598 nor the Fortran @code{allocate} statement provide strong enough
cannam@127 1599 guarantees about data alignment. If you don't use @code{fftw_malloc},
cannam@127 1600 therefore, you probably have to use @code{FFTW_UNALIGNED} (which
cannam@127 1601 disables most SIMD support). If possible, it is probably better for
cannam@127 1602 you to simply create multiple plans (creating a new plan is quick once
cannam@127 1603 one exists for a given size), or better yet re-use the same array for
cannam@127 1604 your transforms.
cannam@127 1605
cannam@127 1606 @findex fftw_alignment_of
cannam@127 1607 For rare circumstances in which you cannot control the alignment of
cannam@127 1608 allocated memory, but wish to determine where a given array is
cannam@127 1609 aligned like the original array for which a plan was created, you can
cannam@127 1610 use the @code{fftw_alignment_of} function:
cannam@127 1611 @example
cannam@127 1612 int fftw_alignment_of(double *p);
cannam@127 1613 @end example
cannam@127 1614 Two arrays have equivalent alignment (for the purposes of applying a
cannam@127 1615 plan) if and only if @code{fftw_alignment_of} returns the same value
cannam@127 1616 for the corresponding pointers to their data (typecast to @code{double*}
cannam@127 1617 if necessary).
cannam@127 1618
cannam@127 1619 If you are tempted to use the new-array execute interface because you
cannam@127 1620 want to transform a known bunch of arrays of the same size, you should
cannam@127 1621 probably go use the advanced interface instead (@pxref{Advanced
cannam@127 1622 Interface})).
cannam@127 1623
cannam@127 1624 The new-array execute functions are:
cannam@127 1625
cannam@127 1626 @example
cannam@127 1627 void fftw_execute_dft(
cannam@127 1628 const fftw_plan p,
cannam@127 1629 fftw_complex *in, fftw_complex *out);
cannam@127 1630
cannam@127 1631 void fftw_execute_split_dft(
cannam@127 1632 const fftw_plan p,
cannam@127 1633 double *ri, double *ii, double *ro, double *io);
cannam@127 1634
cannam@127 1635 void fftw_execute_dft_r2c(
cannam@127 1636 const fftw_plan p,
cannam@127 1637 double *in, fftw_complex *out);
cannam@127 1638
cannam@127 1639 void fftw_execute_split_dft_r2c(
cannam@127 1640 const fftw_plan p,
cannam@127 1641 double *in, double *ro, double *io);
cannam@127 1642
cannam@127 1643 void fftw_execute_dft_c2r(
cannam@127 1644 const fftw_plan p,
cannam@127 1645 fftw_complex *in, double *out);
cannam@127 1646
cannam@127 1647 void fftw_execute_split_dft_c2r(
cannam@127 1648 const fftw_plan p,
cannam@127 1649 double *ri, double *ii, double *out);
cannam@127 1650
cannam@127 1651 void fftw_execute_r2r(
cannam@127 1652 const fftw_plan p,
cannam@127 1653 double *in, double *out);
cannam@127 1654 @end example
cannam@127 1655 @findex fftw_execute_dft
cannam@127 1656 @findex fftw_execute_split_dft
cannam@127 1657 @findex fftw_execute_dft_r2c
cannam@127 1658 @findex fftw_execute_split_dft_r2c
cannam@127 1659 @findex fftw_execute_dft_c2r
cannam@127 1660 @findex fftw_execute_split_dft_c2r
cannam@127 1661 @findex fftw_execute_r2r
cannam@127 1662
cannam@127 1663 These execute the @code{plan} to compute the corresponding transform on
cannam@127 1664 the input/output arrays specified by the subsequent arguments. The
cannam@127 1665 input/output array arguments have the same meanings as the ones passed
cannam@127 1666 to the guru planner routines in the preceding sections. The @code{plan}
cannam@127 1667 is not modified, and these routines can be called as many times as
cannam@127 1668 desired, or intermixed with calls to the ordinary @code{fftw_execute}.
cannam@127 1669
cannam@127 1670 The @code{plan} @emph{must} have been created for the transform type
cannam@127 1671 corresponding to the execute function, e.g. it must be a complex-DFT
cannam@127 1672 plan for @code{fftw_execute_dft}. Any of the planner routines for that
cannam@127 1673 transform type, from the basic to the guru interface, could have been
cannam@127 1674 used to create the plan, however.
cannam@127 1675
cannam@127 1676 @c ------------------------------------------------------------
cannam@127 1677 @node Wisdom, What FFTW Really Computes, New-array Execute Functions, FFTW Reference
cannam@127 1678 @section Wisdom
cannam@127 1679 @cindex wisdom
cannam@127 1680 @cindex saving plans to disk
cannam@127 1681
cannam@127 1682 This section documents the FFTW mechanism for saving and restoring
cannam@127 1683 plans from disk. This mechanism is called @dfn{wisdom}.
cannam@127 1684
cannam@127 1685 @menu
cannam@127 1686 * Wisdom Export::
cannam@127 1687 * Wisdom Import::
cannam@127 1688 * Forgetting Wisdom::
cannam@127 1689 * Wisdom Utilities::
cannam@127 1690 @end menu
cannam@127 1691
cannam@127 1692 @c =========>
cannam@127 1693 @node Wisdom Export, Wisdom Import, Wisdom, Wisdom
cannam@127 1694 @subsection Wisdom Export
cannam@127 1695
cannam@127 1696 @example
cannam@127 1697 int fftw_export_wisdom_to_filename(const char *filename);
cannam@127 1698 void fftw_export_wisdom_to_file(FILE *output_file);
cannam@127 1699 char *fftw_export_wisdom_to_string(void);
cannam@127 1700 void fftw_export_wisdom(void (*write_char)(char c, void *), void *data);
cannam@127 1701 @end example
cannam@127 1702 @findex fftw_export_wisdom
cannam@127 1703 @findex fftw_export_wisdom_to_filename
cannam@127 1704 @findex fftw_export_wisdom_to_file
cannam@127 1705 @findex fftw_export_wisdom_to_string
cannam@127 1706
cannam@127 1707 These functions allow you to export all currently accumulated wisdom
cannam@127 1708 in a form from which it can be later imported and restored, even
cannam@127 1709 during a separate run of the program. (@xref{Words of Wisdom-Saving
cannam@127 1710 Plans}.) The current store of wisdom is not affected by calling any
cannam@127 1711 of these routines.
cannam@127 1712
cannam@127 1713 @code{fftw_export_wisdom} exports the wisdom to any output
cannam@127 1714 medium, as specified by the callback function
cannam@127 1715 @code{write_char}. @code{write_char} is a @code{putc}-like function that
cannam@127 1716 writes the character @code{c} to some output; its second parameter is
cannam@127 1717 the @code{data} pointer passed to @code{fftw_export_wisdom}. For
cannam@127 1718 convenience, the following three ``wrapper'' routines are provided:
cannam@127 1719
cannam@127 1720 @code{fftw_export_wisdom_to_filename} writes wisdom to a file named
cannam@127 1721 @code{filename} (which is created or overwritten), returning @code{1}
cannam@127 1722 on success and @code{0} on failure. A lower-level function, which
cannam@127 1723 requires you to open and close the file yourself (e.g. if you want to
cannam@127 1724 write wisdom to a portion of a larger file) is
cannam@127 1725 @code{fftw_export_wisdom_to_file}. This writes the wisdom to the
cannam@127 1726 current position in @code{output_file}, which should be open with
cannam@127 1727 write permission; upon exit, the file remains open and is positioned
cannam@127 1728 at the end of the wisdom data.
cannam@127 1729
cannam@127 1730 @code{fftw_export_wisdom_to_string} returns a pointer to a
cannam@127 1731 @code{NULL}-terminated string holding the wisdom data. This string is
cannam@127 1732 dynamically allocated, and it is the responsibility of the caller to
cannam@127 1733 deallocate it with @code{free} when it is no longer needed.
cannam@127 1734
cannam@127 1735 All of these routines export the wisdom in the same format, which we
cannam@127 1736 will not document here except to say that it is LISP-like ASCII text
cannam@127 1737 that is insensitive to white space.
cannam@127 1738
cannam@127 1739 @c =========>
cannam@127 1740 @node Wisdom Import, Forgetting Wisdom, Wisdom Export, Wisdom
cannam@127 1741 @subsection Wisdom Import
cannam@127 1742
cannam@127 1743 @example
cannam@127 1744 int fftw_import_system_wisdom(void);
cannam@127 1745 int fftw_import_wisdom_from_filename(const char *filename);
cannam@127 1746 int fftw_import_wisdom_from_string(const char *input_string);
cannam@127 1747 int fftw_import_wisdom(int (*read_char)(void *), void *data);
cannam@127 1748 @end example
cannam@127 1749 @findex fftw_import_wisdom
cannam@127 1750 @findex fftw_import_system_wisdom
cannam@127 1751 @findex fftw_import_wisdom_from_filename
cannam@127 1752 @findex fftw_import_wisdom_from_file
cannam@127 1753 @findex fftw_import_wisdom_from_string
cannam@127 1754
cannam@127 1755 These functions import wisdom into a program from data stored by the
cannam@127 1756 @code{fftw_export_wisdom} functions above. (@xref{Words of
cannam@127 1757 Wisdom-Saving Plans}.) The imported wisdom replaces any wisdom
cannam@127 1758 already accumulated by the running program.
cannam@127 1759
cannam@127 1760 @code{fftw_import_wisdom} imports wisdom from any input medium, as
cannam@127 1761 specified by the callback function @code{read_char}. @code{read_char} is
cannam@127 1762 a @code{getc}-like function that returns the next character in the
cannam@127 1763 input; its parameter is the @code{data} pointer passed to
cannam@127 1764 @code{fftw_import_wisdom}. If the end of the input data is reached
cannam@127 1765 (which should never happen for valid data), @code{read_char} should
cannam@127 1766 return @code{EOF} (as defined in @code{<stdio.h>}). For convenience,
cannam@127 1767 the following three ``wrapper'' routines are provided:
cannam@127 1768
cannam@127 1769 @code{fftw_import_wisdom_from_filename} reads wisdom from a file named
cannam@127 1770 @code{filename}. A lower-level function, which requires you to open
cannam@127 1771 and close the file yourself (e.g. if you want to read wisdom from a
cannam@127 1772 portion of a larger file) is @code{fftw_import_wisdom_from_file}. This
cannam@127 1773 reads wisdom from the current position in @code{input_file} (which
cannam@127 1774 should be open with read permission); upon exit, the file remains
cannam@127 1775 open, but the position of the read pointer is unspecified.
cannam@127 1776
cannam@127 1777 @code{fftw_import_wisdom_from_string} reads wisdom from the
cannam@127 1778 @code{NULL}-terminated string @code{input_string}.
cannam@127 1779
cannam@127 1780 @code{fftw_import_system_wisdom} reads wisdom from an
cannam@127 1781 implementation-defined standard file (@code{/etc/fftw/wisdom} on Unix
cannam@127 1782 and GNU systems).
cannam@127 1783 @cindex wisdom, system-wide
cannam@127 1784
cannam@127 1785
cannam@127 1786 The return value of these import routines is @code{1} if the wisdom was
cannam@127 1787 read successfully and @code{0} otherwise. Note that, in all of these
cannam@127 1788 functions, any data in the input stream past the end of the wisdom data
cannam@127 1789 is simply ignored.
cannam@127 1790
cannam@127 1791 @c =========>
cannam@127 1792 @node Forgetting Wisdom, Wisdom Utilities, Wisdom Import, Wisdom
cannam@127 1793 @subsection Forgetting Wisdom
cannam@127 1794
cannam@127 1795 @example
cannam@127 1796 void fftw_forget_wisdom(void);
cannam@127 1797 @end example
cannam@127 1798 @findex fftw_forget_wisdom
cannam@127 1799
cannam@127 1800 Calling @code{fftw_forget_wisdom} causes all accumulated @code{wisdom}
cannam@127 1801 to be discarded and its associated memory to be freed. (New
cannam@127 1802 @code{wisdom} can still be gathered subsequently, however.)
cannam@127 1803
cannam@127 1804 @c =========>
cannam@127 1805 @node Wisdom Utilities, , Forgetting Wisdom, Wisdom
cannam@127 1806 @subsection Wisdom Utilities
cannam@127 1807
cannam@127 1808 FFTW includes two standalone utility programs that deal with wisdom. We
cannam@127 1809 merely summarize them here, since they come with their own @code{man}
cannam@127 1810 pages for Unix and GNU systems (with HTML versions on our web site).
cannam@127 1811
cannam@127 1812 The first program is @code{fftw-wisdom} (or @code{fftwf-wisdom} in
cannam@127 1813 single precision, etcetera), which can be used to create a wisdom file
cannam@127 1814 containing plans for any of the transform sizes and types supported by
cannam@127 1815 FFTW. It is preferable to create wisdom directly from your executable
cannam@127 1816 (@pxref{Caveats in Using Wisdom}), but this program is useful for
cannam@127 1817 creating global wisdom files for @code{fftw_import_system_wisdom}.
cannam@127 1818 @cindex fftw-wisdom utility
cannam@127 1819
cannam@127 1820
cannam@127 1821 The second program is @code{fftw-wisdom-to-conf}, which takes a wisdom
cannam@127 1822 file as input and produces a @dfn{configuration routine} as output. The
cannam@127 1823 latter is a C subroutine that you can compile and link into your
cannam@127 1824 program, replacing a routine of the same name in the FFTW library, that
cannam@127 1825 determines which parts of FFTW are callable by your program.
cannam@127 1826 @code{fftw-wisdom-to-conf} produces a configuration routine that links
cannam@127 1827 to only those parts of FFTW needed by the saved plans in the wisdom,
cannam@127 1828 greatly reducing the size of statically linked executables (which should
cannam@127 1829 only attempt to create plans corresponding to those in the wisdom,
cannam@127 1830 however).
cannam@127 1831 @cindex fftw-wisdom-to-conf utility
cannam@127 1832 @cindex configuration routines
cannam@127 1833
cannam@127 1834 @c ------------------------------------------------------------
cannam@127 1835 @node What FFTW Really Computes, , Wisdom, FFTW Reference
cannam@127 1836 @section What FFTW Really Computes
cannam@127 1837
cannam@127 1838 In this section, we provide precise mathematical definitions for the
cannam@127 1839 transforms that FFTW computes. These transform definitions are fairly
cannam@127 1840 standard, but some authors follow slightly different conventions for the
cannam@127 1841 normalization of the transform (the constant factor in front) and the
cannam@127 1842 sign of the complex exponent. We begin by presenting the
cannam@127 1843 one-dimensional (1d) transform definitions, and then give the
cannam@127 1844 straightforward extension to multi-dimensional transforms.
cannam@127 1845
cannam@127 1846 @menu
cannam@127 1847 * The 1d Discrete Fourier Transform (DFT)::
cannam@127 1848 * The 1d Real-data DFT::
cannam@127 1849 * 1d Real-even DFTs (DCTs)::
cannam@127 1850 * 1d Real-odd DFTs (DSTs)::
cannam@127 1851 * 1d Discrete Hartley Transforms (DHTs)::
cannam@127 1852 * Multi-dimensional Transforms::
cannam@127 1853 @end menu
cannam@127 1854
cannam@127 1855 @c =========>
cannam@127 1856 @node The 1d Discrete Fourier Transform (DFT), The 1d Real-data DFT, What FFTW Really Computes, What FFTW Really Computes
cannam@127 1857 @subsection The 1d Discrete Fourier Transform (DFT)
cannam@127 1858
cannam@127 1859 @cindex discrete Fourier transform
cannam@127 1860 @cindex DFT
cannam@127 1861 The forward (@code{FFTW_FORWARD}) discrete Fourier transform (DFT) of a
cannam@127 1862 1d complex array @math{X} of size @math{n} computes an array @math{Y},
cannam@127 1863 where:
cannam@127 1864 @tex
cannam@127 1865 $$
cannam@127 1866 Y_k = \sum_{j = 0}^{n - 1} X_j e^{-2\pi j k \sqrt{-1}/n} \ .
cannam@127 1867 $$
cannam@127 1868 @end tex
cannam@127 1869 @ifinfo
cannam@127 1870 @center Y[k] = sum for j = 0 to (n - 1) of X[j] * exp(-2 pi j k sqrt(-1)/n) .
cannam@127 1871 @end ifinfo
cannam@127 1872 @html
cannam@127 1873 <center><img src="equation-dft.png" align="top">.</center>
cannam@127 1874 @end html
cannam@127 1875 The backward (@code{FFTW_BACKWARD}) DFT computes:
cannam@127 1876 @tex
cannam@127 1877 $$
cannam@127 1878 Y_k = \sum_{j = 0}^{n - 1} X_j e^{2\pi j k \sqrt{-1}/n} \ .
cannam@127 1879 $$
cannam@127 1880 @end tex
cannam@127 1881 @ifinfo
cannam@127 1882 @center Y[k] = sum for j = 0 to (n - 1) of X[j] * exp(2 pi j k sqrt(-1)/n) .
cannam@127 1883 @end ifinfo
cannam@127 1884 @html
cannam@127 1885 <center><img src="equation-idft.png" align="top">.</center>
cannam@127 1886 @end html
cannam@127 1887
cannam@127 1888 @cindex normalization
cannam@127 1889 FFTW computes an unnormalized transform, in that there is no coefficient
cannam@127 1890 in front of the summation in the DFT. In other words, applying the
cannam@127 1891 forward and then the backward transform will multiply the input by
cannam@127 1892 @math{n}.
cannam@127 1893
cannam@127 1894 @cindex frequency
cannam@127 1895 From above, an @code{FFTW_FORWARD} transform corresponds to a sign of
cannam@127 1896 @math{-1} in the exponent of the DFT. Note also that we use the
cannam@127 1897 standard ``in-order'' output ordering---the @math{k}-th output
cannam@127 1898 corresponds to the frequency @math{k/n} (or @math{k/T}, where @math{T}
cannam@127 1899 is your total sampling period). For those who like to think in terms of
cannam@127 1900 positive and negative frequencies, this means that the positive
cannam@127 1901 frequencies are stored in the first half of the output and the negative
cannam@127 1902 frequencies are stored in backwards order in the second half of the
cannam@127 1903 output. (The frequency @math{-k/n} is the same as the frequency
cannam@127 1904 @math{(n-k)/n}.)
cannam@127 1905
cannam@127 1906 @c =========>
cannam@127 1907 @node The 1d Real-data DFT, 1d Real-even DFTs (DCTs), The 1d Discrete Fourier Transform (DFT), What FFTW Really Computes
cannam@127 1908 @subsection The 1d Real-data DFT
cannam@127 1909
cannam@127 1910 The real-input (r2c) DFT in FFTW computes the @emph{forward} transform
cannam@127 1911 @math{Y} of the size @code{n} real array @math{X}, exactly as defined
cannam@127 1912 above, i.e.
cannam@127 1913 @tex
cannam@127 1914 $$
cannam@127 1915 Y_k = \sum_{j = 0}^{n - 1} X_j e^{-2\pi j k \sqrt{-1}/n} \ .
cannam@127 1916 $$
cannam@127 1917 @end tex
cannam@127 1918 @ifinfo
cannam@127 1919 @center Y[k] = sum for j = 0 to (n - 1) of X[j] * exp(-2 pi j k sqrt(-1)/n) .
cannam@127 1920 @end ifinfo
cannam@127 1921 @html
cannam@127 1922 <center><img src="equation-dft.png" align="top">.</center>
cannam@127 1923 @end html
cannam@127 1924 This output array @math{Y} can easily be shown to possess the
cannam@127 1925 ``Hermitian'' symmetry
cannam@127 1926 @cindex Hermitian
cannam@127 1927 @tex
cannam@127 1928 $Y_k = Y_{n-k}^*$,
cannam@127 1929 @end tex
cannam@127 1930 @ifinfo
cannam@127 1931 Y[k] = Y[n-k]*,
cannam@127 1932 @end ifinfo
cannam@127 1933 @html
cannam@127 1934 <i>Y<sub>k</sub> = Y<sub>n-k</sub></i><sup>*</sup>,
cannam@127 1935 @end html
cannam@127 1936 where we take @math{Y} to be periodic so that
cannam@127 1937 @tex
cannam@127 1938 $Y_n = Y_0$.
cannam@127 1939 @end tex
cannam@127 1940 @ifinfo
cannam@127 1941 Y[n] = Y[0].
cannam@127 1942 @end ifinfo
cannam@127 1943 @html
cannam@127 1944 <i>Y<sub>n</sub> = Y</i><sub>0</sub>.
cannam@127 1945 @end html
cannam@127 1946
cannam@127 1947 As a result of this symmetry, half of the output @math{Y} is redundant
cannam@127 1948 (being the complex conjugate of the other half), and so the 1d r2c
cannam@127 1949 transforms only output elements @math{0}@dots{}@math{n/2} of @math{Y}
cannam@127 1950 (@math{n/2+1} complex numbers), where the division by @math{2} is
cannam@127 1951 rounded down.
cannam@127 1952
cannam@127 1953 Moreover, the Hermitian symmetry implies that
cannam@127 1954 @tex
cannam@127 1955 $Y_0$
cannam@127 1956 @end tex
cannam@127 1957 @ifinfo
cannam@127 1958 Y[0]
cannam@127 1959 @end ifinfo
cannam@127 1960 @html
cannam@127 1961 <i>Y</i><sub>0</sub>
cannam@127 1962 @end html
cannam@127 1963 and, if @math{n} is even, the
cannam@127 1964 @tex
cannam@127 1965 $Y_{n/2}$
cannam@127 1966 @end tex
cannam@127 1967 @ifinfo
cannam@127 1968 Y[n/2]
cannam@127 1969 @end ifinfo
cannam@127 1970 @html
cannam@127 1971 <i>Y</i><sub><i>n</i>/2</sub>
cannam@127 1972 @end html
cannam@127 1973 element, are purely real. So, for the @code{R2HC} r2r transform, the
cannam@127 1974 halfcomplex format does not store the imaginary parts of these elements.
cannam@127 1975 @cindex r2r
cannam@127 1976 @ctindex R2HC
cannam@127 1977 @cindex halfcomplex format
cannam@127 1978
cannam@127 1979
cannam@127 1980 The c2r and @code{H2RC} r2r transforms compute the backward DFT of the
cannam@127 1981 @emph{complex} array @math{X} with Hermitian symmetry, stored in the
cannam@127 1982 r2c/@code{R2HC} output formats, respectively, where the backward
cannam@127 1983 transform is defined exactly as for the complex case:
cannam@127 1984 @tex
cannam@127 1985 $$
cannam@127 1986 Y_k = \sum_{j = 0}^{n - 1} X_j e^{2\pi j k \sqrt{-1}/n} \ .
cannam@127 1987 $$
cannam@127 1988 @end tex
cannam@127 1989 @ifinfo
cannam@127 1990 @center Y[k] = sum for j = 0 to (n - 1) of X[j] * exp(2 pi j k sqrt(-1)/n) .
cannam@127 1991 @end ifinfo
cannam@127 1992 @html
cannam@127 1993 <center><img src="equation-idft.png" align="top">.</center>
cannam@127 1994 @end html
cannam@127 1995 The outputs @code{Y} of this transform can easily be seen to be purely
cannam@127 1996 real, and are stored as an array of real numbers.
cannam@127 1997
cannam@127 1998 @cindex normalization
cannam@127 1999 Like FFTW's complex DFT, these transforms are unnormalized. In other
cannam@127 2000 words, applying the real-to-complex (forward) and then the
cannam@127 2001 complex-to-real (backward) transform will multiply the input by
cannam@127 2002 @math{n}.
cannam@127 2003
cannam@127 2004 @c =========>
cannam@127 2005 @node 1d Real-even DFTs (DCTs), 1d Real-odd DFTs (DSTs), The 1d Real-data DFT, What FFTW Really Computes
cannam@127 2006 @subsection 1d Real-even DFTs (DCTs)
cannam@127 2007
cannam@127 2008 The Real-even symmetry DFTs in FFTW are exactly equivalent to the unnormalized
cannam@127 2009 forward (and backward) DFTs as defined above, where the input array
cannam@127 2010 @math{X} of length @math{N} is purely real and is also @dfn{even} symmetry. In
cannam@127 2011 this case, the output array is likewise real and even symmetry.
cannam@127 2012 @cindex real-even DFT
cannam@127 2013 @cindex REDFT
cannam@127 2014
cannam@127 2015
cannam@127 2016 @ctindex REDFT00
cannam@127 2017 For the case of @code{REDFT00}, this even symmetry means that
cannam@127 2018 @tex
cannam@127 2019 $X_j = X_{N-j}$,
cannam@127 2020 @end tex
cannam@127 2021 @ifinfo
cannam@127 2022 X[j] = X[N-j],
cannam@127 2023 @end ifinfo
cannam@127 2024 @html
cannam@127 2025 <i>X<sub>j</sub> = X<sub>N-j</sub></i>,
cannam@127 2026 @end html
cannam@127 2027 where we take @math{X} to be periodic so that
cannam@127 2028 @tex
cannam@127 2029 $X_N = X_0$.
cannam@127 2030 @end tex
cannam@127 2031 @ifinfo
cannam@127 2032 X[N] = X[0].
cannam@127 2033 @end ifinfo
cannam@127 2034 @html
cannam@127 2035 <i>X<sub>N</sub> = X</i><sub>0</sub>.
cannam@127 2036 @end html
cannam@127 2037 Because of this redundancy, only the first @math{n} real numbers are
cannam@127 2038 actually stored, where @math{N = 2(n-1)}.
cannam@127 2039
cannam@127 2040 The proper definition of even symmetry for @code{REDFT10},
cannam@127 2041 @code{REDFT01}, and @code{REDFT11} transforms is somewhat more intricate
cannam@127 2042 because of the shifts by @math{1/2} of the input and/or output, although
cannam@127 2043 the corresponding boundary conditions are given in @ref{Real even/odd
cannam@127 2044 DFTs (cosine/sine transforms)}. Because of the even symmetry, however,
cannam@127 2045 the sine terms in the DFT all cancel and the remaining cosine terms are
cannam@127 2046 written explicitly below. This formulation often leads people to call
cannam@127 2047 such a transform a @dfn{discrete cosine transform} (DCT), although it is
cannam@127 2048 really just a special case of the DFT.
cannam@127 2049 @cindex discrete cosine transform
cannam@127 2050 @cindex DCT
cannam@127 2051
cannam@127 2052
cannam@127 2053 In each of the definitions below, we transform a real array @math{X} of
cannam@127 2054 length @math{n} to a real array @math{Y} of length @math{n}:
cannam@127 2055
cannam@127 2056 @subsubheading REDFT00 (DCT-I)
cannam@127 2057 @ctindex REDFT00
cannam@127 2058 An @code{REDFT00} transform (type-I DCT) in FFTW is defined by:
cannam@127 2059 @tex
cannam@127 2060 $$
cannam@127 2061 Y_k = X_0 + (-1)^k X_{n-1}
cannam@127 2062 + 2 \sum_{j=1}^{n-2} X_j \cos [ \pi j k / (n-1)].
cannam@127 2063 $$
cannam@127 2064 @end tex
cannam@127 2065 @ifinfo
cannam@127 2066 Y[k] = X[0] + (-1)^k X[n-1] + 2 (sum for j = 1 to n-2 of X[j] cos(pi jk /(n-1))).
cannam@127 2067 @end ifinfo
cannam@127 2068 @html
cannam@127 2069 <center><img src="equation-redft00.png" align="top">.</center>
cannam@127 2070 @end html
cannam@127 2071 Note that this transform is not defined for @math{n=1}. For @math{n=2},
cannam@127 2072 the summation term above is dropped as you might expect.
cannam@127 2073
cannam@127 2074 @subsubheading REDFT10 (DCT-II)
cannam@127 2075 @ctindex REDFT10
cannam@127 2076 An @code{REDFT10} transform (type-II DCT, sometimes called ``the'' DCT) in FFTW is defined by:
cannam@127 2077 @tex
cannam@127 2078 $$
cannam@127 2079 Y_k = 2 \sum_{j=0}^{n-1} X_j \cos [\pi (j+1/2) k / n].
cannam@127 2080 $$
cannam@127 2081 @end tex
cannam@127 2082 @ifinfo
cannam@127 2083 Y[k] = 2 (sum for j = 0 to n-1 of X[j] cos(pi (j+1/2) k / n)).
cannam@127 2084 @end ifinfo
cannam@127 2085 @html
cannam@127 2086 <center><img src="equation-redft10.png" align="top">.</center>
cannam@127 2087 @end html
cannam@127 2088
cannam@127 2089 @subsubheading REDFT01 (DCT-III)
cannam@127 2090 @ctindex REDFT01
cannam@127 2091 An @code{REDFT01} transform (type-III DCT) in FFTW is defined by:
cannam@127 2092 @tex
cannam@127 2093 $$
cannam@127 2094 Y_k = X_0 + 2 \sum_{j=1}^{n-1} X_j \cos [\pi j (k+1/2) / n].
cannam@127 2095 $$
cannam@127 2096 @end tex
cannam@127 2097 @ifinfo
cannam@127 2098 Y[k] = X[0] + 2 (sum for j = 1 to n-1 of X[j] cos(pi j (k+1/2) / n)).
cannam@127 2099 @end ifinfo
cannam@127 2100 @html
cannam@127 2101 <center><img src="equation-redft01.png" align="top">.</center>
cannam@127 2102 @end html
cannam@127 2103 In the case of @math{n=1}, this reduces to
cannam@127 2104 @tex
cannam@127 2105 $Y_0 = X_0$.
cannam@127 2106 @end tex
cannam@127 2107 @ifinfo
cannam@127 2108 Y[0] = X[0].
cannam@127 2109 @end ifinfo
cannam@127 2110 @html
cannam@127 2111 <i>Y</i><sub>0</sub> = <i>X</i><sub>0</sub>.
cannam@127 2112 @end html
cannam@127 2113 Up to a scale factor (see below), this is the inverse of @code{REDFT10} (``the'' DCT), and so the @code{REDFT01} (DCT-III) is sometimes called the ``IDCT''.
cannam@127 2114 @cindex IDCT
cannam@127 2115
cannam@127 2116 @subsubheading REDFT11 (DCT-IV)
cannam@127 2117 @ctindex REDFT11
cannam@127 2118 An @code{REDFT11} transform (type-IV DCT) in FFTW is defined by:
cannam@127 2119 @tex
cannam@127 2120 $$
cannam@127 2121 Y_k = 2 \sum_{j=0}^{n-1} X_j \cos [\pi (j+1/2) (k+1/2) / n].
cannam@127 2122 $$
cannam@127 2123 @end tex
cannam@127 2124 @ifinfo
cannam@127 2125 Y[k] = 2 (sum for j = 0 to n-1 of X[j] cos(pi (j+1/2) (k+1/2) / n)).
cannam@127 2126 @end ifinfo
cannam@127 2127 @html
cannam@127 2128 <center><img src="equation-redft11.png" align="top">.</center>
cannam@127 2129 @end html
cannam@127 2130
cannam@127 2131 @subsubheading Inverses and Normalization
cannam@127 2132
cannam@127 2133 These definitions correspond directly to the unnormalized DFTs used
cannam@127 2134 elsewhere in FFTW (hence the factors of @math{2} in front of the
cannam@127 2135 summations). The unnormalized inverse of @code{REDFT00} is
cannam@127 2136 @code{REDFT00}, of @code{REDFT10} is @code{REDFT01} and vice versa, and
cannam@127 2137 of @code{REDFT11} is @code{REDFT11}. Each unnormalized inverse results
cannam@127 2138 in the original array multiplied by @math{N}, where @math{N} is the
cannam@127 2139 @emph{logical} DFT size. For @code{REDFT00}, @math{N=2(n-1)} (note that
cannam@127 2140 @math{n=1} is not defined); otherwise, @math{N=2n}.
cannam@127 2141 @cindex normalization
cannam@127 2142
cannam@127 2143
cannam@127 2144 In defining the discrete cosine transform, some authors also include
cannam@127 2145 additional factors of
cannam@127 2146 @ifinfo
cannam@127 2147 sqrt(2)
cannam@127 2148 @end ifinfo
cannam@127 2149 @html
cannam@127 2150 &radic;2
cannam@127 2151 @end html
cannam@127 2152 @tex
cannam@127 2153 $\sqrt{2}$
cannam@127 2154 @end tex
cannam@127 2155 (or its inverse) multiplying selected inputs and/or outputs. This is a
cannam@127 2156 mostly cosmetic change that makes the transform orthogonal, but
cannam@127 2157 sacrifices the direct equivalence to a symmetric DFT.
cannam@127 2158
cannam@127 2159 @c =========>
cannam@127 2160 @node 1d Real-odd DFTs (DSTs), 1d Discrete Hartley Transforms (DHTs), 1d Real-even DFTs (DCTs), What FFTW Really Computes
cannam@127 2161 @subsection 1d Real-odd DFTs (DSTs)
cannam@127 2162
cannam@127 2163 The Real-odd symmetry DFTs in FFTW are exactly equivalent to the unnormalized
cannam@127 2164 forward (and backward) DFTs as defined above, where the input array
cannam@127 2165 @math{X} of length @math{N} is purely real and is also @dfn{odd} symmetry. In
cannam@127 2166 this case, the output is odd symmetry and purely imaginary.
cannam@127 2167 @cindex real-odd DFT
cannam@127 2168 @cindex RODFT
cannam@127 2169
cannam@127 2170
cannam@127 2171 @ctindex RODFT00
cannam@127 2172 For the case of @code{RODFT00}, this odd symmetry means that
cannam@127 2173 @tex
cannam@127 2174 $X_j = -X_{N-j}$,
cannam@127 2175 @end tex
cannam@127 2176 @ifinfo
cannam@127 2177 X[j] = -X[N-j],
cannam@127 2178 @end ifinfo
cannam@127 2179 @html
cannam@127 2180 <i>X<sub>j</sub> = -X<sub>N-j</sub></i>,
cannam@127 2181 @end html
cannam@127 2182 where we take @math{X} to be periodic so that
cannam@127 2183 @tex
cannam@127 2184 $X_N = X_0$.
cannam@127 2185 @end tex
cannam@127 2186 @ifinfo
cannam@127 2187 X[N] = X[0].
cannam@127 2188 @end ifinfo
cannam@127 2189 @html
cannam@127 2190 <i>X<sub>N</sub> = X</i><sub>0</sub>.
cannam@127 2191 @end html
cannam@127 2192 Because of this redundancy, only the first @math{n} real numbers
cannam@127 2193 starting at @math{j=1} are actually stored (the @math{j=0} element is
cannam@127 2194 zero), where @math{N = 2(n+1)}.
cannam@127 2195
cannam@127 2196 The proper definition of odd symmetry for @code{RODFT10},
cannam@127 2197 @code{RODFT01}, and @code{RODFT11} transforms is somewhat more intricate
cannam@127 2198 because of the shifts by @math{1/2} of the input and/or output, although
cannam@127 2199 the corresponding boundary conditions are given in @ref{Real even/odd
cannam@127 2200 DFTs (cosine/sine transforms)}. Because of the odd symmetry, however,
cannam@127 2201 the cosine terms in the DFT all cancel and the remaining sine terms are
cannam@127 2202 written explicitly below. This formulation often leads people to call
cannam@127 2203 such a transform a @dfn{discrete sine transform} (DST), although it is
cannam@127 2204 really just a special case of the DFT.
cannam@127 2205 @cindex discrete sine transform
cannam@127 2206 @cindex DST
cannam@127 2207
cannam@127 2208
cannam@127 2209 In each of the definitions below, we transform a real array @math{X} of
cannam@127 2210 length @math{n} to a real array @math{Y} of length @math{n}:
cannam@127 2211
cannam@127 2212 @subsubheading RODFT00 (DST-I)
cannam@127 2213 @ctindex RODFT00
cannam@127 2214 An @code{RODFT00} transform (type-I DST) in FFTW is defined by:
cannam@127 2215 @tex
cannam@127 2216 $$
cannam@127 2217 Y_k = 2 \sum_{j=0}^{n-1} X_j \sin [ \pi (j+1) (k+1) / (n+1)].
cannam@127 2218 $$
cannam@127 2219 @end tex
cannam@127 2220 @ifinfo
cannam@127 2221 Y[k] = 2 (sum for j = 0 to n-1 of X[j] sin(pi (j+1)(k+1) / (n+1))).
cannam@127 2222 @end ifinfo
cannam@127 2223 @html
cannam@127 2224 <center><img src="equation-rodft00.png" align="top">.</center>
cannam@127 2225 @end html
cannam@127 2226
cannam@127 2227 @subsubheading RODFT10 (DST-II)
cannam@127 2228 @ctindex RODFT10
cannam@127 2229 An @code{RODFT10} transform (type-II DST) in FFTW is defined by:
cannam@127 2230 @tex
cannam@127 2231 $$
cannam@127 2232 Y_k = 2 \sum_{j=0}^{n-1} X_j \sin [\pi (j+1/2) (k+1) / n].
cannam@127 2233 $$
cannam@127 2234 @end tex
cannam@127 2235 @ifinfo
cannam@127 2236 Y[k] = 2 (sum for j = 0 to n-1 of X[j] sin(pi (j+1/2) (k+1) / n)).
cannam@127 2237 @end ifinfo
cannam@127 2238 @html
cannam@127 2239 <center><img src="equation-rodft10.png" align="top">.</center>
cannam@127 2240 @end html
cannam@127 2241
cannam@127 2242 @subsubheading RODFT01 (DST-III)
cannam@127 2243 @ctindex RODFT01
cannam@127 2244 An @code{RODFT01} transform (type-III DST) in FFTW is defined by:
cannam@127 2245 @tex
cannam@127 2246 $$
cannam@127 2247 Y_k = (-1)^k X_{n-1} + 2 \sum_{j=0}^{n-2} X_j \sin [\pi (j+1) (k+1/2) / n].
cannam@127 2248 $$
cannam@127 2249 @end tex
cannam@127 2250 @ifinfo
cannam@127 2251 Y[k] = (-1)^k X[n-1] + 2 (sum for j = 0 to n-2 of X[j] sin(pi (j+1) (k+1/2) / n)).
cannam@127 2252 @end ifinfo
cannam@127 2253 @html
cannam@127 2254 <center><img src="equation-rodft01.png" align="top">.</center>
cannam@127 2255 @end html
cannam@127 2256 In the case of @math{n=1}, this reduces to
cannam@127 2257 @tex
cannam@127 2258 $Y_0 = X_0$.
cannam@127 2259 @end tex
cannam@127 2260 @ifinfo
cannam@127 2261 Y[0] = X[0].
cannam@127 2262 @end ifinfo
cannam@127 2263 @html
cannam@127 2264 <i>Y</i><sub>0</sub> = <i>X</i><sub>0</sub>.
cannam@127 2265 @end html
cannam@127 2266
cannam@127 2267 @subsubheading RODFT11 (DST-IV)
cannam@127 2268 @ctindex RODFT11
cannam@127 2269 An @code{RODFT11} transform (type-IV DST) in FFTW is defined by:
cannam@127 2270 @tex
cannam@127 2271 $$
cannam@127 2272 Y_k = 2 \sum_{j=0}^{n-1} X_j \sin [\pi (j+1/2) (k+1/2) / n].
cannam@127 2273 $$
cannam@127 2274 @end tex
cannam@127 2275 @ifinfo
cannam@127 2276 Y[k] = 2 (sum for j = 0 to n-1 of X[j] sin(pi (j+1/2) (k+1/2) / n)).
cannam@127 2277 @end ifinfo
cannam@127 2278 @html
cannam@127 2279 <center><img src="equation-rodft11.png" align="top">.</center>
cannam@127 2280 @end html
cannam@127 2281
cannam@127 2282 @subsubheading Inverses and Normalization
cannam@127 2283
cannam@127 2284 These definitions correspond directly to the unnormalized DFTs used
cannam@127 2285 elsewhere in FFTW (hence the factors of @math{2} in front of the
cannam@127 2286 summations). The unnormalized inverse of @code{RODFT00} is
cannam@127 2287 @code{RODFT00}, of @code{RODFT10} is @code{RODFT01} and vice versa, and
cannam@127 2288 of @code{RODFT11} is @code{RODFT11}. Each unnormalized inverse results
cannam@127 2289 in the original array multiplied by @math{N}, where @math{N} is the
cannam@127 2290 @emph{logical} DFT size. For @code{RODFT00}, @math{N=2(n+1)};
cannam@127 2291 otherwise, @math{N=2n}.
cannam@127 2292 @cindex normalization
cannam@127 2293
cannam@127 2294
cannam@127 2295 In defining the discrete sine transform, some authors also include
cannam@127 2296 additional factors of
cannam@127 2297 @ifinfo
cannam@127 2298 sqrt(2)
cannam@127 2299 @end ifinfo
cannam@127 2300 @html
cannam@127 2301 &radic;2
cannam@127 2302 @end html
cannam@127 2303 @tex
cannam@127 2304 $\sqrt{2}$
cannam@127 2305 @end tex
cannam@127 2306 (or its inverse) multiplying selected inputs and/or outputs. This is a
cannam@127 2307 mostly cosmetic change that makes the transform orthogonal, but
cannam@127 2308 sacrifices the direct equivalence to an antisymmetric DFT.
cannam@127 2309
cannam@127 2310 @c =========>
cannam@127 2311 @node 1d Discrete Hartley Transforms (DHTs), Multi-dimensional Transforms, 1d Real-odd DFTs (DSTs), What FFTW Really Computes
cannam@127 2312 @subsection 1d Discrete Hartley Transforms (DHTs)
cannam@127 2313
cannam@127 2314 @cindex discrete Hartley transform
cannam@127 2315 @cindex DHT
cannam@127 2316 The discrete Hartley transform (DHT) of a 1d real array @math{X} of size
cannam@127 2317 @math{n} computes a real array @math{Y} of the same size, where:
cannam@127 2318 @tex
cannam@127 2319 $$
cannam@127 2320 Y_k = \sum_{j = 0}^{n - 1} X_j [ \cos(2\pi j k / n) + \sin(2\pi j k / n)].
cannam@127 2321 $$
cannam@127 2322 @end tex
cannam@127 2323 @ifinfo
cannam@127 2324 @center Y[k] = sum for j = 0 to (n - 1) of X[j] * [cos(2 pi j k / n) + sin(2 pi j k / n)].
cannam@127 2325 @end ifinfo
cannam@127 2326 @html
cannam@127 2327 <center><img src="equation-dht.png" align="top">.</center>
cannam@127 2328 @end html
cannam@127 2329
cannam@127 2330 @cindex normalization
cannam@127 2331 FFTW computes an unnormalized transform, in that there is no coefficient
cannam@127 2332 in front of the summation in the DHT. In other words, applying the
cannam@127 2333 transform twice (the DHT is its own inverse) will multiply the input by
cannam@127 2334 @math{n}.
cannam@127 2335
cannam@127 2336 @c =========>
cannam@127 2337 @node Multi-dimensional Transforms, , 1d Discrete Hartley Transforms (DHTs), What FFTW Really Computes
cannam@127 2338 @subsection Multi-dimensional Transforms
cannam@127 2339
cannam@127 2340 The multi-dimensional transforms of FFTW, in general, compute simply the
cannam@127 2341 separable product of the given 1d transform along each dimension of the
cannam@127 2342 array. Since each of these transforms is unnormalized, computing the
cannam@127 2343 forward followed by the backward/inverse multi-dimensional transform
cannam@127 2344 will result in the original array scaled by the product of the
cannam@127 2345 normalization factors for each dimension (e.g. the product of the
cannam@127 2346 dimension sizes, for a multi-dimensional DFT).
cannam@127 2347
cannam@127 2348 @tex
cannam@127 2349 As an explicit example, consider the following exact mathematical
cannam@127 2350 definition of our multi-dimensional DFT. Let $X$ be a $d$-dimensional
cannam@127 2351 complex array whose elements are $X[j_1, j_2, \ldots, j_d]$, where $0
cannam@127 2352 \leq j_s < n_s$ for all~$s \in \{ 1, 2, \ldots, d \}$. Let also
cannam@127 2353 $\omega_s = e^{2\pi \sqrt{-1}/n_s}$, for all ~$s \in \{ 1, 2, \ldots, d
cannam@127 2354 \}$.
cannam@127 2355
cannam@127 2356 The forward transform computes a complex array~$Y$, whose
cannam@127 2357 structure is the same as that of~$X$, defined by
cannam@127 2358
cannam@127 2359 $$
cannam@127 2360 Y[k_1, k_2, \ldots, k_d] =
cannam@127 2361 \sum_{j_1 = 0}^{n_1 - 1}
cannam@127 2362 \sum_{j_2 = 0}^{n_2 - 1}
cannam@127 2363 \cdots
cannam@127 2364 \sum_{j_d = 0}^{n_d - 1}
cannam@127 2365 X[j_1, j_2, \ldots, j_d]
cannam@127 2366 \omega_1^{-j_1 k_1}
cannam@127 2367 \omega_2^{-j_2 k_2}
cannam@127 2368 \cdots
cannam@127 2369 \omega_d^{-j_d k_d} \ .
cannam@127 2370 $$
cannam@127 2371
cannam@127 2372 The backward transform computes
cannam@127 2373 $$
cannam@127 2374 Y[k_1, k_2, \ldots, k_d] =
cannam@127 2375 \sum_{j_1 = 0}^{n_1 - 1}
cannam@127 2376 \sum_{j_2 = 0}^{n_2 - 1}
cannam@127 2377 \cdots
cannam@127 2378 \sum_{j_d = 0}^{n_d - 1}
cannam@127 2379 X[j_1, j_2, \ldots, j_d]
cannam@127 2380 \omega_1^{j_1 k_1}
cannam@127 2381 \omega_2^{j_2 k_2}
cannam@127 2382 \cdots
cannam@127 2383 \omega_d^{j_d k_d} \ .
cannam@127 2384 $$
cannam@127 2385
cannam@127 2386 Computing the forward transform followed by the backward transform
cannam@127 2387 will multiply the array by $\prod_{s=1}^{d} n_d$.
cannam@127 2388 @end tex
cannam@127 2389
cannam@127 2390 @cindex r2c
cannam@127 2391 The definition of FFTW's multi-dimensional DFT of real data (r2c)
cannam@127 2392 deserves special attention. In this case, we logically compute the full
cannam@127 2393 multi-dimensional DFT of the input data; since the input data are purely
cannam@127 2394 real, the output data have the Hermitian symmetry and therefore only one
cannam@127 2395 non-redundant half need be stored. More specifically, for an @ndims multi-dimensional real-input DFT, the full (logical) complex output array
cannam@127 2396 @tex
cannam@127 2397 $Y[k_0, k_1, \ldots, k_{d-1}]$
cannam@127 2398 @end tex
cannam@127 2399 @html
cannam@127 2400 <i>Y</i>[<i>k</i><sub>0</sub>, <i>k</i><sub>1</sub>, ...,
cannam@127 2401 <i>k</i><sub><i>d-1</i></sub>]
cannam@127 2402 @end html
cannam@127 2403 @ifinfo
cannam@127 2404 Y[k[0], k[1], ..., k[d-1]]
cannam@127 2405 @end ifinfo
cannam@127 2406 has the symmetry:
cannam@127 2407 @tex
cannam@127 2408 $$
cannam@127 2409 Y[k_0, k_1, \ldots, k_{d-1}] = Y[n_0 - k_0, n_1 - k_1, \ldots, n_{d-1} - k_{d-1}]^*
cannam@127 2410 $$
cannam@127 2411 @end tex
cannam@127 2412 @html
cannam@127 2413 <i>Y</i>[<i>k</i><sub>0</sub>, <i>k</i><sub>1</sub>, ...,
cannam@127 2414 <i>k</i><sub><i>d-1</i></sub>] = <i>Y</i>[<i>n</i><sub>0</sub> -
cannam@127 2415 <i>k</i><sub>0</sub>, <i>n</i><sub>1</sub> - <i>k</i><sub>1</sub>, ...,
cannam@127 2416 <i>n</i><sub><i>d-1</i></sub> - <i>k</i><sub><i>d-1</i></sub>]<sup>*</sup>
cannam@127 2417 @end html
cannam@127 2418 @ifinfo
cannam@127 2419 Y[k[0], k[1], ..., k[d-1]] = Y[n[0] - k[0], n[1] - k[1], ..., n[d-1] - k[d-1]]*
cannam@127 2420 @end ifinfo
cannam@127 2421 (where each dimension is periodic). Because of this symmetry, we only
cannam@127 2422 store the
cannam@127 2423 @tex
cannam@127 2424 $k_{d-1} = 0 \cdots n_{d-1}/2$
cannam@127 2425 @end tex
cannam@127 2426 @html
cannam@127 2427 <i>k</i><sub><i>d-1</i></sub> = 0...<i>n</i><sub><i>d-1</i></sub>/2+1
cannam@127 2428 @end html
cannam@127 2429 @ifinfo
cannam@127 2430 k[d-1] = 0...n[d-1]/2
cannam@127 2431 @end ifinfo
cannam@127 2432 elements of the @emph{last} dimension (division by @math{2} is rounded
cannam@127 2433 down). (We could instead have cut any other dimension in half, but the
cannam@127 2434 last dimension proved computationally convenient.) This results in the
cannam@127 2435 peculiar array format described in more detail by @ref{Real-data DFT
cannam@127 2436 Array Format}.
cannam@127 2437
cannam@127 2438 The multi-dimensional c2r transform is simply the unnormalized inverse
cannam@127 2439 of the r2c transform. i.e. it is the same as FFTW's complex backward
cannam@127 2440 multi-dimensional DFT, operating on a Hermitian input array in the
cannam@127 2441 peculiar format mentioned above and outputting a real array (since the
cannam@127 2442 DFT output is purely real).
cannam@127 2443
cannam@127 2444 We should remind the user that the separable product of 1d transforms
cannam@127 2445 along each dimension, as computed by FFTW, is not always the same thing
cannam@127 2446 as the usual multi-dimensional transform. A multi-dimensional
cannam@127 2447 @code{R2HC} (or @code{HC2R}) transform is not identical to the
cannam@127 2448 multi-dimensional DFT, requiring some post-processing to combine the
cannam@127 2449 requisite real and imaginary parts, as was described in @ref{The
cannam@127 2450 Halfcomplex-format DFT}. Likewise, FFTW's multidimensional
cannam@127 2451 @code{FFTW_DHT} r2r transform is not the same thing as the logical
cannam@127 2452 multi-dimensional discrete Hartley transform defined in the literature,
cannam@127 2453 as discussed in @ref{The Discrete Hartley Transform}.
cannam@127 2454