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