Chris@19: Chris@19: Chris@19: Reversing array dimensions - FFTW 3.3.4 Chris@19: Chris@19: Chris@19: Chris@19: Chris@19: Chris@19: Chris@19: Chris@19: Chris@19: Chris@19: Chris@19: Chris@19: Chris@19: Chris@19:
Chris@19: Chris@19:

Chris@19: Next: , Chris@19: Previous: Overview of Fortran interface, Chris@19: Up: Calling FFTW from Modern Fortran Chris@19:


Chris@19:
Chris@19: Chris@19:

7.2 Reversing array dimensions

Chris@19: Chris@19:

A minor annoyance in calling FFTW from Fortran is that FFTW's array Chris@19: dimensions are defined in the C convention (row-major order), while Chris@19: Fortran's array dimensions are the opposite convention (column-major Chris@19: order). See Multi-dimensional Array Format. This is just a Chris@19: bookkeeping difference, with no effect on performance. The only Chris@19: consequence of this is that, whenever you create an FFTW plan for a Chris@19: multi-dimensional transform, you must always reverse the Chris@19: ordering of the dimensions. Chris@19: Chris@19:

For example, consider the three-dimensional (L × M × N) arrays: Chris@19: Chris@19:

       complex(C_DOUBLE_COMPLEX), dimension(L,M,N) :: in, out
Chris@19: 
Chris@19:

To plan a DFT for these arrays using fftw_plan_dft_3d, you could do: Chris@19: Chris@19:

Chris@19:

       plan = fftw_plan_dft_3d(N,M,L, in,out, FFTW_FORWARD,FFTW_ESTIMATE)
Chris@19: 
Chris@19:

That is, from FFTW's perspective this is a N × M × L array. Chris@19: No data transposition need occur, as this is only Chris@19: notation. Similarly, to use the more generic routine Chris@19: fftw_plan_dft with the same arrays, you could do: Chris@19: Chris@19:

       integer(C_INT), dimension(3) :: n = [N,M,L]
Chris@19:        plan = fftw_plan_dft_3d(3, n, in,out, FFTW_FORWARD,FFTW_ESTIMATE)
Chris@19: 
Chris@19:

Note, by the way, that this is different from the legacy Fortran Chris@19: interface (see Fortran-interface routines), which automatically Chris@19: reverses the order of the array dimension for you. Here, you are Chris@19: calling the C interface directly, so there is no “translation” layer. Chris@19: Chris@19:

An important thing to keep in mind is the implication of this for Chris@19: multidimensional real-to-complex transforms (see Multi-Dimensional DFTs of Real Data). In C, a multidimensional real-to-complex DFT Chris@19: chops the last dimension roughly in half (N × M × L real input Chris@19: goes to N × M × L/2+1 complex output). In Fortran, because Chris@19: the array dimension notation is reversed, the first dimension of Chris@19: the complex data is chopped roughly in half. For example consider the Chris@19: ‘r2c’ transform of L × M × N real input in Fortran: Chris@19: Chris@19:

Chris@19:

       type(C_PTR) :: plan
Chris@19:        real(C_DOUBLE), dimension(L,M,N) :: in
Chris@19:        complex(C_DOUBLE_COMPLEX), dimension(L/2+1,M,N) :: out
Chris@19:        plan = fftw_plan_dft_r2c_3d(N,M,L, in,out, FFTW_ESTIMATE)
Chris@19:        ...
Chris@19:        call fftw_execute_dft_r2c(plan, in, out)
Chris@19: 
Chris@19:

Alternatively, for an in-place r2c transform, as described in the C Chris@19: documentation we must pad the first dimension of the Chris@19: real input with an extra two entries (which are ignored by FFTW) so as Chris@19: to leave enough space for the complex output. The input is Chris@19: allocated as a 2[L/2+1] × M × N array, even though only Chris@19: L × M × N of it is actually used. In this example, we will Chris@19: allocate the array as a pointer type, using ‘fftw_alloc’ to Chris@19: ensure aligned memory for maximum performance (see Allocating aligned memory in Fortran); this also makes it easy to reference the Chris@19: same memory as both a real array and a complex array. Chris@19: Chris@19:

Chris@19:

       real(C_DOUBLE), pointer :: in(:,:,:)
Chris@19:        complex(C_DOUBLE_COMPLEX), pointer :: out(:,:,:)
Chris@19:        type(C_PTR) :: plan, data
Chris@19:        data = fftw_alloc_complex(int((L/2+1) * M * N, C_SIZE_T))
Chris@19:        call c_f_pointer(data, in, [2*(L/2+1),M,N])
Chris@19:        call c_f_pointer(data, out, [L/2+1,M,N])
Chris@19:        plan = fftw_plan_dft_r2c_3d(N,M,L, in,out, FFTW_ESTIMATE)
Chris@19:        ...
Chris@19:        call fftw_execute_dft_r2c(plan, in, out)
Chris@19:        ...
Chris@19:        call fftw_destroy_plan(plan)
Chris@19:        call fftw_free(data)
Chris@19: 
Chris@19: Chris@19: Chris@19: