Next: , Previous: Guru Interface, Up: FFTW Reference


4.6 New-array Execute Functions

Normally, one executes a plan for the arrays with which the plan was created, by calling fftw_execute(plan) as described in Using Plans. However, it is possible for sophisticated users to apply a given plan to a different array using the “new-array execute” functions detailed below, provided that the following conditions are met:

The alignment issue is especially critical, because if you don't use fftw_malloc then you may have little control over the alignment of arrays in memory. For example, neither the C++ new function nor the Fortran allocate statement provide strong enough guarantees about data alignment. If you don't use fftw_malloc, therefore, you probably have to use FFTW_UNALIGNED (which disables most SIMD support). If possible, it is probably better for you to simply create multiple plans (creating a new plan is quick once one exists for a given size), or better yet re-use the same array for your transforms.

If you are tempted to use the new-array execute interface because you want to transform a known bunch of arrays of the same size, you should probably go use the advanced interface instead (see Advanced Interface)).

The new-array execute functions are:

     void fftw_execute_dft(
          const fftw_plan p,
          fftw_complex *in, fftw_complex *out);
     
     void fftw_execute_split_dft(
          const fftw_plan p,
          double *ri, double *ii, double *ro, double *io);
     
     void fftw_execute_dft_r2c(
          const fftw_plan p,
          double *in, fftw_complex *out);
     
     void fftw_execute_split_dft_r2c(
          const fftw_plan p,
          double *in, double *ro, double *io);
     
     void fftw_execute_dft_c2r(
          const fftw_plan p,
          fftw_complex *in, double *out);
     
     void fftw_execute_split_dft_c2r(
          const fftw_plan p,
          double *ri, double *ii, double *out);
     
     void fftw_execute_r2r(
          const fftw_plan p,
          double *in, double *out);

These execute the plan to compute the corresponding transform on the input/output arrays specified by the subsequent arguments. The input/output array arguments have the same meanings as the ones passed to the guru planner routines in the preceding sections. The plan is not modified, and these routines can be called as many times as desired, or intermixed with calls to the ordinary fftw_execute.

The plan must have been created for the transform type corresponding to the execute function, e.g. it must be a complex-DFT plan for fftw_execute_dft. Any of the planner routines for that transform type, from the basic to the guru interface, could have been used to create the plan, however.