Chris@10: Chris@10: Chris@10: Other Multi-dimensional Real-data MPI Transforms - FFTW 3.3.3 Chris@10: Chris@10: Chris@10: Chris@10: Chris@10: Chris@10: Chris@10: Chris@10: Chris@10: Chris@10: Chris@10: Chris@10: Chris@10: Chris@10:
Chris@10: Chris@10: Chris@10:

Chris@10: Next: , Chris@10: Previous: Multi-dimensional MPI DFTs of Real Data, Chris@10: Up: Distributed-memory FFTW with MPI Chris@10:


Chris@10:
Chris@10: Chris@10:

6.6 Other multi-dimensional Real-Data MPI Transforms

Chris@10: Chris@10:

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

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

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

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

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