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