annotate src/fftw-3.3.3/doc/mpi.texi @ 10:37bf6b4a2645

Add FFTW3
author Chris Cannam
date Wed, 20 Mar 2013 15:35:50 +0000
parents
children
rev   line source
Chris@10 1 @node Distributed-memory FFTW with MPI, Calling FFTW from Modern Fortran, Multi-threaded FFTW, Top
Chris@10 2 @chapter Distributed-memory FFTW with MPI
Chris@10 3 @cindex MPI
Chris@10 4
Chris@10 5 @cindex parallel transform
Chris@10 6 In this chapter we document the parallel FFTW routines for parallel
Chris@10 7 systems supporting the MPI message-passing interface. Unlike the
Chris@10 8 shared-memory threads described in the previous chapter, MPI allows
Chris@10 9 you to use @emph{distributed-memory} parallelism, where each CPU has
Chris@10 10 its own separate memory, and which can scale up to clusters of many
Chris@10 11 thousands of processors. This capability comes at a price, however:
Chris@10 12 each process only stores a @emph{portion} of the data to be
Chris@10 13 transformed, which means that the data structures and
Chris@10 14 programming-interface are quite different from the serial or threads
Chris@10 15 versions of FFTW.
Chris@10 16 @cindex data distribution
Chris@10 17
Chris@10 18
Chris@10 19 Distributed-memory parallelism is especially useful when you are
Chris@10 20 transforming arrays so large that they do not fit into the memory of a
Chris@10 21 single processor. The storage per-process required by FFTW's MPI
Chris@10 22 routines is proportional to the total array size divided by the number
Chris@10 23 of processes. Conversely, distributed-memory parallelism can easily
Chris@10 24 pose an unacceptably high communications overhead for small problems;
Chris@10 25 the threshold problem size for which parallelism becomes advantageous
Chris@10 26 will depend on the precise problem you are interested in, your
Chris@10 27 hardware, and your MPI implementation.
Chris@10 28
Chris@10 29 A note on terminology: in MPI, you divide the data among a set of
Chris@10 30 ``processes'' which each run in their own memory address space.
Chris@10 31 Generally, each process runs on a different physical processor, but
Chris@10 32 this is not required. A set of processes in MPI is described by an
Chris@10 33 opaque data structure called a ``communicator,'' the most common of
Chris@10 34 which is the predefined communicator @code{MPI_COMM_WORLD} which
Chris@10 35 refers to @emph{all} processes. For more information on these and
Chris@10 36 other concepts common to all MPI programs, we refer the reader to the
Chris@10 37 documentation at @uref{http://www.mcs.anl.gov/research/projects/mpi/, the MPI home
Chris@10 38 page}.
Chris@10 39 @cindex MPI communicator
Chris@10 40 @ctindex MPI_COMM_WORLD
Chris@10 41
Chris@10 42
Chris@10 43 We assume in this chapter that the reader is familiar with the usage
Chris@10 44 of the serial (uniprocessor) FFTW, and focus only on the concepts new
Chris@10 45 to the MPI interface.
Chris@10 46
Chris@10 47 @menu
Chris@10 48 * FFTW MPI Installation::
Chris@10 49 * Linking and Initializing MPI FFTW::
Chris@10 50 * 2d MPI example::
Chris@10 51 * MPI Data Distribution::
Chris@10 52 * Multi-dimensional MPI DFTs of Real Data::
Chris@10 53 * Other Multi-dimensional Real-data MPI Transforms::
Chris@10 54 * FFTW MPI Transposes::
Chris@10 55 * FFTW MPI Wisdom::
Chris@10 56 * Avoiding MPI Deadlocks::
Chris@10 57 * FFTW MPI Performance Tips::
Chris@10 58 * Combining MPI and Threads::
Chris@10 59 * FFTW MPI Reference::
Chris@10 60 * FFTW MPI Fortran Interface::
Chris@10 61 @end menu
Chris@10 62
Chris@10 63 @c ------------------------------------------------------------
Chris@10 64 @node FFTW MPI Installation, Linking and Initializing MPI FFTW, Distributed-memory FFTW with MPI, Distributed-memory FFTW with MPI
Chris@10 65 @section FFTW MPI Installation
Chris@10 66
Chris@10 67 All of the FFTW MPI code is located in the @code{mpi} subdirectory of
Chris@10 68 the FFTW package. On Unix systems, the FFTW MPI libraries and header
Chris@10 69 files are automatically configured, compiled, and installed along with
Chris@10 70 the uniprocessor FFTW libraries simply by including
Chris@10 71 @code{--enable-mpi} in the flags to the @code{configure} script
Chris@10 72 (@pxref{Installation on Unix}).
Chris@10 73 @fpindex configure
Chris@10 74
Chris@10 75
Chris@10 76 Any implementation of the MPI standard, version 1 or later, should
Chris@10 77 work with FFTW. The @code{configure} script will attempt to
Chris@10 78 automatically detect how to compile and link code using your MPI
Chris@10 79 implementation. In some cases, especially if you have multiple
Chris@10 80 different MPI implementations installed or have an unusual MPI
Chris@10 81 software package, you may need to provide this information explicitly.
Chris@10 82
Chris@10 83 Most commonly, one compiles MPI code by invoking a special compiler
Chris@10 84 command, typically @code{mpicc} for C code. The @code{configure}
Chris@10 85 script knows the most common names for this command, but you can
Chris@10 86 specify the MPI compilation command explicitly by setting the
Chris@10 87 @code{MPICC} variable, as in @samp{./configure MPICC=mpicc ...}.
Chris@10 88 @fpindex mpicc
Chris@10 89
Chris@10 90
Chris@10 91 If, instead of a special compiler command, you need to link a certain
Chris@10 92 library, you can specify the link command via the @code{MPILIBS}
Chris@10 93 variable, as in @samp{./configure MPILIBS=-lmpi ...}. Note that if
Chris@10 94 your MPI library is installed in a non-standard location (one the
Chris@10 95 compiler does not know about by default), you may also have to specify
Chris@10 96 the location of the library and header files via @code{LDFLAGS} and
Chris@10 97 @code{CPPFLAGS} variables, respectively, as in @samp{./configure
Chris@10 98 LDFLAGS=-L/path/to/mpi/libs CPPFLAGS=-I/path/to/mpi/include ...}.
Chris@10 99
Chris@10 100 @c ------------------------------------------------------------
Chris@10 101 @node Linking and Initializing MPI FFTW, 2d MPI example, FFTW MPI Installation, Distributed-memory FFTW with MPI
Chris@10 102 @section Linking and Initializing MPI FFTW
Chris@10 103
Chris@10 104 Programs using the MPI FFTW routines should be linked with
Chris@10 105 @code{-lfftw3_mpi -lfftw3 -lm} on Unix in double precision,
Chris@10 106 @code{-lfftw3f_mpi -lfftw3f -lm} in single precision, and so on
Chris@10 107 (@pxref{Precision}). You will also need to link with whatever library
Chris@10 108 is responsible for MPI on your system; in most MPI implementations,
Chris@10 109 there is a special compiler alias named @code{mpicc} to compile and
Chris@10 110 link MPI code.
Chris@10 111 @fpindex mpicc
Chris@10 112 @cindex linking on Unix
Chris@10 113 @cindex precision
Chris@10 114
Chris@10 115
Chris@10 116 @findex fftw_init_threads
Chris@10 117 Before calling any FFTW routines except possibly
Chris@10 118 @code{fftw_init_threads} (@pxref{Combining MPI and Threads}), but after calling
Chris@10 119 @code{MPI_Init}, you should call the function:
Chris@10 120
Chris@10 121 @example
Chris@10 122 void fftw_mpi_init(void);
Chris@10 123 @end example
Chris@10 124 @findex fftw_mpi_init
Chris@10 125
Chris@10 126 If, at the end of your program, you want to get rid of all memory and
Chris@10 127 other resources allocated internally by FFTW, for both the serial and
Chris@10 128 MPI routines, you can call:
Chris@10 129
Chris@10 130 @example
Chris@10 131 void fftw_mpi_cleanup(void);
Chris@10 132 @end example
Chris@10 133 @findex fftw_mpi_cleanup
Chris@10 134
Chris@10 135 which is much like the @code{fftw_cleanup()} function except that it
Chris@10 136 also gets rid of FFTW's MPI-related data. You must @emph{not} execute
Chris@10 137 any previously created plans after calling this function.
Chris@10 138
Chris@10 139 @c ------------------------------------------------------------
Chris@10 140 @node 2d MPI example, MPI Data Distribution, Linking and Initializing MPI FFTW, Distributed-memory FFTW with MPI
Chris@10 141 @section 2d MPI example
Chris@10 142
Chris@10 143 Before we document the FFTW MPI interface in detail, we begin with a
Chris@10 144 simple example outlining how one would perform a two-dimensional
Chris@10 145 @code{N0} by @code{N1} complex DFT.
Chris@10 146
Chris@10 147 @example
Chris@10 148 #include <fftw3-mpi.h>
Chris@10 149
Chris@10 150 int main(int argc, char **argv)
Chris@10 151 @{
Chris@10 152 const ptrdiff_t N0 = ..., N1 = ...;
Chris@10 153 fftw_plan plan;
Chris@10 154 fftw_complex *data;
Chris@10 155 ptrdiff_t alloc_local, local_n0, local_0_start, i, j;
Chris@10 156
Chris@10 157 MPI_Init(&argc, &argv);
Chris@10 158 fftw_mpi_init();
Chris@10 159
Chris@10 160 /* @r{get local data size and allocate} */
Chris@10 161 alloc_local = fftw_mpi_local_size_2d(N0, N1, MPI_COMM_WORLD,
Chris@10 162 &local_n0, &local_0_start);
Chris@10 163 data = fftw_alloc_complex(alloc_local);
Chris@10 164
Chris@10 165 /* @r{create plan for in-place forward DFT} */
Chris@10 166 plan = fftw_mpi_plan_dft_2d(N0, N1, data, data, MPI_COMM_WORLD,
Chris@10 167 FFTW_FORWARD, FFTW_ESTIMATE);
Chris@10 168
Chris@10 169 /* @r{initialize data to some function} my_function(x,y) */
Chris@10 170 for (i = 0; i < local_n0; ++i) for (j = 0; j < N1; ++j)
Chris@10 171 data[i*N1 + j] = my_function(local_0_start + i, j);
Chris@10 172
Chris@10 173 /* @r{compute transforms, in-place, as many times as desired} */
Chris@10 174 fftw_execute(plan);
Chris@10 175
Chris@10 176 fftw_destroy_plan(plan);
Chris@10 177
Chris@10 178 MPI_Finalize();
Chris@10 179 @}
Chris@10 180 @end example
Chris@10 181
Chris@10 182 As can be seen above, the MPI interface follows the same basic style
Chris@10 183 of allocate/plan/execute/destroy as the serial FFTW routines. All of
Chris@10 184 the MPI-specific routines are prefixed with @samp{fftw_mpi_} instead
Chris@10 185 of @samp{fftw_}. There are a few important differences, however:
Chris@10 186
Chris@10 187 First, we must call @code{fftw_mpi_init()} after calling
Chris@10 188 @code{MPI_Init} (required in all MPI programs) and before calling any
Chris@10 189 other @samp{fftw_mpi_} routine.
Chris@10 190 @findex MPI_Init
Chris@10 191 @findex fftw_mpi_init
Chris@10 192
Chris@10 193
Chris@10 194 Second, when we create the plan with @code{fftw_mpi_plan_dft_2d},
Chris@10 195 analogous to @code{fftw_plan_dft_2d}, we pass an additional argument:
Chris@10 196 the communicator, indicating which processes will participate in the
Chris@10 197 transform (here @code{MPI_COMM_WORLD}, indicating all processes).
Chris@10 198 Whenever you create, execute, or destroy a plan for an MPI transform,
Chris@10 199 you must call the corresponding FFTW routine on @emph{all} processes
Chris@10 200 in the communicator for that transform. (That is, these are
Chris@10 201 @emph{collective} calls.) Note that the plan for the MPI transform
Chris@10 202 uses the standard @code{fftw_execute} and @code{fftw_destroy} routines
Chris@10 203 (on the other hand, there are MPI-specific new-array execute functions
Chris@10 204 documented below).
Chris@10 205 @cindex collective function
Chris@10 206 @findex fftw_mpi_plan_dft_2d
Chris@10 207 @ctindex MPI_COMM_WORLD
Chris@10 208
Chris@10 209
Chris@10 210 Third, all of the FFTW MPI routines take @code{ptrdiff_t} arguments
Chris@10 211 instead of @code{int} as for the serial FFTW. @code{ptrdiff_t} is a
Chris@10 212 standard C integer type which is (at least) 32 bits wide on a 32-bit
Chris@10 213 machine and 64 bits wide on a 64-bit machine. This is to make it easy
Chris@10 214 to specify very large parallel transforms on a 64-bit machine. (You
Chris@10 215 can specify 64-bit transform sizes in the serial FFTW, too, but only
Chris@10 216 by using the @samp{guru64} planner interface. @xref{64-bit Guru
Chris@10 217 Interface}.)
Chris@10 218 @tindex ptrdiff_t
Chris@10 219 @cindex 64-bit architecture
Chris@10 220
Chris@10 221
Chris@10 222 Fourth, and most importantly, you don't allocate the entire
Chris@10 223 two-dimensional array on each process. Instead, you call
Chris@10 224 @code{fftw_mpi_local_size_2d} to find out what @emph{portion} of the
Chris@10 225 array resides on each processor, and how much space to allocate.
Chris@10 226 Here, the portion of the array on each process is a @code{local_n0} by
Chris@10 227 @code{N1} slice of the total array, starting at index
Chris@10 228 @code{local_0_start}. The total number of @code{fftw_complex} numbers
Chris@10 229 to allocate is given by the @code{alloc_local} return value, which
Chris@10 230 @emph{may} be greater than @code{local_n0 * N1} (in case some
Chris@10 231 intermediate calculations require additional storage). The data
Chris@10 232 distribution in FFTW's MPI interface is described in more detail by
Chris@10 233 the next section.
Chris@10 234 @findex fftw_mpi_local_size_2d
Chris@10 235 @cindex data distribution
Chris@10 236
Chris@10 237
Chris@10 238 Given the portion of the array that resides on the local process, it
Chris@10 239 is straightforward to initialize the data (here to a function
Chris@10 240 @code{myfunction}) and otherwise manipulate it. Of course, at the end
Chris@10 241 of the program you may want to output the data somehow, but
Chris@10 242 synchronizing this output is up to you and is beyond the scope of this
Chris@10 243 manual. (One good way to output a large multi-dimensional distributed
Chris@10 244 array in MPI to a portable binary file is to use the free HDF5
Chris@10 245 library; see the @uref{http://www.hdfgroup.org/, HDF home page}.)
Chris@10 246 @cindex HDF5
Chris@10 247 @cindex MPI I/O
Chris@10 248
Chris@10 249 @c ------------------------------------------------------------
Chris@10 250 @node MPI Data Distribution, Multi-dimensional MPI DFTs of Real Data, 2d MPI example, Distributed-memory FFTW with MPI
Chris@10 251 @section MPI Data Distribution
Chris@10 252 @cindex data distribution
Chris@10 253
Chris@10 254 The most important concept to understand in using FFTW's MPI interface
Chris@10 255 is the data distribution. With a serial or multithreaded FFT, all of
Chris@10 256 the inputs and outputs are stored as a single contiguous chunk of
Chris@10 257 memory. With a distributed-memory FFT, the inputs and outputs are
Chris@10 258 broken into disjoint blocks, one per process.
Chris@10 259
Chris@10 260 In particular, FFTW uses a @emph{1d block distribution} of the data,
Chris@10 261 distributed along the @emph{first dimension}. For example, if you
Chris@10 262 want to perform a @twodims{100,200} complex DFT, distributed over 4
Chris@10 263 processes, each process will get a @twodims{25,200} slice of the data.
Chris@10 264 That is, process 0 will get rows 0 through 24, process 1 will get rows
Chris@10 265 25 through 49, process 2 will get rows 50 through 74, and process 3
Chris@10 266 will get rows 75 through 99. If you take the same array but
Chris@10 267 distribute it over 3 processes, then it is not evenly divisible so the
Chris@10 268 different processes will have unequal chunks. FFTW's default choice
Chris@10 269 in this case is to assign 34 rows to processes 0 and 1, and 32 rows to
Chris@10 270 process 2.
Chris@10 271 @cindex block distribution
Chris@10 272
Chris@10 273
Chris@10 274 FFTW provides several @samp{fftw_mpi_local_size} routines that you can
Chris@10 275 call to find out what portion of an array is stored on the current
Chris@10 276 process. In most cases, you should use the default block sizes picked
Chris@10 277 by FFTW, but it is also possible to specify your own block size. For
Chris@10 278 example, with a @twodims{100,200} array on three processes, you can
Chris@10 279 tell FFTW to use a block size of 40, which would assign 40 rows to
Chris@10 280 processes 0 and 1, and 20 rows to process 2. FFTW's default is to
Chris@10 281 divide the data equally among the processes if possible, and as best
Chris@10 282 it can otherwise. The rows are always assigned in ``rank order,''
Chris@10 283 i.e. process 0 gets the first block of rows, then process 1, and so
Chris@10 284 on. (You can change this by using @code{MPI_Comm_split} to create a
Chris@10 285 new communicator with re-ordered processes.) However, you should
Chris@10 286 always call the @samp{fftw_mpi_local_size} routines, if possible,
Chris@10 287 rather than trying to predict FFTW's distribution choices.
Chris@10 288
Chris@10 289 In particular, it is critical that you allocate the storage size that
Chris@10 290 is returned by @samp{fftw_mpi_local_size}, which is @emph{not}
Chris@10 291 necessarily the size of the local slice of the array. The reason is
Chris@10 292 that intermediate steps of FFTW's algorithms involve transposing the
Chris@10 293 array and redistributing the data, so at these intermediate steps FFTW
Chris@10 294 may require more local storage space (albeit always proportional to
Chris@10 295 the total size divided by the number of processes). The
Chris@10 296 @samp{fftw_mpi_local_size} functions know how much storage is required
Chris@10 297 for these intermediate steps and tell you the correct amount to
Chris@10 298 allocate.
Chris@10 299
Chris@10 300 @menu
Chris@10 301 * Basic and advanced distribution interfaces::
Chris@10 302 * Load balancing::
Chris@10 303 * Transposed distributions::
Chris@10 304 * One-dimensional distributions::
Chris@10 305 @end menu
Chris@10 306
Chris@10 307 @node Basic and advanced distribution interfaces, Load balancing, MPI Data Distribution, MPI Data Distribution
Chris@10 308 @subsection Basic and advanced distribution interfaces
Chris@10 309
Chris@10 310 As with the planner interface, the @samp{fftw_mpi_local_size}
Chris@10 311 distribution interface is broken into basic and advanced
Chris@10 312 (@samp{_many}) interfaces, where the latter allows you to specify the
Chris@10 313 block size manually and also to request block sizes when computing
Chris@10 314 multiple transforms simultaneously. These functions are documented
Chris@10 315 more exhaustively by the FFTW MPI Reference, but we summarize the
Chris@10 316 basic ideas here using a couple of two-dimensional examples.
Chris@10 317
Chris@10 318 For the @twodims{100,200} complex-DFT example, above, we would find
Chris@10 319 the distribution by calling the following function in the basic
Chris@10 320 interface:
Chris@10 321
Chris@10 322 @example
Chris@10 323 ptrdiff_t fftw_mpi_local_size_2d(ptrdiff_t n0, ptrdiff_t n1, MPI_Comm comm,
Chris@10 324 ptrdiff_t *local_n0, ptrdiff_t *local_0_start);
Chris@10 325 @end example
Chris@10 326 @findex fftw_mpi_local_size_2d
Chris@10 327
Chris@10 328 Given the total size of the data to be transformed (here, @code{n0 =
Chris@10 329 100} and @code{n1 = 200}) and an MPI communicator (@code{comm}), this
Chris@10 330 function provides three numbers.
Chris@10 331
Chris@10 332 First, it describes the shape of the local data: the current process
Chris@10 333 should store a @code{local_n0} by @code{n1} slice of the overall
Chris@10 334 dataset, in row-major order (@code{n1} dimension contiguous), starting
Chris@10 335 at index @code{local_0_start}. That is, if the total dataset is
Chris@10 336 viewed as a @code{n0} by @code{n1} matrix, the current process should
Chris@10 337 store the rows @code{local_0_start} to
Chris@10 338 @code{local_0_start+local_n0-1}. Obviously, if you are running with
Chris@10 339 only a single MPI process, that process will store the entire array:
Chris@10 340 @code{local_0_start} will be zero and @code{local_n0} will be
Chris@10 341 @code{n0}. @xref{Row-major Format}.
Chris@10 342 @cindex row-major
Chris@10 343
Chris@10 344
Chris@10 345 Second, the return value is the total number of data elements (e.g.,
Chris@10 346 complex numbers for a complex DFT) that should be allocated for the
Chris@10 347 input and output arrays on the current process (ideally with
Chris@10 348 @code{fftw_malloc} or an @samp{fftw_alloc} function, to ensure optimal
Chris@10 349 alignment). It might seem that this should always be equal to
Chris@10 350 @code{local_n0 * n1}, but this is @emph{not} the case. FFTW's
Chris@10 351 distributed FFT algorithms require data redistributions at
Chris@10 352 intermediate stages of the transform, and in some circumstances this
Chris@10 353 may require slightly larger local storage. This is discussed in more
Chris@10 354 detail below, under @ref{Load balancing}.
Chris@10 355 @findex fftw_malloc
Chris@10 356 @findex fftw_alloc_complex
Chris@10 357
Chris@10 358
Chris@10 359 @cindex advanced interface
Chris@10 360 The advanced-interface @samp{local_size} function for multidimensional
Chris@10 361 transforms returns the same three things (@code{local_n0},
Chris@10 362 @code{local_0_start}, and the total number of elements to allocate),
Chris@10 363 but takes more inputs:
Chris@10 364
Chris@10 365 @example
Chris@10 366 ptrdiff_t fftw_mpi_local_size_many(int rnk, const ptrdiff_t *n,
Chris@10 367 ptrdiff_t howmany,
Chris@10 368 ptrdiff_t block0,
Chris@10 369 MPI_Comm comm,
Chris@10 370 ptrdiff_t *local_n0,
Chris@10 371 ptrdiff_t *local_0_start);
Chris@10 372 @end example
Chris@10 373 @findex fftw_mpi_local_size_many
Chris@10 374
Chris@10 375 The two-dimensional case above corresponds to @code{rnk = 2} and an
Chris@10 376 array @code{n} of length 2 with @code{n[0] = n0} and @code{n[1] = n1}.
Chris@10 377 This routine is for any @code{rnk > 1}; one-dimensional transforms
Chris@10 378 have their own interface because they work slightly differently, as
Chris@10 379 discussed below.
Chris@10 380
Chris@10 381 First, the advanced interface allows you to perform multiple
Chris@10 382 transforms at once, of interleaved data, as specified by the
Chris@10 383 @code{howmany} parameter. (@code{hoamany} is 1 for a single
Chris@10 384 transform.)
Chris@10 385
Chris@10 386 Second, here you can specify your desired block size in the @code{n0}
Chris@10 387 dimension, @code{block0}. To use FFTW's default block size, pass
Chris@10 388 @code{FFTW_MPI_DEFAULT_BLOCK} (0) for @code{block0}. Otherwise, on
Chris@10 389 @code{P} processes, FFTW will return @code{local_n0} equal to
Chris@10 390 @code{block0} on the first @code{P / block0} processes (rounded down),
Chris@10 391 return @code{local_n0} equal to @code{n0 - block0 * (P / block0)} on
Chris@10 392 the next process, and @code{local_n0} equal to zero on any remaining
Chris@10 393 processes. In general, we recommend using the default block size
Chris@10 394 (which corresponds to @code{n0 / P}, rounded up).
Chris@10 395 @ctindex FFTW_MPI_DEFAULT_BLOCK
Chris@10 396 @cindex block distribution
Chris@10 397
Chris@10 398
Chris@10 399 For example, suppose you have @code{P = 4} processes and @code{n0 =
Chris@10 400 21}. The default will be a block size of @code{6}, which will give
Chris@10 401 @code{local_n0 = 6} on the first three processes and @code{local_n0 =
Chris@10 402 3} on the last process. Instead, however, you could specify
Chris@10 403 @code{block0 = 5} if you wanted, which would give @code{local_n0 = 5}
Chris@10 404 on processes 0 to 2, @code{local_n0 = 6} on process 3. (This choice,
Chris@10 405 while it may look superficially more ``balanced,'' has the same
Chris@10 406 critical path as FFTW's default but requires more communications.)
Chris@10 407
Chris@10 408 @node Load balancing, Transposed distributions, Basic and advanced distribution interfaces, MPI Data Distribution
Chris@10 409 @subsection Load balancing
Chris@10 410 @cindex load balancing
Chris@10 411
Chris@10 412 Ideally, when you parallelize a transform over some @math{P}
Chris@10 413 processes, each process should end up with work that takes equal time.
Chris@10 414 Otherwise, all of the processes end up waiting on whichever process is
Chris@10 415 slowest. This goal is known as ``load balancing.'' In this section,
Chris@10 416 we describe the circumstances under which FFTW is able to load-balance
Chris@10 417 well, and in particular how you should choose your transform size in
Chris@10 418 order to load balance.
Chris@10 419
Chris@10 420 Load balancing is especially difficult when you are parallelizing over
Chris@10 421 heterogeneous machines; for example, if one of your processors is a
Chris@10 422 old 486 and another is a Pentium IV, obviously you should give the
Chris@10 423 Pentium more work to do than the 486 since the latter is much slower.
Chris@10 424 FFTW does not deal with this problem, however---it assumes that your
Chris@10 425 processes run on hardware of comparable speed, and that the goal is
Chris@10 426 therefore to divide the problem as equally as possible.
Chris@10 427
Chris@10 428 For a multi-dimensional complex DFT, FFTW can divide the problem
Chris@10 429 equally among the processes if: (i) the @emph{first} dimension
Chris@10 430 @code{n0} is divisible by @math{P}; and (ii), the @emph{product} of
Chris@10 431 the subsequent dimensions is divisible by @math{P}. (For the advanced
Chris@10 432 interface, where you can specify multiple simultaneous transforms via
Chris@10 433 some ``vector'' length @code{howmany}, a factor of @code{howmany} is
Chris@10 434 included in the product of the subsequent dimensions.)
Chris@10 435
Chris@10 436 For a one-dimensional complex DFT, the length @code{N} of the data
Chris@10 437 should be divisible by @math{P} @emph{squared} to be able to divide
Chris@10 438 the problem equally among the processes.
Chris@10 439
Chris@10 440 @node Transposed distributions, One-dimensional distributions, Load balancing, MPI Data Distribution
Chris@10 441 @subsection Transposed distributions
Chris@10 442
Chris@10 443 Internally, FFTW's MPI transform algorithms work by first computing
Chris@10 444 transforms of the data local to each process, then by globally
Chris@10 445 @emph{transposing} the data in some fashion to redistribute the data
Chris@10 446 among the processes, transforming the new data local to each process,
Chris@10 447 and transposing back. For example, a two-dimensional @code{n0} by
Chris@10 448 @code{n1} array, distributed across the @code{n0} dimension, is
Chris@10 449 transformd by: (i) transforming the @code{n1} dimension, which are
Chris@10 450 local to each process; (ii) transposing to an @code{n1} by @code{n0}
Chris@10 451 array, distributed across the @code{n1} dimension; (iii) transforming
Chris@10 452 the @code{n0} dimension, which is now local to each process; (iv)
Chris@10 453 transposing back.
Chris@10 454 @cindex transpose
Chris@10 455
Chris@10 456
Chris@10 457 However, in many applications it is acceptable to compute a
Chris@10 458 multidimensional DFT whose results are produced in transposed order
Chris@10 459 (e.g., @code{n1} by @code{n0} in two dimensions). This provides a
Chris@10 460 significant performance advantage, because it means that the final
Chris@10 461 transposition step can be omitted. FFTW supports this optimization,
Chris@10 462 which you specify by passing the flag @code{FFTW_MPI_TRANSPOSED_OUT}
Chris@10 463 to the planner routines. To compute the inverse transform of
Chris@10 464 transposed output, you specify @code{FFTW_MPI_TRANSPOSED_IN} to tell
Chris@10 465 it that the input is transposed. In this section, we explain how to
Chris@10 466 interpret the output format of such a transform.
Chris@10 467 @ctindex FFTW_MPI_TRANSPOSED_OUT
Chris@10 468 @ctindex FFTW_MPI_TRANSPOSED_IN
Chris@10 469
Chris@10 470
Chris@10 471 Suppose you have are transforming multi-dimensional data with (at
Chris@10 472 least two) dimensions @ndims{}. As always, it is distributed along
Chris@10 473 the first dimension @dimk{0}. Now, if we compute its DFT with the
Chris@10 474 @code{FFTW_MPI_TRANSPOSED_OUT} flag, the resulting output data are stored
Chris@10 475 with the first @emph{two} dimensions transposed: @ndimstrans{},
Chris@10 476 distributed along the @dimk{1} dimension. Conversely, if we take the
Chris@10 477 @ndimstrans{} data and transform it with the
Chris@10 478 @code{FFTW_MPI_TRANSPOSED_IN} flag, then the format goes back to the
Chris@10 479 original @ndims{} array.
Chris@10 480
Chris@10 481 There are two ways to find the portion of the transposed array that
Chris@10 482 resides on the current process. First, you can simply call the
Chris@10 483 appropriate @samp{local_size} function, passing @ndimstrans{} (the
Chris@10 484 transposed dimensions). This would mean calling the @samp{local_size}
Chris@10 485 function twice, once for the transposed and once for the
Chris@10 486 non-transposed dimensions. Alternatively, you can call one of the
Chris@10 487 @samp{local_size_transposed} functions, which returns both the
Chris@10 488 non-transposed and transposed data distribution from a single call.
Chris@10 489 For example, for a 3d transform with transposed output (or input), you
Chris@10 490 might call:
Chris@10 491
Chris@10 492 @example
Chris@10 493 ptrdiff_t fftw_mpi_local_size_3d_transposed(
Chris@10 494 ptrdiff_t n0, ptrdiff_t n1, ptrdiff_t n2, MPI_Comm comm,
Chris@10 495 ptrdiff_t *local_n0, ptrdiff_t *local_0_start,
Chris@10 496 ptrdiff_t *local_n1, ptrdiff_t *local_1_start);
Chris@10 497 @end example
Chris@10 498 @findex fftw_mpi_local_size_3d_transposed
Chris@10 499
Chris@10 500 Here, @code{local_n0} and @code{local_0_start} give the size and
Chris@10 501 starting index of the @code{n0} dimension for the
Chris@10 502 @emph{non}-transposed data, as in the previous sections. For
Chris@10 503 @emph{transposed} data (e.g. the output for
Chris@10 504 @code{FFTW_MPI_TRANSPOSED_OUT}), @code{local_n1} and
Chris@10 505 @code{local_1_start} give the size and starting index of the @code{n1}
Chris@10 506 dimension, which is the first dimension of the transposed data
Chris@10 507 (@code{n1} by @code{n0} by @code{n2}).
Chris@10 508
Chris@10 509 (Note that @code{FFTW_MPI_TRANSPOSED_IN} is completely equivalent to
Chris@10 510 performing @code{FFTW_MPI_TRANSPOSED_OUT} and passing the first two
Chris@10 511 dimensions to the planner in reverse order, or vice versa. If you
Chris@10 512 pass @emph{both} the @code{FFTW_MPI_TRANSPOSED_IN} and
Chris@10 513 @code{FFTW_MPI_TRANSPOSED_OUT} flags, it is equivalent to swapping the
Chris@10 514 first two dimensions passed to the planner and passing @emph{neither}
Chris@10 515 flag.)
Chris@10 516
Chris@10 517 @node One-dimensional distributions, , Transposed distributions, MPI Data Distribution
Chris@10 518 @subsection One-dimensional distributions
Chris@10 519
Chris@10 520 For one-dimensional distributed DFTs using FFTW, matters are slightly
Chris@10 521 more complicated because the data distribution is more closely tied to
Chris@10 522 how the algorithm works. In particular, you can no longer pass an
Chris@10 523 arbitrary block size and must accept FFTW's default; also, the block
Chris@10 524 sizes may be different for input and output. Also, the data
Chris@10 525 distribution depends on the flags and transform direction, in order
Chris@10 526 for forward and backward transforms to work correctly.
Chris@10 527
Chris@10 528 @example
Chris@10 529 ptrdiff_t fftw_mpi_local_size_1d(ptrdiff_t n0, MPI_Comm comm,
Chris@10 530 int sign, unsigned flags,
Chris@10 531 ptrdiff_t *local_ni, ptrdiff_t *local_i_start,
Chris@10 532 ptrdiff_t *local_no, ptrdiff_t *local_o_start);
Chris@10 533 @end example
Chris@10 534 @findex fftw_mpi_local_size_1d
Chris@10 535
Chris@10 536 This function computes the data distribution for a 1d transform of
Chris@10 537 size @code{n0} with the given transform @code{sign} and @code{flags}.
Chris@10 538 Both input and output data use block distributions. The input on the
Chris@10 539 current process will consist of @code{local_ni} numbers starting at
Chris@10 540 index @code{local_i_start}; e.g. if only a single process is used,
Chris@10 541 then @code{local_ni} will be @code{n0} and @code{local_i_start} will
Chris@10 542 be @code{0}. Similarly for the output, with @code{local_no} numbers
Chris@10 543 starting at index @code{local_o_start}. The return value of
Chris@10 544 @code{fftw_mpi_local_size_1d} will be the total number of elements to
Chris@10 545 allocate on the current process (which might be slightly larger than
Chris@10 546 the local size due to intermediate steps in the algorithm).
Chris@10 547
Chris@10 548 As mentioned above (@pxref{Load balancing}), the data will be divided
Chris@10 549 equally among the processes if @code{n0} is divisible by the
Chris@10 550 @emph{square} of the number of processes. In this case,
Chris@10 551 @code{local_ni} will equal @code{local_no}. Otherwise, they may be
Chris@10 552 different.
Chris@10 553
Chris@10 554 For some applications, such as convolutions, the order of the output
Chris@10 555 data is irrelevant. In this case, performance can be improved by
Chris@10 556 specifying that the output data be stored in an FFTW-defined
Chris@10 557 ``scrambled'' format. (In particular, this is the analogue of
Chris@10 558 transposed output in the multidimensional case: scrambled output saves
Chris@10 559 a communications step.) If you pass @code{FFTW_MPI_SCRAMBLED_OUT} in
Chris@10 560 the flags, then the output is stored in this (undocumented) scrambled
Chris@10 561 order. Conversely, to perform the inverse transform of data in
Chris@10 562 scrambled order, pass the @code{FFTW_MPI_SCRAMBLED_IN} flag.
Chris@10 563 @ctindex FFTW_MPI_SCRAMBLED_OUT
Chris@10 564 @ctindex FFTW_MPI_SCRAMBLED_IN
Chris@10 565
Chris@10 566
Chris@10 567 In MPI FFTW, only composite sizes @code{n0} can be parallelized; we
Chris@10 568 have not yet implemented a parallel algorithm for large prime sizes.
Chris@10 569
Chris@10 570 @c ------------------------------------------------------------
Chris@10 571 @node Multi-dimensional MPI DFTs of Real Data, Other Multi-dimensional Real-data MPI Transforms, MPI Data Distribution, Distributed-memory FFTW with MPI
Chris@10 572 @section Multi-dimensional MPI DFTs of Real Data
Chris@10 573
Chris@10 574 FFTW's MPI interface also supports multi-dimensional DFTs of real
Chris@10 575 data, similar to the serial r2c and c2r interfaces. (Parallel
Chris@10 576 one-dimensional real-data DFTs are not currently supported; you must
Chris@10 577 use a complex transform and set the imaginary parts of the inputs to
Chris@10 578 zero.)
Chris@10 579
Chris@10 580 The key points to understand for r2c and c2r MPI transforms (compared
Chris@10 581 to the MPI complex DFTs or the serial r2c/c2r transforms), are:
Chris@10 582
Chris@10 583 @itemize @bullet
Chris@10 584
Chris@10 585 @item
Chris@10 586 Just as for serial transforms, r2c/c2r DFTs transform @ndims{} real
Chris@10 587 data to/from @ndimshalf{} complex data: the last dimension of the
Chris@10 588 complex data is cut in half (rounded down), plus one. As for the
Chris@10 589 serial transforms, the sizes you pass to the @samp{plan_dft_r2c} and
Chris@10 590 @samp{plan_dft_c2r} are the @ndims{} dimensions of the real data.
Chris@10 591
Chris@10 592 @item
Chris@10 593 @cindex padding
Chris@10 594 Although the real data is @emph{conceptually} @ndims{}, it is
Chris@10 595 @emph{physically} stored as an @ndimspad{} array, where the last
Chris@10 596 dimension has been @emph{padded} to make it the same size as the
Chris@10 597 complex output. This is much like the in-place serial r2c/c2r
Chris@10 598 interface (@pxref{Multi-Dimensional DFTs of Real Data}), except that
Chris@10 599 in MPI the padding is required even for out-of-place data. The extra
Chris@10 600 padding numbers are ignored by FFTW (they are @emph{not} like
Chris@10 601 zero-padding the transform to a larger size); they are only used to
Chris@10 602 determine the data layout.
Chris@10 603
Chris@10 604 @item
Chris@10 605 @cindex data distribution
Chris@10 606 The data distribution in MPI for @emph{both} the real and complex data
Chris@10 607 is determined by the shape of the @emph{complex} data. That is, you
Chris@10 608 call the appropriate @samp{local size} function for the @ndimshalf{}
Chris@10 609 complex data, and then use the @emph{same} distribution for the real
Chris@10 610 data except that the last complex dimension is replaced by a (padded)
Chris@10 611 real dimension of twice the length.
Chris@10 612
Chris@10 613 @end itemize
Chris@10 614
Chris@10 615 For example suppose we are performing an out-of-place r2c transform of
Chris@10 616 @threedims{L,M,N} real data [padded to @threedims{L,M,2(N/2+1)}],
Chris@10 617 resulting in @threedims{L,M,N/2+1} complex data. Similar to the
Chris@10 618 example in @ref{2d MPI example}, we might do something like:
Chris@10 619
Chris@10 620 @example
Chris@10 621 #include <fftw3-mpi.h>
Chris@10 622
Chris@10 623 int main(int argc, char **argv)
Chris@10 624 @{
Chris@10 625 const ptrdiff_t L = ..., M = ..., N = ...;
Chris@10 626 fftw_plan plan;
Chris@10 627 double *rin;
Chris@10 628 fftw_complex *cout;
Chris@10 629 ptrdiff_t alloc_local, local_n0, local_0_start, i, j, k;
Chris@10 630
Chris@10 631 MPI_Init(&argc, &argv);
Chris@10 632 fftw_mpi_init();
Chris@10 633
Chris@10 634 /* @r{get local data size and allocate} */
Chris@10 635 alloc_local = fftw_mpi_local_size_3d(L, M, N/2+1, MPI_COMM_WORLD,
Chris@10 636 &local_n0, &local_0_start);
Chris@10 637 rin = fftw_alloc_real(2 * alloc_local);
Chris@10 638 cout = fftw_alloc_complex(alloc_local);
Chris@10 639
Chris@10 640 /* @r{create plan for out-of-place r2c DFT} */
Chris@10 641 plan = fftw_mpi_plan_dft_r2c_3d(L, M, N, rin, cout, MPI_COMM_WORLD,
Chris@10 642 FFTW_MEASURE);
Chris@10 643
Chris@10 644 /* @r{initialize rin to some function} my_func(x,y,z) */
Chris@10 645 for (i = 0; i < local_n0; ++i)
Chris@10 646 for (j = 0; j < M; ++j)
Chris@10 647 for (k = 0; k < N; ++k)
Chris@10 648 rin[(i*M + j) * (2*(N/2+1)) + k] = my_func(local_0_start+i, j, k);
Chris@10 649
Chris@10 650 /* @r{compute transforms as many times as desired} */
Chris@10 651 fftw_execute(plan);
Chris@10 652
Chris@10 653 fftw_destroy_plan(plan);
Chris@10 654
Chris@10 655 MPI_Finalize();
Chris@10 656 @}
Chris@10 657 @end example
Chris@10 658
Chris@10 659 @findex fftw_alloc_real
Chris@10 660 @cindex row-major
Chris@10 661 Note that we allocated @code{rin} using @code{fftw_alloc_real} with an
Chris@10 662 argument of @code{2 * alloc_local}: since @code{alloc_local} is the
Chris@10 663 number of @emph{complex} values to allocate, the number of @emph{real}
Chris@10 664 values is twice as many. The @code{rin} array is then
Chris@10 665 @threedims{local_n0,M,2(N/2+1)} in row-major order, so its
Chris@10 666 @code{(i,j,k)} element is at the index @code{(i*M + j) * (2*(N/2+1)) +
Chris@10 667 k} (@pxref{Multi-dimensional Array Format }).
Chris@10 668
Chris@10 669 @cindex transpose
Chris@10 670 @ctindex FFTW_TRANSPOSED_OUT
Chris@10 671 @ctindex FFTW_TRANSPOSED_IN
Chris@10 672 As for the complex transforms, improved performance can be obtained by
Chris@10 673 specifying that the output is the transpose of the input or vice versa
Chris@10 674 (@pxref{Transposed distributions}). In our @threedims{L,M,N} r2c
Chris@10 675 example, including @code{FFTW_TRANSPOSED_OUT} in the flags means that
Chris@10 676 the input would be a padded @threedims{L,M,2(N/2+1)} real array
Chris@10 677 distributed over the @code{L} dimension, while the output would be a
Chris@10 678 @threedims{M,L,N/2+1} complex array distributed over the @code{M}
Chris@10 679 dimension. To perform the inverse c2r transform with the same data
Chris@10 680 distributions, you would use the @code{FFTW_TRANSPOSED_IN} flag.
Chris@10 681
Chris@10 682 @c ------------------------------------------------------------
Chris@10 683 @node Other Multi-dimensional Real-data MPI Transforms, FFTW MPI Transposes, Multi-dimensional MPI DFTs of Real Data, Distributed-memory FFTW with MPI
Chris@10 684 @section Other multi-dimensional Real-Data MPI Transforms
Chris@10 685
Chris@10 686 @cindex r2r
Chris@10 687 FFTW's MPI interface also supports multi-dimensional @samp{r2r}
Chris@10 688 transforms of all kinds supported by the serial interface
Chris@10 689 (e.g. discrete cosine and sine transforms, discrete Hartley
Chris@10 690 transforms, etc.). Only multi-dimensional @samp{r2r} transforms, not
Chris@10 691 one-dimensional transforms, are currently parallelized.
Chris@10 692
Chris@10 693 @tindex fftw_r2r_kind
Chris@10 694 These are used much like the multidimensional complex DFTs discussed
Chris@10 695 above, except that the data is real rather than complex, and one needs
Chris@10 696 to pass an r2r transform kind (@code{fftw_r2r_kind}) for each
Chris@10 697 dimension as in the serial FFTW (@pxref{More DFTs of Real Data}).
Chris@10 698
Chris@10 699 For example, one might perform a two-dimensional @twodims{L,M} that is
Chris@10 700 an REDFT10 (DCT-II) in the first dimension and an RODFT10 (DST-II) in
Chris@10 701 the second dimension with code like:
Chris@10 702
Chris@10 703 @example
Chris@10 704 const ptrdiff_t L = ..., M = ...;
Chris@10 705 fftw_plan plan;
Chris@10 706 double *data;
Chris@10 707 ptrdiff_t alloc_local, local_n0, local_0_start, i, j;
Chris@10 708
Chris@10 709 /* @r{get local data size and allocate} */
Chris@10 710 alloc_local = fftw_mpi_local_size_2d(L, M, MPI_COMM_WORLD,
Chris@10 711 &local_n0, &local_0_start);
Chris@10 712 data = fftw_alloc_real(alloc_local);
Chris@10 713
Chris@10 714 /* @r{create plan for in-place REDFT10 x RODFT10} */
Chris@10 715 plan = fftw_mpi_plan_r2r_2d(L, M, data, data, MPI_COMM_WORLD,
Chris@10 716 FFTW_REDFT10, FFTW_RODFT10, FFTW_MEASURE);
Chris@10 717
Chris@10 718 /* @r{initialize data to some function} my_function(x,y) */
Chris@10 719 for (i = 0; i < local_n0; ++i) for (j = 0; j < M; ++j)
Chris@10 720 data[i*M + j] = my_function(local_0_start + i, j);
Chris@10 721
Chris@10 722 /* @r{compute transforms, in-place, as many times as desired} */
Chris@10 723 fftw_execute(plan);
Chris@10 724
Chris@10 725 fftw_destroy_plan(plan);
Chris@10 726 @end example
Chris@10 727
Chris@10 728 @findex fftw_alloc_real
Chris@10 729 Notice that we use the same @samp{local_size} functions as we did for
Chris@10 730 complex data, only now we interpret the sizes in terms of real rather
Chris@10 731 than complex values, and correspondingly use @code{fftw_alloc_real}.
Chris@10 732
Chris@10 733 @c ------------------------------------------------------------
Chris@10 734 @node FFTW MPI Transposes, FFTW MPI Wisdom, Other Multi-dimensional Real-data MPI Transforms, Distributed-memory FFTW with MPI
Chris@10 735 @section FFTW MPI Transposes
Chris@10 736 @cindex transpose
Chris@10 737
Chris@10 738 The FFTW's MPI Fourier transforms rely on one or more @emph{global
Chris@10 739 transposition} step for their communications. For example, the
Chris@10 740 multidimensional transforms work by transforming along some
Chris@10 741 dimensions, then transposing to make the first dimension local and
Chris@10 742 transforming that, then transposing back. Because global
Chris@10 743 transposition of a block-distributed matrix has many other potential
Chris@10 744 uses besides FFTs, FFTW's transpose routines can be called directly,
Chris@10 745 as documented in this section.
Chris@10 746
Chris@10 747 @menu
Chris@10 748 * Basic distributed-transpose interface::
Chris@10 749 * Advanced distributed-transpose interface::
Chris@10 750 * An improved replacement for MPI_Alltoall::
Chris@10 751 @end menu
Chris@10 752
Chris@10 753 @node Basic distributed-transpose interface, Advanced distributed-transpose interface, FFTW MPI Transposes, FFTW MPI Transposes
Chris@10 754 @subsection Basic distributed-transpose interface
Chris@10 755
Chris@10 756 In particular, suppose that we have an @code{n0} by @code{n1} array in
Chris@10 757 row-major order, block-distributed across the @code{n0} dimension. To
Chris@10 758 transpose this into an @code{n1} by @code{n0} array block-distributed
Chris@10 759 across the @code{n1} dimension, we would create a plan by calling the
Chris@10 760 following function:
Chris@10 761
Chris@10 762 @example
Chris@10 763 fftw_plan fftw_mpi_plan_transpose(ptrdiff_t n0, ptrdiff_t n1,
Chris@10 764 double *in, double *out,
Chris@10 765 MPI_Comm comm, unsigned flags);
Chris@10 766 @end example
Chris@10 767 @findex fftw_mpi_plan_transpose
Chris@10 768
Chris@10 769 The input and output arrays (@code{in} and @code{out}) can be the
Chris@10 770 same. The transpose is actually executed by calling
Chris@10 771 @code{fftw_execute} on the plan, as usual.
Chris@10 772 @findex fftw_execute
Chris@10 773
Chris@10 774
Chris@10 775 The @code{flags} are the usual FFTW planner flags, but support
Chris@10 776 two additional flags: @code{FFTW_MPI_TRANSPOSED_OUT} and/or
Chris@10 777 @code{FFTW_MPI_TRANSPOSED_IN}. What these flags indicate, for
Chris@10 778 transpose plans, is that the output and/or input, respectively, are
Chris@10 779 @emph{locally} transposed. That is, on each process input data is
Chris@10 780 normally stored as a @code{local_n0} by @code{n1} array in row-major
Chris@10 781 order, but for an @code{FFTW_MPI_TRANSPOSED_IN} plan the input data is
Chris@10 782 stored as @code{n1} by @code{local_n0} in row-major order. Similarly,
Chris@10 783 @code{FFTW_MPI_TRANSPOSED_OUT} means that the output is @code{n0} by
Chris@10 784 @code{local_n1} instead of @code{local_n1} by @code{n0}.
Chris@10 785 @ctindex FFTW_MPI_TRANSPOSED_OUT
Chris@10 786 @ctindex FFTW_MPI_TRANSPOSED_IN
Chris@10 787
Chris@10 788
Chris@10 789 To determine the local size of the array on each process before and
Chris@10 790 after the transpose, as well as the amount of storage that must be
Chris@10 791 allocated, one should call @code{fftw_mpi_local_size_2d_transposed},
Chris@10 792 just as for a 2d DFT as described in the previous section:
Chris@10 793 @cindex data distribution
Chris@10 794
Chris@10 795 @example
Chris@10 796 ptrdiff_t fftw_mpi_local_size_2d_transposed
Chris@10 797 (ptrdiff_t n0, ptrdiff_t n1, MPI_Comm comm,
Chris@10 798 ptrdiff_t *local_n0, ptrdiff_t *local_0_start,
Chris@10 799 ptrdiff_t *local_n1, ptrdiff_t *local_1_start);
Chris@10 800 @end example
Chris@10 801 @findex fftw_mpi_local_size_2d_transposed
Chris@10 802
Chris@10 803 Again, the return value is the local storage to allocate, which in
Chris@10 804 this case is the number of @emph{real} (@code{double}) values rather
Chris@10 805 than complex numbers as in the previous examples.
Chris@10 806
Chris@10 807 @node Advanced distributed-transpose interface, An improved replacement for MPI_Alltoall, Basic distributed-transpose interface, FFTW MPI Transposes
Chris@10 808 @subsection Advanced distributed-transpose interface
Chris@10 809
Chris@10 810 The above routines are for a transpose of a matrix of numbers (of type
Chris@10 811 @code{double}), using FFTW's default block sizes. More generally, one
Chris@10 812 can perform transposes of @emph{tuples} of numbers, with
Chris@10 813 user-specified block sizes for the input and output:
Chris@10 814
Chris@10 815 @example
Chris@10 816 fftw_plan fftw_mpi_plan_many_transpose
Chris@10 817 (ptrdiff_t n0, ptrdiff_t n1, ptrdiff_t howmany,
Chris@10 818 ptrdiff_t block0, ptrdiff_t block1,
Chris@10 819 double *in, double *out, MPI_Comm comm, unsigned flags);
Chris@10 820 @end example
Chris@10 821 @findex fftw_mpi_plan_many_transpose
Chris@10 822
Chris@10 823 In this case, one is transposing an @code{n0} by @code{n1} matrix of
Chris@10 824 @code{howmany}-tuples (e.g. @code{howmany = 2} for complex numbers).
Chris@10 825 The input is distributed along the @code{n0} dimension with block size
Chris@10 826 @code{block0}, and the @code{n1} by @code{n0} output is distributed
Chris@10 827 along the @code{n1} dimension with block size @code{block1}. If
Chris@10 828 @code{FFTW_MPI_DEFAULT_BLOCK} (0) is passed for a block size then FFTW
Chris@10 829 uses its default block size. To get the local size of the data on
Chris@10 830 each process, you should then call @code{fftw_mpi_local_size_many_transposed}.
Chris@10 831 @ctindex FFTW_MPI_DEFAULT_BLOCK
Chris@10 832 @findex fftw_mpi_local_size_many_transposed
Chris@10 833
Chris@10 834 @node An improved replacement for MPI_Alltoall, , Advanced distributed-transpose interface, FFTW MPI Transposes
Chris@10 835 @subsection An improved replacement for MPI_Alltoall
Chris@10 836
Chris@10 837 We close this section by noting that FFTW's MPI transpose routines can
Chris@10 838 be thought of as a generalization for the @code{MPI_Alltoall} function
Chris@10 839 (albeit only for floating-point types), and in some circumstances can
Chris@10 840 function as an improved replacement.
Chris@10 841 @findex MPI_Alltoall
Chris@10 842
Chris@10 843
Chris@10 844 @code{MPI_Alltoall} is defined by the MPI standard as:
Chris@10 845
Chris@10 846 @example
Chris@10 847 int MPI_Alltoall(void *sendbuf, int sendcount, MPI_Datatype sendtype,
Chris@10 848 void *recvbuf, int recvcnt, MPI_Datatype recvtype,
Chris@10 849 MPI_Comm comm);
Chris@10 850 @end example
Chris@10 851
Chris@10 852 In particular, for @code{double*} arrays @code{in} and @code{out},
Chris@10 853 consider the call:
Chris@10 854
Chris@10 855 @example
Chris@10 856 MPI_Alltoall(in, howmany, MPI_DOUBLE, out, howmany MPI_DOUBLE, comm);
Chris@10 857 @end example
Chris@10 858
Chris@10 859 This is completely equivalent to:
Chris@10 860
Chris@10 861 @example
Chris@10 862 MPI_Comm_size(comm, &P);
Chris@10 863 plan = fftw_mpi_plan_many_transpose(P, P, howmany, 1, 1, in, out, comm, FFTW_ESTIMATE);
Chris@10 864 fftw_execute(plan);
Chris@10 865 fftw_destroy_plan(plan);
Chris@10 866 @end example
Chris@10 867
Chris@10 868 That is, computing a @twodims{P,P} transpose on @code{P} processes,
Chris@10 869 with a block size of 1, is just a standard all-to-all communication.
Chris@10 870
Chris@10 871 However, using the FFTW routine instead of @code{MPI_Alltoall} may
Chris@10 872 have certain advantages. First of all, FFTW's routine can operate
Chris@10 873 in-place (@code{in == out}) whereas @code{MPI_Alltoall} can only
Chris@10 874 operate out-of-place.
Chris@10 875 @cindex in-place
Chris@10 876
Chris@10 877
Chris@10 878 Second, even for out-of-place plans, FFTW's routine may be faster,
Chris@10 879 especially if you need to perform the all-to-all communication many
Chris@10 880 times and can afford to use @code{FFTW_MEASURE} or
Chris@10 881 @code{FFTW_PATIENT}. It should certainly be no slower, not including
Chris@10 882 the time to create the plan, since one of the possible algorithms that
Chris@10 883 FFTW uses for an out-of-place transpose @emph{is} simply to call
Chris@10 884 @code{MPI_Alltoall}. However, FFTW also considers several other
Chris@10 885 possible algorithms that, depending on your MPI implementation and
Chris@10 886 your hardware, may be faster.
Chris@10 887 @ctindex FFTW_MEASURE
Chris@10 888 @ctindex FFTW_PATIENT
Chris@10 889
Chris@10 890 @c ------------------------------------------------------------
Chris@10 891 @node FFTW MPI Wisdom, Avoiding MPI Deadlocks, FFTW MPI Transposes, Distributed-memory FFTW with MPI
Chris@10 892 @section FFTW MPI Wisdom
Chris@10 893 @cindex wisdom
Chris@10 894 @cindex saving plans to disk
Chris@10 895
Chris@10 896 FFTW's ``wisdom'' facility (@pxref{Words of Wisdom-Saving Plans}) can
Chris@10 897 be used to save MPI plans as well as to save uniprocessor plans.
Chris@10 898 However, for MPI there are several unavoidable complications.
Chris@10 899
Chris@10 900 @cindex MPI I/O
Chris@10 901 First, the MPI standard does not guarantee that every process can
Chris@10 902 perform file I/O (at least, not using C stdio routines)---in general,
Chris@10 903 we may only assume that process 0 is capable of I/O.@footnote{In fact,
Chris@10 904 even this assumption is not technically guaranteed by the standard,
Chris@10 905 although it seems to be universal in actual MPI implementations and is
Chris@10 906 widely assumed by MPI-using software. Technically, you need to query
Chris@10 907 the @code{MPI_IO} attribute of @code{MPI_COMM_WORLD} with
Chris@10 908 @code{MPI_Attr_get}. If this attribute is @code{MPI_PROC_NULL}, no
Chris@10 909 I/O is possible. If it is @code{MPI_ANY_SOURCE}, any process can
Chris@10 910 perform I/O. Otherwise, it is the rank of a process that can perform
Chris@10 911 I/O ... but since it is not guaranteed to yield the @emph{same} rank
Chris@10 912 on all processes, you have to do an @code{MPI_Allreduce} of some kind
Chris@10 913 if you want all processes to agree about which is going to do I/O.
Chris@10 914 And even then, the standard only guarantees that this process can
Chris@10 915 perform output, but not input. See e.g. @cite{Parallel Programming
Chris@10 916 with MPI} by P. S. Pacheco, section 8.1.3. Needless to say, in our
Chris@10 917 experience virtually no MPI programmers worry about this.} So, if we
Chris@10 918 want to export the wisdom from a single process to a file, we must
Chris@10 919 first export the wisdom to a string, then send it to process 0, then
Chris@10 920 write it to a file.
Chris@10 921
Chris@10 922 Second, in principle we may want to have separate wisdom for every
Chris@10 923 process, since in general the processes may run on different hardware
Chris@10 924 even for a single MPI program. However, in practice FFTW's MPI code
Chris@10 925 is designed for the case of homogeneous hardware (@pxref{Load
Chris@10 926 balancing}), and in this case it is convenient to use the same wisdom
Chris@10 927 for every process. Thus, we need a mechanism to synchronize the wisdom.
Chris@10 928
Chris@10 929 To address both of these problems, FFTW provides the following two
Chris@10 930 functions:
Chris@10 931
Chris@10 932 @example
Chris@10 933 void fftw_mpi_broadcast_wisdom(MPI_Comm comm);
Chris@10 934 void fftw_mpi_gather_wisdom(MPI_Comm comm);
Chris@10 935 @end example
Chris@10 936 @findex fftw_mpi_gather_wisdom
Chris@10 937 @findex fftw_mpi_broadcast_wisdom
Chris@10 938
Chris@10 939 Given a communicator @code{comm}, @code{fftw_mpi_broadcast_wisdom}
Chris@10 940 will broadcast the wisdom from process 0 to all other processes.
Chris@10 941 Conversely, @code{fftw_mpi_gather_wisdom} will collect wisdom from all
Chris@10 942 processes onto process 0. (If the plans created for the same problem
Chris@10 943 by different processes are not the same, @code{fftw_mpi_gather_wisdom}
Chris@10 944 will arbitrarily choose one of the plans.) Both of these functions
Chris@10 945 may result in suboptimal plans for different processes if the
Chris@10 946 processes are running on non-identical hardware. Both of these
Chris@10 947 functions are @emph{collective} calls, which means that they must be
Chris@10 948 executed by all processes in the communicator.
Chris@10 949 @cindex collective function
Chris@10 950
Chris@10 951
Chris@10 952 So, for example, a typical code snippet to import wisdom from a file
Chris@10 953 and use it on all processes would be:
Chris@10 954
Chris@10 955 @example
Chris@10 956 @{
Chris@10 957 int rank;
Chris@10 958
Chris@10 959 fftw_mpi_init();
Chris@10 960 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
Chris@10 961 if (rank == 0) fftw_import_wisdom_from_filename("mywisdom");
Chris@10 962 fftw_mpi_broadcast_wisdom(MPI_COMM_WORLD);
Chris@10 963 @}
Chris@10 964 @end example
Chris@10 965
Chris@10 966 (Note that we must call @code{fftw_mpi_init} before importing any
Chris@10 967 wisdom that might contain MPI plans.) Similarly, a typical code
Chris@10 968 snippet to export wisdom from all processes to a file is:
Chris@10 969 @findex fftw_mpi_init
Chris@10 970
Chris@10 971 @example
Chris@10 972 @{
Chris@10 973 int rank;
Chris@10 974
Chris@10 975 fftw_mpi_gather_wisdom(MPI_COMM_WORLD);
Chris@10 976 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
Chris@10 977 if (rank == 0) fftw_export_wisdom_to_filename("mywisdom");
Chris@10 978 @}
Chris@10 979 @end example
Chris@10 980
Chris@10 981 @c ------------------------------------------------------------
Chris@10 982 @node Avoiding MPI Deadlocks, FFTW MPI Performance Tips, FFTW MPI Wisdom, Distributed-memory FFTW with MPI
Chris@10 983 @section Avoiding MPI Deadlocks
Chris@10 984 @cindex deadlock
Chris@10 985
Chris@10 986 An MPI program can @emph{deadlock} if one process is waiting for a
Chris@10 987 message from another process that never gets sent. To avoid deadlocks
Chris@10 988 when using FFTW's MPI routines, it is important to know which
Chris@10 989 functions are @emph{collective}: that is, which functions must
Chris@10 990 @emph{always} be called in the @emph{same order} from @emph{every}
Chris@10 991 process in a given communicator. (For example, @code{MPI_Barrier} is
Chris@10 992 the canonical example of a collective function in the MPI standard.)
Chris@10 993 @cindex collective function
Chris@10 994 @findex MPI_Barrier
Chris@10 995
Chris@10 996
Chris@10 997 The functions in FFTW that are @emph{always} collective are: every
Chris@10 998 function beginning with @samp{fftw_mpi_plan}, as well as
Chris@10 999 @code{fftw_mpi_broadcast_wisdom} and @code{fftw_mpi_gather_wisdom}.
Chris@10 1000 Also, the following functions from the ordinary FFTW interface are
Chris@10 1001 collective when they are applied to a plan created by an
Chris@10 1002 @samp{fftw_mpi_plan} function: @code{fftw_execute},
Chris@10 1003 @code{fftw_destroy_plan}, and @code{fftw_flops}.
Chris@10 1004 @findex fftw_execute
Chris@10 1005 @findex fftw_destroy_plan
Chris@10 1006 @findex fftw_flops
Chris@10 1007
Chris@10 1008 @c ------------------------------------------------------------
Chris@10 1009 @node FFTW MPI Performance Tips, Combining MPI and Threads, Avoiding MPI Deadlocks, Distributed-memory FFTW with MPI
Chris@10 1010 @section FFTW MPI Performance Tips
Chris@10 1011
Chris@10 1012 In this section, we collect a few tips on getting the best performance
Chris@10 1013 out of FFTW's MPI transforms.
Chris@10 1014
Chris@10 1015 First, because of the 1d block distribution, FFTW's parallelization is
Chris@10 1016 currently limited by the size of the first dimension.
Chris@10 1017 (Multidimensional block distributions may be supported by a future
Chris@10 1018 version.) More generally, you should ideally arrange the dimensions so
Chris@10 1019 that FFTW can divide them equally among the processes. @xref{Load
Chris@10 1020 balancing}.
Chris@10 1021 @cindex block distribution
Chris@10 1022 @cindex load balancing
Chris@10 1023
Chris@10 1024
Chris@10 1025 Second, if it is not too inconvenient, you should consider working
Chris@10 1026 with transposed output for multidimensional plans, as this saves a
Chris@10 1027 considerable amount of communications. @xref{Transposed distributions}.
Chris@10 1028 @cindex transpose
Chris@10 1029
Chris@10 1030
Chris@10 1031 Third, the fastest choices are generally either an in-place transform
Chris@10 1032 or an out-of-place transform with the @code{FFTW_DESTROY_INPUT} flag
Chris@10 1033 (which allows the input array to be used as scratch space). In-place
Chris@10 1034 is especially beneficial if the amount of data per process is large.
Chris@10 1035 @ctindex FFTW_DESTROY_INPUT
Chris@10 1036
Chris@10 1037
Chris@10 1038 Fourth, if you have multiple arrays to transform at once, rather than
Chris@10 1039 calling FFTW's MPI transforms several times it usually seems to be
Chris@10 1040 faster to interleave the data and use the advanced interface. (This
Chris@10 1041 groups the communications together instead of requiring separate
Chris@10 1042 messages for each transform.)
Chris@10 1043
Chris@10 1044 @c ------------------------------------------------------------
Chris@10 1045 @node Combining MPI and Threads, FFTW MPI Reference, FFTW MPI Performance Tips, Distributed-memory FFTW with MPI
Chris@10 1046 @section Combining MPI and Threads
Chris@10 1047 @cindex threads
Chris@10 1048
Chris@10 1049 In certain cases, it may be advantageous to combine MPI
Chris@10 1050 (distributed-memory) and threads (shared-memory) parallelization.
Chris@10 1051 FFTW supports this, with certain caveats. For example, if you have a
Chris@10 1052 cluster of 4-processor shared-memory nodes, you may want to use
Chris@10 1053 threads within the nodes and MPI between the nodes, instead of MPI for
Chris@10 1054 all parallelization.
Chris@10 1055
Chris@10 1056 In particular, it is possible to seamlessly combine the MPI FFTW
Chris@10 1057 routines with the multi-threaded FFTW routines (@pxref{Multi-threaded
Chris@10 1058 FFTW}). However, some care must be taken in the initialization code,
Chris@10 1059 which should look something like this:
Chris@10 1060
Chris@10 1061 @example
Chris@10 1062 int threads_ok;
Chris@10 1063
Chris@10 1064 int main(int argc, char **argv)
Chris@10 1065 @{
Chris@10 1066 int provided;
Chris@10 1067 MPI_Init_thread(&argc, &argv, MPI_THREAD_FUNNELED, &provided);
Chris@10 1068 threads_ok = provided >= MPI_THREAD_FUNNELED;
Chris@10 1069
Chris@10 1070 if (threads_ok) threads_ok = fftw_init_threads();
Chris@10 1071 fftw_mpi_init();
Chris@10 1072
Chris@10 1073 ...
Chris@10 1074 if (threads_ok) fftw_plan_with_nthreads(...);
Chris@10 1075 ...
Chris@10 1076
Chris@10 1077 MPI_Finalize();
Chris@10 1078 @}
Chris@10 1079 @end example
Chris@10 1080 @findex fftw_mpi_init
Chris@10 1081 @findex fftw_init_threads
Chris@10 1082 @findex fftw_plan_with_nthreads
Chris@10 1083
Chris@10 1084 First, note that instead of calling @code{MPI_Init}, you should call
Chris@10 1085 @code{MPI_Init_threads}, which is the initialization routine defined
Chris@10 1086 by the MPI-2 standard to indicate to MPI that your program will be
Chris@10 1087 multithreaded. We pass @code{MPI_THREAD_FUNNELED}, which indicates
Chris@10 1088 that we will only call MPI routines from the main thread. (FFTW will
Chris@10 1089 launch additional threads internally, but the extra threads will not
Chris@10 1090 call MPI code.) (You may also pass @code{MPI_THREAD_SERIALIZED} or
Chris@10 1091 @code{MPI_THREAD_MULTIPLE}, which requests additional multithreading
Chris@10 1092 support from the MPI implementation, but this is not required by
Chris@10 1093 FFTW.) The @code{provided} parameter returns what level of threads
Chris@10 1094 support is actually supported by your MPI implementation; this
Chris@10 1095 @emph{must} be at least @code{MPI_THREAD_FUNNELED} if you want to call
Chris@10 1096 the FFTW threads routines, so we define a global variable
Chris@10 1097 @code{threads_ok} to record this. You should only call
Chris@10 1098 @code{fftw_init_threads} or @code{fftw_plan_with_nthreads} if
Chris@10 1099 @code{threads_ok} is true. For more information on thread safety in
Chris@10 1100 MPI, see the
Chris@10 1101 @uref{http://www.mpi-forum.org/docs/mpi-20-html/node162.htm, MPI and
Chris@10 1102 Threads} section of the MPI-2 standard.
Chris@10 1103 @cindex thread safety
Chris@10 1104
Chris@10 1105
Chris@10 1106 Second, we must call @code{fftw_init_threads} @emph{before}
Chris@10 1107 @code{fftw_mpi_init}. This is critical for technical reasons having
Chris@10 1108 to do with how FFTW initializes its list of algorithms.
Chris@10 1109
Chris@10 1110 Then, if you call @code{fftw_plan_with_nthreads(N)}, @emph{every} MPI
Chris@10 1111 process will launch (up to) @code{N} threads to parallelize its transforms.
Chris@10 1112
Chris@10 1113 For example, in the hypothetical cluster of 4-processor nodes, you
Chris@10 1114 might wish to launch only a single MPI process per node, and then call
Chris@10 1115 @code{fftw_plan_with_nthreads(4)} on each process to use all
Chris@10 1116 processors in the nodes.
Chris@10 1117
Chris@10 1118 This may or may not be faster than simply using as many MPI processes
Chris@10 1119 as you have processors, however. On the one hand, using threads
Chris@10 1120 within a node eliminates the need for explicit message passing within
Chris@10 1121 the node. On the other hand, FFTW's transpose routines are not
Chris@10 1122 multi-threaded, and this means that the communications that do take
Chris@10 1123 place will not benefit from parallelization within the node.
Chris@10 1124 Moreover, many MPI implementations already have optimizations to
Chris@10 1125 exploit shared memory when it is available, so adding the
Chris@10 1126 multithreaded FFTW on top of this may be superfluous.
Chris@10 1127 @cindex transpose
Chris@10 1128
Chris@10 1129 @c ------------------------------------------------------------
Chris@10 1130 @node FFTW MPI Reference, FFTW MPI Fortran Interface, Combining MPI and Threads, Distributed-memory FFTW with MPI
Chris@10 1131 @section FFTW MPI Reference
Chris@10 1132
Chris@10 1133 This chapter provides a complete reference to all FFTW MPI functions,
Chris@10 1134 datatypes, and constants. See also @ref{FFTW Reference} for information
Chris@10 1135 on functions and types in common with the serial interface.
Chris@10 1136
Chris@10 1137 @menu
Chris@10 1138 * MPI Files and Data Types::
Chris@10 1139 * MPI Initialization::
Chris@10 1140 * Using MPI Plans::
Chris@10 1141 * MPI Data Distribution Functions::
Chris@10 1142 * MPI Plan Creation::
Chris@10 1143 * MPI Wisdom Communication::
Chris@10 1144 @end menu
Chris@10 1145
Chris@10 1146 @node MPI Files and Data Types, MPI Initialization, FFTW MPI Reference, FFTW MPI Reference
Chris@10 1147 @subsection MPI Files and Data Types
Chris@10 1148
Chris@10 1149 All programs using FFTW's MPI support should include its header file:
Chris@10 1150
Chris@10 1151 @example
Chris@10 1152 #include <fftw3-mpi.h>
Chris@10 1153 @end example
Chris@10 1154
Chris@10 1155 Note that this header file includes the serial-FFTW @code{fftw3.h}
Chris@10 1156 header file, and also the @code{mpi.h} header file for MPI, so you
Chris@10 1157 need not include those files separately.
Chris@10 1158
Chris@10 1159 You must also link to @emph{both} the FFTW MPI library and to the
Chris@10 1160 serial FFTW library. On Unix, this means adding @code{-lfftw3_mpi
Chris@10 1161 -lfftw3 -lm} at the end of the link command.
Chris@10 1162
Chris@10 1163 @cindex precision
Chris@10 1164 Different precisions are handled as in the serial interface:
Chris@10 1165 @xref{Precision}. That is, @samp{fftw_} functions become
Chris@10 1166 @code{fftwf_} (in single precision) etcetera, and the libraries become
Chris@10 1167 @code{-lfftw3f_mpi -lfftw3f -lm} etcetera on Unix. Long-double
Chris@10 1168 precision is supported in MPI, but quad precision (@samp{fftwq_}) is
Chris@10 1169 not due to the lack of MPI support for this type.
Chris@10 1170
Chris@10 1171 @node MPI Initialization, Using MPI Plans, MPI Files and Data Types, FFTW MPI Reference
Chris@10 1172 @subsection MPI Initialization
Chris@10 1173
Chris@10 1174 Before calling any other FFTW MPI (@samp{fftw_mpi_}) function, and
Chris@10 1175 before importing any wisdom for MPI problems, you must call:
Chris@10 1176
Chris@10 1177 @findex fftw_mpi_init
Chris@10 1178 @example
Chris@10 1179 void fftw_mpi_init(void);
Chris@10 1180 @end example
Chris@10 1181
Chris@10 1182 @findex fftw_init_threads
Chris@10 1183 If FFTW threads support is used, however, @code{fftw_mpi_init} should
Chris@10 1184 be called @emph{after} @code{fftw_init_threads} (@pxref{Combining MPI
Chris@10 1185 and Threads}). Calling @code{fftw_mpi_init} additional times (before
Chris@10 1186 @code{fftw_mpi_cleanup}) has no effect.
Chris@10 1187
Chris@10 1188
Chris@10 1189 If you want to deallocate all persistent data and reset FFTW to the
Chris@10 1190 pristine state it was in when you started your program, you can call:
Chris@10 1191
Chris@10 1192 @findex fftw_mpi_cleanup
Chris@10 1193 @example
Chris@10 1194 void fftw_mpi_cleanup(void);
Chris@10 1195 @end example
Chris@10 1196
Chris@10 1197 @findex fftw_cleanup
Chris@10 1198 (This calls @code{fftw_cleanup}, so you need not call the serial
Chris@10 1199 cleanup routine too, although it is safe to do so.) After calling
Chris@10 1200 @code{fftw_mpi_cleanup}, all existing plans become undefined, and you
Chris@10 1201 should not attempt to execute or destroy them. You must call
Chris@10 1202 @code{fftw_mpi_init} again after @code{fftw_mpi_cleanup} if you want
Chris@10 1203 to resume using the MPI FFTW routines.
Chris@10 1204
Chris@10 1205 @node Using MPI Plans, MPI Data Distribution Functions, MPI Initialization, FFTW MPI Reference
Chris@10 1206 @subsection Using MPI Plans
Chris@10 1207
Chris@10 1208 Once an MPI plan is created, you can execute and destroy it using
Chris@10 1209 @code{fftw_execute}, @code{fftw_destroy_plan}, and the other functions
Chris@10 1210 in the serial interface that operate on generic plans (@pxref{Using
Chris@10 1211 Plans}).
Chris@10 1212
Chris@10 1213 @cindex collective function
Chris@10 1214 @cindex MPI communicator
Chris@10 1215 The @code{fftw_execute} and @code{fftw_destroy_plan} functions, applied to
Chris@10 1216 MPI plans, are @emph{collective} calls: they must be called for all processes
Chris@10 1217 in the communicator that was used to create the plan.
Chris@10 1218
Chris@10 1219 @cindex new-array execution
Chris@10 1220 You must @emph{not} use the serial new-array plan-execution functions
Chris@10 1221 @code{fftw_execute_dft} and so on (@pxref{New-array Execute
Chris@10 1222 Functions}) with MPI plans. Such functions are specialized to the
Chris@10 1223 problem type, and there are specific new-array execute functions for MPI plans:
Chris@10 1224
Chris@10 1225 @findex fftw_mpi_execute_dft
Chris@10 1226 @findex fftw_mpi_execute_dft_r2c
Chris@10 1227 @findex fftw_mpi_execute_dft_c2r
Chris@10 1228 @findex fftw_mpi_execute_r2r
Chris@10 1229 @example
Chris@10 1230 void fftw_mpi_execute_dft(fftw_plan p, fftw_complex *in, fftw_complex *out);
Chris@10 1231 void fftw_mpi_execute_dft_r2c(fftw_plan p, double *in, fftw_complex *out);
Chris@10 1232 void fftw_mpi_execute_dft_c2r(fftw_plan p, fftw_complex *in, double *out);
Chris@10 1233 void fftw_mpi_execute_r2r(fftw_plan p, double *in, double *out);
Chris@10 1234 @end example
Chris@10 1235
Chris@10 1236 @cindex alignment
Chris@10 1237 @findex fftw_malloc
Chris@10 1238 These functions have the same restrictions as those of the serial
Chris@10 1239 new-array execute functions. They are @emph{always} safe to apply to
Chris@10 1240 the @emph{same} @code{in} and @code{out} arrays that were used to
Chris@10 1241 create the plan. They can only be applied to new arrarys if those
Chris@10 1242 arrays have the same types, dimensions, in-placeness, and alignment as
Chris@10 1243 the original arrays, where the best way to ensure the same alignment
Chris@10 1244 is to use FFTW's @code{fftw_malloc} and related allocation functions
Chris@10 1245 for all arrays (@pxref{Memory Allocation}). Note that distributed
Chris@10 1246 transposes (@pxref{FFTW MPI Transposes}) use
Chris@10 1247 @code{fftw_mpi_execute_r2r}, since they count as rank-zero r2r plans
Chris@10 1248 from FFTW's perspective.
Chris@10 1249
Chris@10 1250 @node MPI Data Distribution Functions, MPI Plan Creation, Using MPI Plans, FFTW MPI Reference
Chris@10 1251 @subsection MPI Data Distribution Functions
Chris@10 1252
Chris@10 1253 @cindex data distribution
Chris@10 1254 As described above (@pxref{MPI Data Distribution}), in order to
Chris@10 1255 allocate your arrays, @emph{before} creating a plan, you must first
Chris@10 1256 call one of the following routines to determine the required
Chris@10 1257 allocation size and the portion of the array locally stored on a given
Chris@10 1258 process. The @code{MPI_Comm} communicator passed here must be
Chris@10 1259 equivalent to the communicator used below for plan creation.
Chris@10 1260
Chris@10 1261 The basic interface for multidimensional transforms consists of the
Chris@10 1262 functions:
Chris@10 1263
Chris@10 1264 @findex fftw_mpi_local_size_2d
Chris@10 1265 @findex fftw_mpi_local_size_3d
Chris@10 1266 @findex fftw_mpi_local_size
Chris@10 1267 @findex fftw_mpi_local_size_2d_transposed
Chris@10 1268 @findex fftw_mpi_local_size_3d_transposed
Chris@10 1269 @findex fftw_mpi_local_size_transposed
Chris@10 1270 @example
Chris@10 1271 ptrdiff_t fftw_mpi_local_size_2d(ptrdiff_t n0, ptrdiff_t n1, MPI_Comm comm,
Chris@10 1272 ptrdiff_t *local_n0, ptrdiff_t *local_0_start);
Chris@10 1273 ptrdiff_t fftw_mpi_local_size_3d(ptrdiff_t n0, ptrdiff_t n1, ptrdiff_t n2,
Chris@10 1274 MPI_Comm comm,
Chris@10 1275 ptrdiff_t *local_n0, ptrdiff_t *local_0_start);
Chris@10 1276 ptrdiff_t fftw_mpi_local_size(int rnk, const ptrdiff_t *n, MPI_Comm comm,
Chris@10 1277 ptrdiff_t *local_n0, ptrdiff_t *local_0_start);
Chris@10 1278
Chris@10 1279 ptrdiff_t fftw_mpi_local_size_2d_transposed(ptrdiff_t n0, ptrdiff_t n1, MPI_Comm comm,
Chris@10 1280 ptrdiff_t *local_n0, ptrdiff_t *local_0_start,
Chris@10 1281 ptrdiff_t *local_n1, ptrdiff_t *local_1_start);
Chris@10 1282 ptrdiff_t fftw_mpi_local_size_3d_transposed(ptrdiff_t n0, ptrdiff_t n1, ptrdiff_t n2,
Chris@10 1283 MPI_Comm comm,
Chris@10 1284 ptrdiff_t *local_n0, ptrdiff_t *local_0_start,
Chris@10 1285 ptrdiff_t *local_n1, ptrdiff_t *local_1_start);
Chris@10 1286 ptrdiff_t fftw_mpi_local_size_transposed(int rnk, const ptrdiff_t *n, MPI_Comm comm,
Chris@10 1287 ptrdiff_t *local_n0, ptrdiff_t *local_0_start,
Chris@10 1288 ptrdiff_t *local_n1, ptrdiff_t *local_1_start);
Chris@10 1289 @end example
Chris@10 1290
Chris@10 1291 These functions return the number of elements to allocate (complex
Chris@10 1292 numbers for DFT/r2c/c2r plans, real numbers for r2r plans), whereas
Chris@10 1293 the @code{local_n0} and @code{local_0_start} return the portion
Chris@10 1294 (@code{local_0_start} to @code{local_0_start + local_n0 - 1}) of the
Chris@10 1295 first dimension of an @ndims{} array that is stored on the local
Chris@10 1296 process. @xref{Basic and advanced distribution interfaces}. For
Chris@10 1297 @code{FFTW_MPI_TRANSPOSED_OUT} plans, the @samp{_transposed} variants
Chris@10 1298 are useful in order to also return the local portion of the first
Chris@10 1299 dimension in the @ndimstrans{} transposed output. @xref{Transposed
Chris@10 1300 distributions}. The advanced interface for multidimensional
Chris@10 1301 transforms is:
Chris@10 1302
Chris@10 1303 @cindex advanced interface
Chris@10 1304 @findex fftw_mpi_local_size_many
Chris@10 1305 @findex fftw_mpi_local_size_many_transposed
Chris@10 1306 @example
Chris@10 1307 ptrdiff_t fftw_mpi_local_size_many(int rnk, const ptrdiff_t *n, ptrdiff_t howmany,
Chris@10 1308 ptrdiff_t block0, MPI_Comm comm,
Chris@10 1309 ptrdiff_t *local_n0, ptrdiff_t *local_0_start);
Chris@10 1310 ptrdiff_t fftw_mpi_local_size_many_transposed(int rnk, const ptrdiff_t *n, ptrdiff_t howmany,
Chris@10 1311 ptrdiff_t block0, ptrdiff_t block1, MPI_Comm comm,
Chris@10 1312 ptrdiff_t *local_n0, ptrdiff_t *local_0_start,
Chris@10 1313 ptrdiff_t *local_n1, ptrdiff_t *local_1_start);
Chris@10 1314 @end example
Chris@10 1315
Chris@10 1316 These differ from the basic interface in only two ways. First, they
Chris@10 1317 allow you to specify block sizes @code{block0} and @code{block1} (the
Chris@10 1318 latter for the transposed output); you can pass
Chris@10 1319 @code{FFTW_MPI_DEFAULT_BLOCK} to use FFTW's default block size as in
Chris@10 1320 the basic interface. Second, you can pass a @code{howmany} parameter,
Chris@10 1321 corresponding to the advanced planning interface below: this is for
Chris@10 1322 transforms of contiguous @code{howmany}-tuples of numbers
Chris@10 1323 (@code{howmany = 1} in the basic interface).
Chris@10 1324
Chris@10 1325 The corresponding basic and advanced routines for one-dimensional
Chris@10 1326 transforms (currently only complex DFTs) are:
Chris@10 1327
Chris@10 1328 @findex fftw_mpi_local_size_1d
Chris@10 1329 @findex fftw_mpi_local_size_many_1d
Chris@10 1330 @example
Chris@10 1331 ptrdiff_t fftw_mpi_local_size_1d(
Chris@10 1332 ptrdiff_t n0, MPI_Comm comm, int sign, unsigned flags,
Chris@10 1333 ptrdiff_t *local_ni, ptrdiff_t *local_i_start,
Chris@10 1334 ptrdiff_t *local_no, ptrdiff_t *local_o_start);
Chris@10 1335 ptrdiff_t fftw_mpi_local_size_many_1d(
Chris@10 1336 ptrdiff_t n0, ptrdiff_t howmany,
Chris@10 1337 MPI_Comm comm, int sign, unsigned flags,
Chris@10 1338 ptrdiff_t *local_ni, ptrdiff_t *local_i_start,
Chris@10 1339 ptrdiff_t *local_no, ptrdiff_t *local_o_start);
Chris@10 1340 @end example
Chris@10 1341
Chris@10 1342 @ctindex FFTW_MPI_SCRAMBLED_OUT
Chris@10 1343 @ctindex FFTW_MPI_SCRAMBLED_IN
Chris@10 1344 As above, the return value is the number of elements to allocate
Chris@10 1345 (complex numbers, for complex DFTs). The @code{local_ni} and
Chris@10 1346 @code{local_i_start} arguments return the portion
Chris@10 1347 (@code{local_i_start} to @code{local_i_start + local_ni - 1}) of the
Chris@10 1348 1d array that is stored on this process for the transform
Chris@10 1349 @emph{input}, and @code{local_no} and @code{local_o_start} are the
Chris@10 1350 corresponding quantities for the input. The @code{sign}
Chris@10 1351 (@code{FFTW_FORWARD} or @code{FFTW_BACKWARD}) and @code{flags} must
Chris@10 1352 match the arguments passed when creating a plan. Although the inputs
Chris@10 1353 and outputs have different data distributions in general, it is
Chris@10 1354 guaranteed that the @emph{output} data distribution of an
Chris@10 1355 @code{FFTW_FORWARD} plan will match the @emph{input} data distribution
Chris@10 1356 of an @code{FFTW_BACKWARD} plan and vice versa; similarly for the
Chris@10 1357 @code{FFTW_MPI_SCRAMBLED_OUT} and @code{FFTW_MPI_SCRAMBLED_IN} flags.
Chris@10 1358 @xref{One-dimensional distributions}.
Chris@10 1359
Chris@10 1360 @node MPI Plan Creation, MPI Wisdom Communication, MPI Data Distribution Functions, FFTW MPI Reference
Chris@10 1361 @subsection MPI Plan Creation
Chris@10 1362
Chris@10 1363 @subsubheading Complex-data MPI DFTs
Chris@10 1364
Chris@10 1365 Plans for complex-data DFTs (@pxref{2d MPI example}) are created by:
Chris@10 1366
Chris@10 1367 @findex fftw_mpi_plan_dft_1d
Chris@10 1368 @findex fftw_mpi_plan_dft_2d
Chris@10 1369 @findex fftw_mpi_plan_dft_3d
Chris@10 1370 @findex fftw_mpi_plan_dft
Chris@10 1371 @findex fftw_mpi_plan_many_dft
Chris@10 1372 @example
Chris@10 1373 fftw_plan fftw_mpi_plan_dft_1d(ptrdiff_t n0, fftw_complex *in, fftw_complex *out,
Chris@10 1374 MPI_Comm comm, int sign, unsigned flags);
Chris@10 1375 fftw_plan fftw_mpi_plan_dft_2d(ptrdiff_t n0, ptrdiff_t n1,
Chris@10 1376 fftw_complex *in, fftw_complex *out,
Chris@10 1377 MPI_Comm comm, int sign, unsigned flags);
Chris@10 1378 fftw_plan fftw_mpi_plan_dft_3d(ptrdiff_t n0, ptrdiff_t n1, ptrdiff_t n2,
Chris@10 1379 fftw_complex *in, fftw_complex *out,
Chris@10 1380 MPI_Comm comm, int sign, unsigned flags);
Chris@10 1381 fftw_plan fftw_mpi_plan_dft(int rnk, const ptrdiff_t *n,
Chris@10 1382 fftw_complex *in, fftw_complex *out,
Chris@10 1383 MPI_Comm comm, int sign, unsigned flags);
Chris@10 1384 fftw_plan fftw_mpi_plan_many_dft(int rnk, const ptrdiff_t *n,
Chris@10 1385 ptrdiff_t howmany, ptrdiff_t block, ptrdiff_t tblock,
Chris@10 1386 fftw_complex *in, fftw_complex *out,
Chris@10 1387 MPI_Comm comm, int sign, unsigned flags);
Chris@10 1388 @end example
Chris@10 1389
Chris@10 1390 @cindex MPI communicator
Chris@10 1391 @cindex collective function
Chris@10 1392 These are similar to their serial counterparts (@pxref{Complex DFTs})
Chris@10 1393 in specifying the dimensions, sign, and flags of the transform. The
Chris@10 1394 @code{comm} argument gives an MPI communicator that specifies the set
Chris@10 1395 of processes to participate in the transform; plan creation is a
Chris@10 1396 collective function that must be called for all processes in the
Chris@10 1397 communicator. The @code{in} and @code{out} pointers refer only to a
Chris@10 1398 portion of the overall transform data (@pxref{MPI Data Distribution})
Chris@10 1399 as specified by the @samp{local_size} functions in the previous
Chris@10 1400 section. Unless @code{flags} contains @code{FFTW_ESTIMATE}, these
Chris@10 1401 arrays are overwritten during plan creation as for the serial
Chris@10 1402 interface. For multi-dimensional transforms, any dimensions @code{>
Chris@10 1403 1} are supported; for one-dimensional transforms, only composite
Chris@10 1404 (non-prime) @code{n0} are currently supported (unlike the serial
Chris@10 1405 FFTW). Requesting an unsupported transform size will yield a
Chris@10 1406 @code{NULL} plan. (As in the serial interface, highly composite sizes
Chris@10 1407 generally yield the best performance.)
Chris@10 1408
Chris@10 1409 @cindex advanced interface
Chris@10 1410 @ctindex FFTW_MPI_DEFAULT_BLOCK
Chris@10 1411 @cindex stride
Chris@10 1412 The advanced-interface @code{fftw_mpi_plan_many_dft} additionally
Chris@10 1413 allows you to specify the block sizes for the first dimension
Chris@10 1414 (@code{block}) of the @ndims{} input data and the first dimension
Chris@10 1415 (@code{tblock}) of the @ndimstrans{} transposed data (at intermediate
Chris@10 1416 steps of the transform, and for the output if
Chris@10 1417 @code{FFTW_TRANSPOSED_OUT} is specified in @code{flags}). These must
Chris@10 1418 be the same block sizes as were passed to the corresponding
Chris@10 1419 @samp{local_size} function; you can pass @code{FFTW_MPI_DEFAULT_BLOCK}
Chris@10 1420 to use FFTW's default block size as in the basic interface. Also, the
Chris@10 1421 @code{howmany} parameter specifies that the transform is of contiguous
Chris@10 1422 @code{howmany}-tuples rather than individual complex numbers; this
Chris@10 1423 corresponds to the same parameter in the serial advanced interface
Chris@10 1424 (@pxref{Advanced Complex DFTs}) with @code{stride = howmany} and
Chris@10 1425 @code{dist = 1}.
Chris@10 1426
Chris@10 1427 @subsubheading MPI flags
Chris@10 1428
Chris@10 1429 The @code{flags} can be any of those for the serial FFTW
Chris@10 1430 (@pxref{Planner Flags}), and in addition may include one or more of
Chris@10 1431 the following MPI-specific flags, which improve performance at the
Chris@10 1432 cost of changing the output or input data formats.
Chris@10 1433
Chris@10 1434 @itemize @bullet
Chris@10 1435
Chris@10 1436 @item
Chris@10 1437 @ctindex FFTW_MPI_SCRAMBLED_OUT
Chris@10 1438 @ctindex FFTW_MPI_SCRAMBLED_IN
Chris@10 1439 @code{FFTW_MPI_SCRAMBLED_OUT}, @code{FFTW_MPI_SCRAMBLED_IN}: valid for
Chris@10 1440 1d transforms only, these flags indicate that the output/input of the
Chris@10 1441 transform are in an undocumented ``scrambled'' order. A forward
Chris@10 1442 @code{FFTW_MPI_SCRAMBLED_OUT} transform can be inverted by a backward
Chris@10 1443 @code{FFTW_MPI_SCRAMBLED_IN} (times the usual 1/@i{N} normalization).
Chris@10 1444 @xref{One-dimensional distributions}.
Chris@10 1445
Chris@10 1446 @item
Chris@10 1447 @ctindex FFTW_MPI_TRANSPOSED_OUT
Chris@10 1448 @ctindex FFTW_MPI_TRANSPOSED_IN
Chris@10 1449 @code{FFTW_MPI_TRANSPOSED_OUT}, @code{FFTW_MPI_TRANSPOSED_IN}: valid
Chris@10 1450 for multidimensional (@code{rnk > 1}) transforms only, these flags
Chris@10 1451 specify that the output or input of an @ndims{} transform is
Chris@10 1452 transposed to @ndimstrans{}. @xref{Transposed distributions}.
Chris@10 1453
Chris@10 1454 @end itemize
Chris@10 1455
Chris@10 1456 @subsubheading Real-data MPI DFTs
Chris@10 1457
Chris@10 1458 @cindex r2c
Chris@10 1459 Plans for real-input/output (r2c/c2r) DFTs (@pxref{Multi-dimensional
Chris@10 1460 MPI DFTs of Real Data}) are created by:
Chris@10 1461
Chris@10 1462 @findex fftw_mpi_plan_dft_r2c_2d
Chris@10 1463 @findex fftw_mpi_plan_dft_r2c_2d
Chris@10 1464 @findex fftw_mpi_plan_dft_r2c_3d
Chris@10 1465 @findex fftw_mpi_plan_dft_r2c
Chris@10 1466 @findex fftw_mpi_plan_dft_c2r_2d
Chris@10 1467 @findex fftw_mpi_plan_dft_c2r_2d
Chris@10 1468 @findex fftw_mpi_plan_dft_c2r_3d
Chris@10 1469 @findex fftw_mpi_plan_dft_c2r
Chris@10 1470 @example
Chris@10 1471 fftw_plan fftw_mpi_plan_dft_r2c_2d(ptrdiff_t n0, ptrdiff_t n1,
Chris@10 1472 double *in, fftw_complex *out,
Chris@10 1473 MPI_Comm comm, unsigned flags);
Chris@10 1474 fftw_plan fftw_mpi_plan_dft_r2c_2d(ptrdiff_t n0, ptrdiff_t n1,
Chris@10 1475 double *in, fftw_complex *out,
Chris@10 1476 MPI_Comm comm, unsigned flags);
Chris@10 1477 fftw_plan fftw_mpi_plan_dft_r2c_3d(ptrdiff_t n0, ptrdiff_t n1, ptrdiff_t n2,
Chris@10 1478 double *in, fftw_complex *out,
Chris@10 1479 MPI_Comm comm, unsigned flags);
Chris@10 1480 fftw_plan fftw_mpi_plan_dft_r2c(int rnk, const ptrdiff_t *n,
Chris@10 1481 double *in, fftw_complex *out,
Chris@10 1482 MPI_Comm comm, unsigned flags);
Chris@10 1483 fftw_plan fftw_mpi_plan_dft_c2r_2d(ptrdiff_t n0, ptrdiff_t n1,
Chris@10 1484 fftw_complex *in, double *out,
Chris@10 1485 MPI_Comm comm, unsigned flags);
Chris@10 1486 fftw_plan fftw_mpi_plan_dft_c2r_2d(ptrdiff_t n0, ptrdiff_t n1,
Chris@10 1487 fftw_complex *in, double *out,
Chris@10 1488 MPI_Comm comm, unsigned flags);
Chris@10 1489 fftw_plan fftw_mpi_plan_dft_c2r_3d(ptrdiff_t n0, ptrdiff_t n1, ptrdiff_t n2,
Chris@10 1490 fftw_complex *in, double *out,
Chris@10 1491 MPI_Comm comm, unsigned flags);
Chris@10 1492 fftw_plan fftw_mpi_plan_dft_c2r(int rnk, const ptrdiff_t *n,
Chris@10 1493 fftw_complex *in, double *out,
Chris@10 1494 MPI_Comm comm, unsigned flags);
Chris@10 1495 @end example
Chris@10 1496
Chris@10 1497 Similar to the serial interface (@pxref{Real-data DFTs}), these
Chris@10 1498 transform logically @ndims{} real data to/from @ndimshalf{} complex
Chris@10 1499 data, representing the non-redundant half of the conjugate-symmetry
Chris@10 1500 output of a real-input DFT (@pxref{Multi-dimensional Transforms}).
Chris@10 1501 However, the real array must be stored within a padded @ndimspad{}
Chris@10 1502 array (much like the in-place serial r2c transforms, but here for
Chris@10 1503 out-of-place transforms as well). Currently, only multi-dimensional
Chris@10 1504 (@code{rnk > 1}) r2c/c2r transforms are supported (requesting a plan
Chris@10 1505 for @code{rnk = 1} will yield @code{NULL}). As explained above
Chris@10 1506 (@pxref{Multi-dimensional MPI DFTs of Real Data}), the data
Chris@10 1507 distribution of both the real and complex arrays is given by the
Chris@10 1508 @samp{local_size} function called for the dimensions of the
Chris@10 1509 @emph{complex} array. Similar to the other planning functions, the
Chris@10 1510 input and output arrays are overwritten when the plan is created
Chris@10 1511 except in @code{FFTW_ESTIMATE} mode.
Chris@10 1512
Chris@10 1513 As for the complex DFTs above, there is an advance interface that
Chris@10 1514 allows you to manually specify block sizes and to transform contiguous
Chris@10 1515 @code{howmany}-tuples of real/complex numbers:
Chris@10 1516
Chris@10 1517 @findex fftw_mpi_plan_many_dft_r2c
Chris@10 1518 @findex fftw_mpi_plan_many_dft_c2r
Chris@10 1519 @example
Chris@10 1520 fftw_plan fftw_mpi_plan_many_dft_r2c
Chris@10 1521 (int rnk, const ptrdiff_t *n, ptrdiff_t howmany,
Chris@10 1522 ptrdiff_t iblock, ptrdiff_t oblock,
Chris@10 1523 double *in, fftw_complex *out,
Chris@10 1524 MPI_Comm comm, unsigned flags);
Chris@10 1525 fftw_plan fftw_mpi_plan_many_dft_c2r
Chris@10 1526 (int rnk, const ptrdiff_t *n, ptrdiff_t howmany,
Chris@10 1527 ptrdiff_t iblock, ptrdiff_t oblock,
Chris@10 1528 fftw_complex *in, double *out,
Chris@10 1529 MPI_Comm comm, unsigned flags);
Chris@10 1530 @end example
Chris@10 1531
Chris@10 1532 @subsubheading MPI r2r transforms
Chris@10 1533
Chris@10 1534 @cindex r2r
Chris@10 1535 There are corresponding plan-creation routines for r2r
Chris@10 1536 transforms (@pxref{More DFTs of Real Data}), currently supporting
Chris@10 1537 multidimensional (@code{rnk > 1}) transforms only (@code{rnk = 1} will
Chris@10 1538 yield a @code{NULL} plan):
Chris@10 1539
Chris@10 1540 @example
Chris@10 1541 fftw_plan fftw_mpi_plan_r2r_2d(ptrdiff_t n0, ptrdiff_t n1,
Chris@10 1542 double *in, double *out,
Chris@10 1543 MPI_Comm comm,
Chris@10 1544 fftw_r2r_kind kind0, fftw_r2r_kind kind1,
Chris@10 1545 unsigned flags);
Chris@10 1546 fftw_plan fftw_mpi_plan_r2r_3d(ptrdiff_t n0, ptrdiff_t n1, ptrdiff_t n2,
Chris@10 1547 double *in, double *out,
Chris@10 1548 MPI_Comm comm,
Chris@10 1549 fftw_r2r_kind kind0, fftw_r2r_kind kind1, fftw_r2r_kind kind2,
Chris@10 1550 unsigned flags);
Chris@10 1551 fftw_plan fftw_mpi_plan_r2r(int rnk, const ptrdiff_t *n,
Chris@10 1552 double *in, double *out,
Chris@10 1553 MPI_Comm comm, const fftw_r2r_kind *kind,
Chris@10 1554 unsigned flags);
Chris@10 1555 fftw_plan fftw_mpi_plan_many_r2r(int rnk, const ptrdiff_t *n,
Chris@10 1556 ptrdiff_t iblock, ptrdiff_t oblock,
Chris@10 1557 double *in, double *out,
Chris@10 1558 MPI_Comm comm, const fftw_r2r_kind *kind,
Chris@10 1559 unsigned flags);
Chris@10 1560 @end example
Chris@10 1561
Chris@10 1562 The parameters are much the same as for the complex DFTs above, except
Chris@10 1563 that the arrays are of real numbers (and hence the outputs of the
Chris@10 1564 @samp{local_size} data-distribution functions should be interpreted as
Chris@10 1565 counts of real rather than complex numbers). Also, the @code{kind}
Chris@10 1566 parameters specify the r2r kinds along each dimension as for the
Chris@10 1567 serial interface (@pxref{Real-to-Real Transform Kinds}). @xref{Other
Chris@10 1568 Multi-dimensional Real-data MPI Transforms}.
Chris@10 1569
Chris@10 1570 @subsubheading MPI transposition
Chris@10 1571 @cindex transpose
Chris@10 1572
Chris@10 1573 FFTW also provides routines to plan a transpose of a distributed
Chris@10 1574 @code{n0} by @code{n1} array of real numbers, or an array of
Chris@10 1575 @code{howmany}-tuples of real numbers with specified block sizes
Chris@10 1576 (@pxref{FFTW MPI Transposes}):
Chris@10 1577
Chris@10 1578 @findex fftw_mpi_plan_transpose
Chris@10 1579 @findex fftw_mpi_plan_many_transpose
Chris@10 1580 @example
Chris@10 1581 fftw_plan fftw_mpi_plan_transpose(ptrdiff_t n0, ptrdiff_t n1,
Chris@10 1582 double *in, double *out,
Chris@10 1583 MPI_Comm comm, unsigned flags);
Chris@10 1584 fftw_plan fftw_mpi_plan_many_transpose
Chris@10 1585 (ptrdiff_t n0, ptrdiff_t n1, ptrdiff_t howmany,
Chris@10 1586 ptrdiff_t block0, ptrdiff_t block1,
Chris@10 1587 double *in, double *out, MPI_Comm comm, unsigned flags);
Chris@10 1588 @end example
Chris@10 1589
Chris@10 1590 @cindex new-array execution
Chris@10 1591 @findex fftw_mpi_execute_r2r
Chris@10 1592 These plans are used with the @code{fftw_mpi_execute_r2r} new-array
Chris@10 1593 execute function (@pxref{Using MPI Plans }), since they count as (rank
Chris@10 1594 zero) r2r plans from FFTW's perspective.
Chris@10 1595
Chris@10 1596 @node MPI Wisdom Communication, , MPI Plan Creation, FFTW MPI Reference
Chris@10 1597 @subsection MPI Wisdom Communication
Chris@10 1598
Chris@10 1599 To facilitate synchronizing wisdom among the different MPI processes,
Chris@10 1600 we provide two functions:
Chris@10 1601
Chris@10 1602 @findex fftw_mpi_gather_wisdom
Chris@10 1603 @findex fftw_mpi_broadcast_wisdom
Chris@10 1604 @example
Chris@10 1605 void fftw_mpi_gather_wisdom(MPI_Comm comm);
Chris@10 1606 void fftw_mpi_broadcast_wisdom(MPI_Comm comm);
Chris@10 1607 @end example
Chris@10 1608
Chris@10 1609 The @code{fftw_mpi_gather_wisdom} function gathers all wisdom in the
Chris@10 1610 given communicator @code{comm} to the process of rank 0 in the
Chris@10 1611 communicator: that process obtains the union of all wisdom on all the
Chris@10 1612 processes. As a side effect, some other processes will gain
Chris@10 1613 additional wisdom from other processes, but only process 0 will gain
Chris@10 1614 the complete union.
Chris@10 1615
Chris@10 1616 The @code{fftw_mpi_broadcast_wisdom} does the reverse: it exports
Chris@10 1617 wisdom from process 0 in @code{comm} to all other processes in the
Chris@10 1618 communicator, replacing any wisdom they currently have.
Chris@10 1619
Chris@10 1620 @xref{FFTW MPI Wisdom}.
Chris@10 1621
Chris@10 1622 @c ------------------------------------------------------------
Chris@10 1623 @node FFTW MPI Fortran Interface, , FFTW MPI Reference, Distributed-memory FFTW with MPI
Chris@10 1624 @section FFTW MPI Fortran Interface
Chris@10 1625 @cindex Fortran interface
Chris@10 1626
Chris@10 1627 @cindex iso_c_binding
Chris@10 1628 The FFTW MPI interface is callable from modern Fortran compilers
Chris@10 1629 supporting the Fortran 2003 @code{iso_c_binding} standard for calling
Chris@10 1630 C functions. As described in @ref{Calling FFTW from Modern Fortran},
Chris@10 1631 this means that you can directly call FFTW's C interface from Fortran
Chris@10 1632 with only minor changes in syntax. There are, however, a few things
Chris@10 1633 specific to the MPI interface to keep in mind:
Chris@10 1634
Chris@10 1635 @itemize @bullet
Chris@10 1636
Chris@10 1637 @item
Chris@10 1638 Instead of including @code{fftw3.f03} as in @ref{Overview of Fortran
Chris@10 1639 interface }, you should @code{include 'fftw3-mpi.f03'} (after
Chris@10 1640 @code{use, intrinsic :: iso_c_binding} as before). The
Chris@10 1641 @code{fftw3-mpi.f03} file includes @code{fftw3.f03}, so you should
Chris@10 1642 @emph{not} @code{include} them both yourself. (You will also want to
Chris@10 1643 include the MPI header file, usually via @code{include 'mpif.h'} or
Chris@10 1644 similar, although though this is not needed by @code{fftw3-mpi.f03}
Chris@10 1645 @i{per se}.) (To use the @samp{fftwl_} @code{long double} extended-precision routines in supporting compilers, you should include @code{fftw3f-mpi.f03} in @emph{addition} to @code{fftw3-mpi.f03}. @xref{Extended and quadruple precision in Fortran}.)
Chris@10 1646
Chris@10 1647 @item
Chris@10 1648 Because of the different storage conventions between C and Fortran,
Chris@10 1649 you reverse the order of your array dimensions when passing them to
Chris@10 1650 FFTW (@pxref{Reversing array dimensions}). This is merely a
Chris@10 1651 difference in notation and incurs no performance overhead. However,
Chris@10 1652 it means that, whereas in C the @emph{first} dimension is distributed,
Chris@10 1653 in Fortran the @emph{last} dimension of your array is distributed.
Chris@10 1654
Chris@10 1655 @item
Chris@10 1656 @cindex MPI communicator
Chris@10 1657 In Fortran, communicators are stored as @code{integer} types; there is
Chris@10 1658 no @code{MPI_Comm} type, nor is there any way to access a C
Chris@10 1659 @code{MPI_Comm}. Fortunately, this is taken care of for you by the
Chris@10 1660 FFTW Fortran interface: whenever the C interface expects an
Chris@10 1661 @code{MPI_Comm} type, you should pass the Fortran communicator as an
Chris@10 1662 @code{integer}.@footnote{Technically, this is because you aren't
Chris@10 1663 actually calling the C functions directly. You are calling wrapper
Chris@10 1664 functions that translate the communicator with @code{MPI_Comm_f2c}
Chris@10 1665 before calling the ordinary C interface. This is all done
Chris@10 1666 transparently, however, since the @code{fftw3-mpi.f03} interface file
Chris@10 1667 renames the wrappers so that they are called in Fortran with the same
Chris@10 1668 names as the C interface functions.}
Chris@10 1669
Chris@10 1670 @item
Chris@10 1671 Because you need to call the @samp{local_size} function to find out
Chris@10 1672 how much space to allocate, and this may be @emph{larger} than the
Chris@10 1673 local portion of the array (@pxref{MPI Data Distribution}), you should
Chris@10 1674 @emph{always} allocate your arrays dynamically using FFTW's allocation
Chris@10 1675 routines as described in @ref{Allocating aligned memory in Fortran}.
Chris@10 1676 (Coincidentally, this also provides the best performance by
Chris@10 1677 guaranteeding proper data alignment.)
Chris@10 1678
Chris@10 1679 @item
Chris@10 1680 Because all sizes in the MPI FFTW interface are declared as
Chris@10 1681 @code{ptrdiff_t} in C, you should use @code{integer(C_INTPTR_T)} in
Chris@10 1682 Fortran (@pxref{FFTW Fortran type reference}).
Chris@10 1683
Chris@10 1684 @item
Chris@10 1685 @findex fftw_execute_dft
Chris@10 1686 @findex fftw_mpi_execute_dft
Chris@10 1687 @cindex new-array execution
Chris@10 1688 In Fortran, because of the language semantics, we generally recommend
Chris@10 1689 using the new-array execute functions for all plans, even in the
Chris@10 1690 common case where you are executing the plan on the same arrays for
Chris@10 1691 which the plan was created (@pxref{Plan execution in Fortran}).
Chris@10 1692 However, note that in the MPI interface these functions are changed:
Chris@10 1693 @code{fftw_execute_dft} becomes @code{fftw_mpi_execute_dft},
Chris@10 1694 etcetera. @xref{Using MPI Plans}.
Chris@10 1695
Chris@10 1696 @end itemize
Chris@10 1697
Chris@10 1698 For example, here is a Fortran code snippet to perform a distributed
Chris@10 1699 @twodims{L,M} complex DFT in-place. (This assumes you have already
Chris@10 1700 initialized MPI with @code{MPI_init} and have also performed
Chris@10 1701 @code{call fftw_mpi_init}.)
Chris@10 1702
Chris@10 1703 @example
Chris@10 1704 use, intrinsic :: iso_c_binding
Chris@10 1705 include 'fftw3-mpi.f03'
Chris@10 1706 integer(C_INTPTR_T), parameter :: L = ...
Chris@10 1707 integer(C_INTPTR_T), parameter :: M = ...
Chris@10 1708 type(C_PTR) :: plan, cdata
Chris@10 1709 complex(C_DOUBLE_COMPLEX), pointer :: data(:,:)
Chris@10 1710 integer(C_INTPTR_T) :: i, j, alloc_local, local_M, local_j_offset
Chris@10 1711
Chris@10 1712 ! @r{get local data size and allocate (note dimension reversal)}
Chris@10 1713 alloc_local = fftw_mpi_local_size_2d(M, L, MPI_COMM_WORLD, &
Chris@10 1714 local_M, local_j_offset)
Chris@10 1715 cdata = fftw_alloc_complex(alloc_local)
Chris@10 1716 call c_f_pointer(cdata, data, [L,local_M])
Chris@10 1717
Chris@10 1718 ! @r{create MPI plan for in-place forward DFT (note dimension reversal)}
Chris@10 1719 plan = fftw_mpi_plan_dft_2d(M, L, data, data, MPI_COMM_WORLD, &
Chris@10 1720 FFTW_FORWARD, FFTW_MEASURE)
Chris@10 1721
Chris@10 1722 ! @r{initialize data to some function} my_function(i,j)
Chris@10 1723 do j = 1, local_M
Chris@10 1724 do i = 1, L
Chris@10 1725 data(i, j) = my_function(i, j + local_j_offset)
Chris@10 1726 end do
Chris@10 1727 end do
Chris@10 1728
Chris@10 1729 ! @r{compute transform (as many times as desired)}
Chris@10 1730 call fftw_mpi_execute_dft(plan, data, data)
Chris@10 1731
Chris@10 1732 call fftw_destroy_plan(plan)
Chris@10 1733 call fftw_free(cdata)
Chris@10 1734 @end example
Chris@10 1735
Chris@10 1736 Note that when we called @code{fftw_mpi_local_size_2d} and
Chris@10 1737 @code{fftw_mpi_plan_dft_2d} with the dimensions in reversed order,
Chris@10 1738 since a @twodims{L,M} Fortran array is viewed by FFTW in C as a
Chris@10 1739 @twodims{M, L} array. This means that the array was distributed over
Chris@10 1740 the @code{M} dimension, the local portion of which is a
Chris@10 1741 @twodims{L,local_M} array in Fortran. (You must @emph{not} use an
Chris@10 1742 @code{allocate} statement to allocate an @twodims{L,local_M} array,
Chris@10 1743 however; you must allocate @code{alloc_local} complex numbers, which
Chris@10 1744 may be greater than @code{L * local_M}, in order to reserve space for
Chris@10 1745 intermediate steps of the transform.) Finally, we mention that
Chris@10 1746 because C's array indices are zero-based, the @code{local_j_offset}
Chris@10 1747 argument can conveniently be interpreted as an offset in the 1-based
Chris@10 1748 @code{j} index (rather than as a starting index as in C).
Chris@10 1749
Chris@10 1750 If instead you had used the @code{ior(FFTW_MEASURE,
Chris@10 1751 FFTW_MPI_TRANSPOSED_OUT)} flag, the output of the transform would be a
Chris@10 1752 transposed @twodims{M,local_L} array, associated with the @emph{same}
Chris@10 1753 @code{cdata} allocation (since the transform is in-place), and which
Chris@10 1754 you could declare with:
Chris@10 1755
Chris@10 1756 @example
Chris@10 1757 complex(C_DOUBLE_COMPLEX), pointer :: tdata(:,:)
Chris@10 1758 ...
Chris@10 1759 call c_f_pointer(cdata, tdata, [M,local_L])
Chris@10 1760 @end example
Chris@10 1761
Chris@10 1762 where @code{local_L} would have been obtained by changing the
Chris@10 1763 @code{fftw_mpi_local_size_2d} call to:
Chris@10 1764
Chris@10 1765 @example
Chris@10 1766 alloc_local = fftw_mpi_local_size_2d_transposed(M, L, MPI_COMM_WORLD, &
Chris@10 1767 local_M, local_j_offset, local_L, local_i_offset)
Chris@10 1768 @end example