Chris@82: Chris@82: Chris@82: Chris@82: Chris@82: Chris@82: FFTW 3.3.8: Reversing array dimensions Chris@82: Chris@82: Chris@82: Chris@82: Chris@82: Chris@82: Chris@82: Chris@82: Chris@82: Chris@82: Chris@82: Chris@82: Chris@82: Chris@82: Chris@82: Chris@82: Chris@82: Chris@82: Chris@82: Chris@82: Chris@82:
Chris@82:

Chris@82: Next: , Previous: , Up: Calling FFTW from Modern Fortran   [Contents][Index]

Chris@82:
Chris@82:
Chris@82: Chris@82:

7.2 Reversing array dimensions

Chris@82: Chris@82: Chris@82: Chris@82:

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

Chris@82:

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

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

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

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

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

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

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

Chris@82: Chris@82:

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

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

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

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

Chris@82: Next: , Previous: , Up: Calling FFTW from Modern Fortran   [Contents][Index]

Chris@82:
Chris@82: Chris@82: Chris@82: Chris@82: Chris@82: