Chris@82: Chris@82: Chris@82: Chris@82: Chris@82: Chris@82: FFTW 3.3.8: Other Multi-dimensional Real-data MPI Transforms Chris@82: Chris@82: Chris@82: Chris@82: Chris@82: Chris@82: Chris@82: Chris@82: Chris@82: Chris@82: Chris@82: Chris@82: Chris@82: Chris@82: Chris@82: Chris@82: Chris@82: Chris@82: Chris@82: Chris@82: Chris@82:
Chris@82:

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

Chris@82:
Chris@82:
Chris@82: Chris@82:

6.6 Other multi-dimensional Real-Data MPI Transforms

Chris@82: Chris@82: Chris@82:

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

Chris@82: Chris@82:

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

Chris@82:

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

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

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

Chris@82: Chris@82: Chris@82: Chris@82: Chris@82: