cannam@95: cannam@95: cannam@95: New-array Execute Functions - FFTW 3.3.3 cannam@95: cannam@95: cannam@95: cannam@95: cannam@95: cannam@95: cannam@95: cannam@95: cannam@95: cannam@95: cannam@95: cannam@95: cannam@95: cannam@95:
cannam@95: cannam@95: cannam@95:

cannam@95: Next: , cannam@95: Previous: Guru Interface, cannam@95: Up: FFTW Reference cannam@95:


cannam@95:
cannam@95: cannam@95:

4.6 New-array Execute Functions

cannam@95: cannam@95:

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

cannam@95: cannam@95:

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

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

The new-array execute functions are: cannam@95: cannam@95:

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

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

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