Chris@42: Chris@42: Chris@42: Chris@42: Chris@42:
Chris@42:Chris@42: Next: Wisdom, Previous: Guru Interface, Up: FFTW Reference [Contents][Index]
Chris@42:Normally, one executes a plan for the arrays with which the plan was
Chris@42: created, by calling fftw_execute(plan)
as described in Using Plans.
Chris@42:
Chris@42: However, it is possible for sophisticated users to apply a given plan
Chris@42: to a different array using the “new-array execute” functions
Chris@42: detailed below, provided that the following conditions are met:
Chris@42:
ii-ri
and io-ro
, are the same as they were for
Chris@42: the input and output arrays when the plan was created. (This
Chris@42: condition is automatically satisfied for interleaved arrays.)
Chris@42:
Chris@42: FFTW_UNALIGNED
flag.
Chris@42:
Chris@42: Here, the alignment is a platform-dependent quantity (for example, it is
Chris@42: the address modulo 16 if SSE SIMD instructions are used, but the address
Chris@42: modulo 4 for non-SIMD single-precision FFTW on the same machine). In
Chris@42: general, only arrays allocated with fftw_malloc
are guaranteed to
Chris@42: be equally aligned (see SIMD alignment and fftw_malloc).
Chris@42:
Chris@42: The alignment issue is especially critical, because if you don’t use
Chris@42: fftw_malloc
then you may have little control over the alignment
Chris@42: of arrays in memory. For example, neither the C++ new
function
Chris@42: nor the Fortran allocate
statement provide strong enough
Chris@42: guarantees about data alignment. If you don’t use fftw_malloc
,
Chris@42: therefore, you probably have to use FFTW_UNALIGNED
(which
Chris@42: disables most SIMD support). If possible, it is probably better for
Chris@42: you to simply create multiple plans (creating a new plan is quick once
Chris@42: one exists for a given size), or better yet re-use the same array for
Chris@42: your transforms.
Chris@42:
For rare circumstances in which you cannot control the alignment of
Chris@42: allocated memory, but wish to determine where a given array is
Chris@42: aligned like the original array for which a plan was created, you can
Chris@42: use the fftw_alignment_of
function:
Chris@42:
int fftw_alignment_of(double *p); Chris@42:
Two arrays have equivalent alignment (for the purposes of applying a
Chris@42: plan) if and only if fftw_alignment_of
returns the same value
Chris@42: for the corresponding pointers to their data (typecast to double*
Chris@42: if necessary).
Chris@42:
If you are tempted to use the new-array execute interface because you Chris@42: want to transform a known bunch of arrays of the same size, you should Chris@42: probably go use the advanced interface instead (see Advanced Interface)). Chris@42:
Chris@42:The new-array execute functions are: Chris@42:
Chris@42:void fftw_execute_dft( Chris@42: const fftw_plan p, Chris@42: fftw_complex *in, fftw_complex *out); Chris@42: Chris@42: void fftw_execute_split_dft( Chris@42: const fftw_plan p, Chris@42: double *ri, double *ii, double *ro, double *io); Chris@42: Chris@42: void fftw_execute_dft_r2c( Chris@42: const fftw_plan p, Chris@42: double *in, fftw_complex *out); Chris@42: Chris@42: void fftw_execute_split_dft_r2c( Chris@42: const fftw_plan p, Chris@42: double *in, double *ro, double *io); Chris@42: Chris@42: void fftw_execute_dft_c2r( Chris@42: const fftw_plan p, Chris@42: fftw_complex *in, double *out); Chris@42: Chris@42: void fftw_execute_split_dft_c2r( Chris@42: const fftw_plan p, Chris@42: double *ri, double *ii, double *out); Chris@42: Chris@42: void fftw_execute_r2r( Chris@42: const fftw_plan p, Chris@42: double *in, double *out); Chris@42:
These execute the plan
to compute the corresponding transform on
Chris@42: the input/output arrays specified by the subsequent arguments. The
Chris@42: input/output array arguments have the same meanings as the ones passed
Chris@42: to the guru planner routines in the preceding sections. The plan
Chris@42: is not modified, and these routines can be called as many times as
Chris@42: desired, or intermixed with calls to the ordinary fftw_execute
.
Chris@42:
The plan
must have been created for the transform type
Chris@42: corresponding to the execute function, e.g. it must be a complex-DFT
Chris@42: plan for fftw_execute_dft
. Any of the planner routines for that
Chris@42: transform type, from the basic to the guru interface, could have been
Chris@42: used to create the plan, however.
Chris@42:
Chris@42: Next: Wisdom, Previous: Guru Interface, Up: FFTW Reference [Contents][Index]
Chris@42: