Chris@42: Chris@42: Chris@42: Chris@42: Chris@42:
Chris@42:Chris@42: Previous: FFTW MPI Reference, Up: Distributed-memory FFTW with MPI [Contents][Index]
Chris@42:The FFTW MPI interface is callable from modern Fortran compilers
Chris@42: supporting the Fortran 2003 iso_c_binding
standard for calling
Chris@42: C functions. As described in Calling FFTW from Modern Fortran,
Chris@42: this means that you can directly call FFTW’s C interface from Fortran
Chris@42: with only minor changes in syntax. There are, however, a few things
Chris@42: specific to the MPI interface to keep in mind:
Chris@42:
fftw3.f03
as in Overview of Fortran interface, you should include 'fftw3-mpi.f03'
(after
Chris@42: use, intrinsic :: iso_c_binding
as before). The
Chris@42: fftw3-mpi.f03
file includes fftw3.f03
, so you should
Chris@42: not include
them both yourself. (You will also want to
Chris@42: include the MPI header file, usually via include 'mpif.h'
or
Chris@42: similar, although though this is not needed by fftw3-mpi.f03
Chris@42: per se.) (To use the ‘fftwl_’ long double
extended-precision routines in supporting compilers, you should include fftw3f-mpi.f03
in addition to fftw3-mpi.f03
. See Extended and quadruple precision in Fortran.)
Chris@42:
Chris@42: integer
types; there is
Chris@42: no MPI_Comm
type, nor is there any way to access a C
Chris@42: MPI_Comm
. Fortunately, this is taken care of for you by the
Chris@42: FFTW Fortran interface: whenever the C interface expects an
Chris@42: MPI_Comm
type, you should pass the Fortran communicator as an
Chris@42: integer
.8
Chris@42:
Chris@42: ptrdiff_t
in C, you should use integer(C_INTPTR_T)
in
Chris@42: Fortran (see FFTW Fortran type reference).
Chris@42:
Chris@42: fftw_execute_dft
becomes fftw_mpi_execute_dft
,
Chris@42: etcetera. See Using MPI Plans.
Chris@42:
Chris@42: For example, here is a Fortran code snippet to perform a distributed
Chris@42: L × M complex DFT in-place. (This assumes you have already
Chris@42: initialized MPI with MPI_init
and have also performed
Chris@42: call fftw_mpi_init
.)
Chris@42:
use, intrinsic :: iso_c_binding Chris@42: include 'fftw3-mpi.f03' Chris@42: integer(C_INTPTR_T), parameter :: L = ... Chris@42: integer(C_INTPTR_T), parameter :: M = ... Chris@42: type(C_PTR) :: plan, cdata Chris@42: complex(C_DOUBLE_COMPLEX), pointer :: data(:,:) Chris@42: integer(C_INTPTR_T) :: i, j, alloc_local, local_M, local_j_offset Chris@42: Chris@42: ! get local data size and allocate (note dimension reversal) Chris@42: alloc_local = fftw_mpi_local_size_2d(M, L, MPI_COMM_WORLD, & Chris@42: local_M, local_j_offset) Chris@42: cdata = fftw_alloc_complex(alloc_local) Chris@42: call c_f_pointer(cdata, data, [L,local_M]) Chris@42: Chris@42: ! create MPI plan for in-place forward DFT (note dimension reversal) Chris@42: plan = fftw_mpi_plan_dft_2d(M, L, data, data, MPI_COMM_WORLD, & Chris@42: FFTW_FORWARD, FFTW_MEASURE) Chris@42: Chris@42: ! initialize data to some function my_function(i,j) Chris@42: do j = 1, local_M Chris@42: do i = 1, L Chris@42: data(i, j) = my_function(i, j + local_j_offset) Chris@42: end do Chris@42: end do Chris@42: Chris@42: ! compute transform (as many times as desired) Chris@42: call fftw_mpi_execute_dft(plan, data, data) Chris@42: Chris@42: call fftw_destroy_plan(plan) Chris@42: call fftw_free(cdata) Chris@42:
Note that when we called fftw_mpi_local_size_2d
and
Chris@42: fftw_mpi_plan_dft_2d
with the dimensions in reversed order,
Chris@42: since a L × M Fortran array is viewed by FFTW in C as a
Chris@42: M × L array. This means that the array was distributed over
Chris@42: the M
dimension, the local portion of which is a
Chris@42: L × local_M array in Fortran. (You must not use an
Chris@42: allocate
statement to allocate an L × local_M array,
Chris@42: however; you must allocate alloc_local
complex numbers, which
Chris@42: may be greater than L * local_M
, in order to reserve space for
Chris@42: intermediate steps of the transform.) Finally, we mention that
Chris@42: because C’s array indices are zero-based, the local_j_offset
Chris@42: argument can conveniently be interpreted as an offset in the 1-based
Chris@42: j
index (rather than as a starting index as in C).
Chris@42:
If instead you had used the ior(FFTW_MEASURE,
Chris@42: FFTW_MPI_TRANSPOSED_OUT)
flag, the output of the transform would be a
Chris@42: transposed M × local_L array, associated with the same
Chris@42: cdata
allocation (since the transform is in-place), and which
Chris@42: you could declare with:
Chris@42:
complex(C_DOUBLE_COMPLEX), pointer :: tdata(:,:) Chris@42: ... Chris@42: call c_f_pointer(cdata, tdata, [M,local_L]) Chris@42:
where local_L
would have been obtained by changing the
Chris@42: fftw_mpi_local_size_2d
call to:
Chris@42:
alloc_local = fftw_mpi_local_size_2d_transposed(M, L, MPI_COMM_WORLD, & Chris@42: local_M, local_j_offset, local_L, local_i_offset) Chris@42:
Technically, this is because you aren’t
Chris@42: actually calling the C functions directly. You are calling wrapper
Chris@42: functions that translate the communicator with MPI_Comm_f2c
Chris@42: before calling the ordinary C interface. This is all done
Chris@42: transparently, however, since the fftw3-mpi.f03
interface file
Chris@42: renames the wrappers so that they are called in Fortran with the same
Chris@42: names as the C interface functions.
Chris@42: Previous: FFTW MPI Reference, Up: Distributed-memory FFTW with MPI [Contents][Index]
Chris@42: