cannam@167: cannam@167: cannam@167: cannam@167: cannam@167:
cannam@167:cannam@167: Next: FFTW MPI Transposes, Previous: Multi-dimensional MPI DFTs of Real Data, Up: Distributed-memory FFTW with MPI [Contents][Index]
cannam@167:FFTW’s MPI interface also supports multi-dimensional ‘r2r’ cannam@167: transforms of all kinds supported by the serial interface cannam@167: (e.g. discrete cosine and sine transforms, discrete Hartley cannam@167: transforms, etc.). Only multi-dimensional ‘r2r’ transforms, not cannam@167: one-dimensional transforms, are currently parallelized. cannam@167:
cannam@167: cannam@167:These are used much like the multidimensional complex DFTs discussed
cannam@167: above, except that the data is real rather than complex, and one needs
cannam@167: to pass an r2r transform kind (fftw_r2r_kind
) for each
cannam@167: dimension as in the serial FFTW (see More DFTs of Real Data).
cannam@167:
For example, one might perform a two-dimensional L × M cannam@167: that is cannam@167: an REDFT10 (DCT-II) in the first dimension and an RODFT10 (DST-II) in cannam@167: the second dimension with code like: cannam@167:
cannam@167:const ptrdiff_t L = ..., M = ...; cannam@167: fftw_plan plan; cannam@167: double *data; cannam@167: ptrdiff_t alloc_local, local_n0, local_0_start, i, j; cannam@167: cannam@167: /* get local data size and allocate */ cannam@167: alloc_local = fftw_mpi_local_size_2d(L, M, MPI_COMM_WORLD, cannam@167: &local_n0, &local_0_start); cannam@167: data = fftw_alloc_real(alloc_local); cannam@167: cannam@167: /* create plan for in-place REDFT10 x RODFT10 */ cannam@167: plan = fftw_mpi_plan_r2r_2d(L, M, data, data, MPI_COMM_WORLD, cannam@167: FFTW_REDFT10, FFTW_RODFT10, FFTW_MEASURE); cannam@167: cannam@167: /* initialize data to some function my_function(x,y) */ cannam@167: for (i = 0; i < local_n0; ++i) for (j = 0; j < M; ++j) cannam@167: data[i*M + j] = my_function(local_0_start + i, j); cannam@167: cannam@167: /* compute transforms, in-place, as many times as desired */ cannam@167: fftw_execute(plan); cannam@167: cannam@167: fftw_destroy_plan(plan); cannam@167:
Notice that we use the same ‘local_size’ functions as we did for
cannam@167: complex data, only now we interpret the sizes in terms of real rather
cannam@167: than complex values, and correspondingly use fftw_alloc_real
.
cannam@167: