Next: , Previous: , Up: FFTW MPI Reference   [Contents][Index]


6.12.3 Using MPI Plans

Once an MPI plan is created, you can execute and destroy it using fftw_execute, fftw_destroy_plan, and the other functions in the serial interface that operate on generic plans (see Using Plans).

The fftw_execute and fftw_destroy_plan functions, applied to MPI plans, are collective calls: they must be called for all processes in the communicator that was used to create the plan.

You must not use the serial new-array plan-execution functions fftw_execute_dft and so on (see New-array Execute Functions) with MPI plans. Such functions are specialized to the problem type, and there are specific new-array execute functions for MPI plans:

void fftw_mpi_execute_dft(fftw_plan p, fftw_complex *in, fftw_complex *out);
void fftw_mpi_execute_dft_r2c(fftw_plan p, double *in, fftw_complex *out);
void fftw_mpi_execute_dft_c2r(fftw_plan p, fftw_complex *in, double *out);
void fftw_mpi_execute_r2r(fftw_plan p, double *in, double *out);

These functions have the same restrictions as those of the serial new-array execute functions. They are always safe to apply to the same in and out arrays that were used to create the plan. They can only be applied to new arrarys if those arrays have the same types, dimensions, in-placeness, and alignment as the original arrays, where the best way to ensure the same alignment is to use FFTW’s fftw_malloc and related allocation functions for all arrays (see Memory Allocation). Note that distributed transposes (see FFTW MPI Transposes) use fftw_mpi_execute_r2r, since they count as rank-zero r2r plans from FFTW’s perspective.


Next: MPI Data Distribution Functions, Previous: MPI Initialization, Up: FFTW MPI Reference   [Contents][Index]