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