cannam@127: cannam@127: cannam@127: cannam@127: cannam@127: cannam@127: FFTW 3.3.5: Multi-dimensional MPI DFTs of Real Data cannam@127: cannam@127: cannam@127: cannam@127: cannam@127: cannam@127: cannam@127: cannam@127: cannam@127: cannam@127: cannam@127: cannam@127: cannam@127: cannam@127: cannam@127: cannam@127: cannam@127: cannam@127: cannam@127: cannam@127: cannam@127:
cannam@127:

cannam@127: Next: , Previous: , Up: Distributed-memory FFTW with MPI   [Contents][Index]

cannam@127:
cannam@127:
cannam@127: cannam@127:

6.5 Multi-dimensional MPI DFTs of Real Data

cannam@127: cannam@127:

FFTW’s MPI interface also supports multi-dimensional DFTs of real cannam@127: data, similar to the serial r2c and c2r interfaces. (Parallel cannam@127: one-dimensional real-data DFTs are not currently supported; you must cannam@127: use a complex transform and set the imaginary parts of the inputs to cannam@127: zero.) cannam@127:

cannam@127:

The key points to understand for r2c and c2r MPI transforms (compared cannam@127: to the MPI complex DFTs or the serial r2c/c2r transforms), are: cannam@127:

cannam@127: cannam@127: cannam@127:

For example suppose we are performing an out-of-place r2c transform of cannam@127: L × M × N real data [padded to L × M × 2(N/2+1)], cannam@127: resulting in L × M × N/2+1 complex data. Similar to the cannam@127: example in 2d MPI example, we might do something like: cannam@127:

cannam@127:
cannam@127:
#include <fftw3-mpi.h>
cannam@127: 
cannam@127: int main(int argc, char **argv)
cannam@127: {
cannam@127:     const ptrdiff_t L = ..., M = ..., N = ...;
cannam@127:     fftw_plan plan;
cannam@127:     double *rin;
cannam@127:     fftw_complex *cout;
cannam@127:     ptrdiff_t alloc_local, local_n0, local_0_start, i, j, k;
cannam@127: 
cannam@127:     MPI_Init(&argc, &argv);
cannam@127:     fftw_mpi_init();
cannam@127: 
cannam@127:     /* get local data size and allocate */
cannam@127:     alloc_local = fftw_mpi_local_size_3d(L, M, N/2+1, MPI_COMM_WORLD,
cannam@127:                                          &local_n0, &local_0_start);
cannam@127:     rin = fftw_alloc_real(2 * alloc_local);
cannam@127:     cout = fftw_alloc_complex(alloc_local);
cannam@127: 
cannam@127:     /* create plan for out-of-place r2c DFT */
cannam@127:     plan = fftw_mpi_plan_dft_r2c_3d(L, M, N, rin, cout, MPI_COMM_WORLD,
cannam@127:                                     FFTW_MEASURE);
cannam@127: 
cannam@127:     /* initialize rin to some function my_func(x,y,z) */
cannam@127:     for (i = 0; i < local_n0; ++i)
cannam@127:        for (j = 0; j < M; ++j)
cannam@127:          for (k = 0; k < N; ++k)
cannam@127:        rin[(i*M + j) * (2*(N/2+1)) + k] = my_func(local_0_start+i, j, k);
cannam@127: 
cannam@127:     /* compute transforms as many times as desired */
cannam@127:     fftw_execute(plan);
cannam@127: 
cannam@127:     fftw_destroy_plan(plan);
cannam@127: 
cannam@127:     MPI_Finalize();
cannam@127: }
cannam@127: 
cannam@127: cannam@127: cannam@127: cannam@127:

Note that we allocated rin using fftw_alloc_real with an cannam@127: argument of 2 * alloc_local: since alloc_local is the cannam@127: number of complex values to allocate, the number of real cannam@127: values is twice as many. The rin array is then cannam@127: local_n0 × M × 2(N/2+1) in row-major order, so its cannam@127: (i,j,k) element is at the index (i*M + j) * (2*(N/2+1)) + cannam@127: k (see Multi-dimensional Array Format). cannam@127:

cannam@127: cannam@127: cannam@127: cannam@127:

As for the complex transforms, improved performance can be obtained by cannam@127: specifying that the output is the transpose of the input or vice versa cannam@127: (see Transposed distributions). In our L × M × N r2c cannam@127: example, including FFTW_TRANSPOSED_OUT in the flags means that cannam@127: the input would be a padded L × M × 2(N/2+1) real array cannam@127: distributed over the L dimension, while the output would be a cannam@127: M × L × N/2+1 complex array distributed over the M cannam@127: dimension. To perform the inverse c2r transform with the same data cannam@127: distributions, you would use the FFTW_TRANSPOSED_IN flag. cannam@127:

cannam@127:
cannam@127:
cannam@127:

cannam@127: Next: , Previous: , Up: Distributed-memory FFTW with MPI   [Contents][Index]

cannam@127:
cannam@127: cannam@127: cannam@127: cannam@127: cannam@127: