cannam@127: cannam@127: cannam@127: cannam@127: cannam@127: cannam@127: FFTW 3.3.5: Other Multi-dimensional Real-data MPI Transforms 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.6 Other multi-dimensional Real-Data MPI Transforms

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

FFTW’s MPI interface also supports multi-dimensional ‘r2r’ cannam@127: transforms of all kinds supported by the serial interface cannam@127: (e.g. discrete cosine and sine transforms, discrete Hartley cannam@127: transforms, etc.). Only multi-dimensional ‘r2r’ transforms, not cannam@127: one-dimensional transforms, are currently parallelized. cannam@127:

cannam@127: cannam@127:

These are used much like the multidimensional complex DFTs discussed cannam@127: above, except that the data is real rather than complex, and one needs cannam@127: to pass an r2r transform kind (fftw_r2r_kind) for each cannam@127: dimension as in the serial FFTW (see More DFTs of Real Data). cannam@127:

cannam@127:

For example, one might perform a two-dimensional L × M that is cannam@127: an REDFT10 (DCT-II) in the first dimension and an RODFT10 (DST-II) in cannam@127: the second dimension with code like: cannam@127:

cannam@127:
cannam@127:
    const ptrdiff_t L = ..., M = ...;
cannam@127:     fftw_plan plan;
cannam@127:     double *data;
cannam@127:     ptrdiff_t alloc_local, local_n0, local_0_start, i, j;
cannam@127: 
cannam@127:     /* get local data size and allocate */
cannam@127:     alloc_local = fftw_mpi_local_size_2d(L, M, MPI_COMM_WORLD,
cannam@127:                                          &local_n0, &local_0_start);
cannam@127:     data = fftw_alloc_real(alloc_local);
cannam@127: 
cannam@127:     /* create plan for in-place REDFT10 x RODFT10 */
cannam@127:     plan = fftw_mpi_plan_r2r_2d(L, M, data, data, MPI_COMM_WORLD,
cannam@127:                                 FFTW_REDFT10, FFTW_RODFT10, FFTW_MEASURE);
cannam@127: 
cannam@127:     /* initialize data to some function my_function(x,y) */
cannam@127:     for (i = 0; i < local_n0; ++i) for (j = 0; j < M; ++j)
cannam@127:        data[i*M + j] = my_function(local_0_start + i, j);
cannam@127: 
cannam@127:     /* compute transforms, in-place, as many times as desired */
cannam@127:     fftw_execute(plan);
cannam@127: 
cannam@127:     fftw_destroy_plan(plan);
cannam@127: 
cannam@127: cannam@127: cannam@127:

Notice that we use the same ‘local_size’ functions as we did for cannam@127: complex data, only now we interpret the sizes in terms of real rather cannam@127: than complex values, and correspondingly use fftw_alloc_real. cannam@127:

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