Chris@82: Chris@82: Chris@82: Chris@82: Chris@82:
Chris@82:Chris@82: Next: MPI Data Distribution, Previous: Linking and Initializing MPI FFTW, Up: Distributed-memory FFTW with MPI [Contents][Index]
Chris@82:Before we document the FFTW MPI interface in detail, we begin with a
Chris@82: simple example outlining how one would perform a two-dimensional
Chris@82: N0
by N1
complex DFT.
Chris@82:
#include <fftw3-mpi.h> Chris@82: Chris@82: int main(int argc, char **argv) Chris@82: { Chris@82: const ptrdiff_t N0 = ..., N1 = ...; Chris@82: fftw_plan plan; Chris@82: fftw_complex *data; Chris@82: ptrdiff_t alloc_local, local_n0, local_0_start, i, j; Chris@82: Chris@82: MPI_Init(&argc, &argv); Chris@82: fftw_mpi_init(); Chris@82: Chris@82: /* get local data size and allocate */ Chris@82: alloc_local = fftw_mpi_local_size_2d(N0, N1, MPI_COMM_WORLD, Chris@82: &local_n0, &local_0_start); Chris@82: data = fftw_alloc_complex(alloc_local); Chris@82: Chris@82: /* create plan for in-place forward DFT */ Chris@82: plan = fftw_mpi_plan_dft_2d(N0, N1, data, data, MPI_COMM_WORLD, Chris@82: FFTW_FORWARD, FFTW_ESTIMATE); 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 < N1; ++j) Chris@82: data[i*N1 + 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: MPI_Finalize(); Chris@82: } Chris@82:
As can be seen above, the MPI interface follows the same basic style Chris@82: of allocate/plan/execute/destroy as the serial FFTW routines. All of Chris@82: the MPI-specific routines are prefixed with ‘fftw_mpi_’ instead Chris@82: of ‘fftw_’. There are a few important differences, however: Chris@82:
Chris@82:First, we must call fftw_mpi_init()
after calling
Chris@82: MPI_Init
(required in all MPI programs) and before calling any
Chris@82: other ‘fftw_mpi_’ routine.
Chris@82:
Chris@82:
Chris@82:
Second, when we create the plan with fftw_mpi_plan_dft_2d
,
Chris@82: analogous to fftw_plan_dft_2d
, we pass an additional argument:
Chris@82: the communicator, indicating which processes will participate in the
Chris@82: transform (here MPI_COMM_WORLD
, indicating all processes).
Chris@82: Whenever you create, execute, or destroy a plan for an MPI transform,
Chris@82: you must call the corresponding FFTW routine on all processes
Chris@82: in the communicator for that transform. (That is, these are
Chris@82: collective calls.) Note that the plan for the MPI transform
Chris@82: uses the standard fftw_execute
and fftw_destroy
routines
Chris@82: (on the other hand, there are MPI-specific new-array execute functions
Chris@82: documented below).
Chris@82:
Chris@82:
Chris@82:
Chris@82:
Third, all of the FFTW MPI routines take ptrdiff_t
arguments
Chris@82: instead of int
as for the serial FFTW. ptrdiff_t
is a
Chris@82: standard C integer type which is (at least) 32 bits wide on a 32-bit
Chris@82: machine and 64 bits wide on a 64-bit machine. This is to make it easy
Chris@82: to specify very large parallel transforms on a 64-bit machine. (You
Chris@82: can specify 64-bit transform sizes in the serial FFTW, too, but only
Chris@82: by using the ‘guru64’ planner interface. See 64-bit Guru Interface.)
Chris@82:
Chris@82:
Chris@82:
Fourth, and most importantly, you don’t allocate the entire
Chris@82: two-dimensional array on each process. Instead, you call
Chris@82: fftw_mpi_local_size_2d
to find out what portion of the
Chris@82: array resides on each processor, and how much space to allocate.
Chris@82: Here, the portion of the array on each process is a local_n0
by
Chris@82: N1
slice of the total array, starting at index
Chris@82: local_0_start
. The total number of fftw_complex
numbers
Chris@82: to allocate is given by the alloc_local
return value, which
Chris@82: may be greater than local_n0 * N1
(in case some
Chris@82: intermediate calculations require additional storage). The data
Chris@82: distribution in FFTW’s MPI interface is described in more detail by
Chris@82: the next section.
Chris@82:
Chris@82:
Chris@82:
Given the portion of the array that resides on the local process, it
Chris@82: is straightforward to initialize the data (here to a function
Chris@82: myfunction
) and otherwise manipulate it. Of course, at the end
Chris@82: of the program you may want to output the data somehow, but
Chris@82: synchronizing this output is up to you and is beyond the scope of this
Chris@82: manual. (One good way to output a large multi-dimensional distributed
Chris@82: array in MPI to a portable binary file is to use the free HDF5
Chris@82: library; see the HDF home page.)
Chris@82:
Chris@82:
Chris@82:
Chris@82: Next: MPI Data Distribution, Previous: Linking and Initializing MPI FFTW, Up: Distributed-memory FFTW with MPI [Contents][Index]
Chris@82: