Whole document tree
    

Whole document tree

FFTW - Calling FFTW from Fortran Go to the first, previous, next, last section, table of contents.


Calling FFTW from Fortran

The standard FFTW libraries include special wrapper functions that allow Fortran programs to call FFTW subroutines. This chapter describes how those functions may be employed to use FFTW from Fortran. We assume here that the reader is already familiar with the usage of FFTW in C, as described elsewhere in this manual.

In general, it is not possible to call C functions directly from Fortran, due to Fortran's inability to pass arguments by value and also because Fortran compilers typically expect identifiers to be mangled somehow for linking. However, if C functions are written in a special way, they are callable from Fortran, and we have employed this technique to create Fortran-callable "wrapper" functions around the main FFTW routines. These wrapper functions are included in the FFTW libraries by default, unless a Fortran compiler isn't found on your system or --disable-fortran is included in the configure flags.

As a result, calling FFTW from Fortran requires little more than appending `_f77' to the function names and then linking normally with the FFTW libraries. There are a few wrinkles, however, as we shall discuss below.

Wrapper Routines

All of the uniprocessor and multi-threaded transform routines have Fortran-callable wrappers, except for the wisdom import/export functions (since it is not possible to exchange string and file arguments portably with Fortran) and the specific planner routines (see Section Discussion on Specific Plans). The name of the wrapper routine is the same as that of the corresponding C routine, but with fftw/fftwnd/rfftw/rfftwnd replaced by fftw_f77/fftwnd_f77/rfftw_f77/rfftwnd_f77. For example, in Fortran, instead of calling fftw_one you would call fftw_f77_one.(8) For the most part, all of the arguments to the functions are the same, with the following exceptions:

  • plan variables (what would be of type fftw_plan, rfftwnd_plan, etcetera, in C), must be declared as a type that is the same size as a pointer (address) on your machine. (Fortran has no generic pointer type.) The Fortran integer type is usually the same size as a pointer, but you need to be wary (especially on 64-bit machines). (You could also use integer*4 on a 32-bit machine and integer*8 on a 64-bit machine.) Ugh. (g77 has a special type, integer(kind=7), that is defined to be the same size as a pointer.)
  • Any function that returns a value (e.g. fftw_create_plan) is converted into a subroutine. The return value is converted into an additional (first) parameter of the wrapper subroutine. (The reason for this is that some Fortran implementations seem to have trouble with C function return values.)
  • When performing one-dimensional FFTW_IN_PLACE transforms, you don't have the option of passing NULL for the out argument (since there is no way to pass NULL from Fortran). Therefore, when performing such transforms, you must allocate and pass a contiguous scratch array of the same size as the transform. Note that for in-place multi-dimensional ((r)fftwnd) transforms, the out argument is ignored, so you can pass anything for that parameter.
  • The wrapper routines expect multi-dimensional arrays to be in column-major order, which is the ordinary format of Fortran arrays. They do this transparently and costlessly simply by reversing the order of the dimensions passed to FFTW, but this has one important consequence for multi-dimensional real-complex transforms, discussed below.

In general, you should take care to use Fortran data types that correspond to (i.e. are the same size as) the C types used by FFTW. If your C and Fortran compilers are made by the same vendor, the correspondence is usually straightforward (i.e. integer corresponds to int, real corresponds to float, etcetera). Such simple correspondences are assumed in the examples below. The examples also assume that FFTW was compiled in double precision (the default).

FFTW Constants in Fortran

When creating plans in FFTW, a number of constants are used to specify options, such as FFTW_FORWARD or FFTW_USE_WISDOM. The same constants must be used with the wrapper routines, but of course the C header files where the constants are defined can't be incorporated directly into Fortran code.

Instead, we have placed Fortran equivalents of the FFTW constant definitions in the file fortran/fftw_f77.i of the FFTW package. If your Fortran compiler supports a preprocessor, you can use that to incorporate this file into your code whenever you need to call FFTW. Otherwise, you will have to paste the constant definitions in directly. They are:

      integer FFTW_FORWARD,FFTW_BACKWARD
      parameter (FFTW_FORWARD=-1,FFTW_BACKWARD=1)

      integer FFTW_REAL_TO_COMPLEX,FFTW_COMPLEX_TO_REAL
      parameter (FFTW_REAL_TO_COMPLEX=-1,FFTW_COMPLEX_TO_REAL=1)

      integer FFTW_ESTIMATE,FFTW_MEASURE
      parameter (FFTW_ESTIMATE=0,FFTW_MEASURE=1)

      integer FFTW_OUT_OF_PLACE,FFTW_IN_PLACE,FFTW_USE_WISDOM
      parameter (FFTW_OUT_OF_PLACE=0)
      parameter (FFTW_IN_PLACE=8,FFTW_USE_WISDOM=16)

      integer FFTW_THREADSAFE
      parameter (FFTW_THREADSAFE=128)

In C, you combine different flags (like FFTW_USE_WISDOM and FFTW_MEASURE) using the `|' operator; in Fortran you should just use `+'.

Fortran Examples

In C you might have something like the following to transform a one-dimensional complex array:

        fftw_complex in[N], *out[N];
        fftw_plan plan;

        plan = fftw_create_plan(N,FFTW_FORWARD,FFTW_ESTIMATE);
        fftw_one(plan,in,out);
        fftw_destroy_plan(plan);

In Fortran, you use the following to accomplish the same thing:

        double complex in, out
        dimension in(N), out(N)
        integer plan

        call fftw_f77_create_plan(plan,N,FFTW_FORWARD,FFTW_ESTIMATE)
        call fftw_f77_one(plan,in,out)
        call fftw_f77_destroy_plan(plan)

Notice how all routines are called as Fortran subroutines, and the plan is returned via the first argument to fftw_f77_create_plan. Important: these examples assume that integer is the same size as a pointer, and may need modification on a 64-bit machine. See Section Wrapper Routines, above. To do the same thing, but using 8 threads in parallel (see Section Multi-threaded FFTW), you would simply replace the call to fftw_f77_one with:

        call fftw_f77_threads_one(8,plan,in,out)

To transform a three-dimensional array in-place with C, you might do:

        fftw_complex arr[L][M][N];
        fftwnd_plan plan;
        int n[3] = {L,M,N};

        plan = fftwnd_create_plan(3,n,FFTW_FORWARD,
                                  FFTW_ESTIMATE | FFTW_IN_PLACE);
        fftwnd_one(plan, arr, 0);
        fftwnd_destroy_plan(plan);

In Fortran, you would use this instead:

        double complex arr
        dimension arr(L,M,N)
        integer n
        dimension n(3)
        integer plan

        n(1) = L
        n(2) = M
        n(3) = N
        call fftwnd_f77_create_plan(plan,3,n,FFTW_FORWARD,
       +                            FFTW_ESTIMATE + FFTW_IN_PLACE)
        call fftwnd_f77_one(plan, arr, 0)
        call fftwnd_f77_destroy_plan(plan)

Instead of calling fftwnd_f77_create_plan(plan,3,n,...), we could also have called fftw3d_f77_create_plan(plan,L,M,N,...).

Note that we pass the array dimensions in the "natural" order; also note that the last argument to fftwnd_f77 is ignored since the transform is FFTW_IN_PLACE.

To transform a one-dimensional real array in Fortran, you might do:

        double precision in, out
        dimension in(N), out(N)
        integer plan

        call rfftw_f77_create_plan(plan,N,FFTW_REAL_TO_COMPLEX,
       +                           FFTW_ESTIMATE)
        call rfftw_f77_one(plan,in,out)
        call rfftw_f77_destroy_plan(plan)

To transform a two-dimensional real array, out of place, you might use the following:

        double precision in
        double complex out
        dimension in(M,N), out(M/2 + 1, N)
        integer plan

        call rfftw2d_f77_create_plan(plan,M,N,FFTW_REAL_TO_COMPLEX,
       +                             FFTW_ESTIMATE)
        call rfftwnd_f77_one_real_to_complex(plan, in, out)
        call rfftwnd_f77_destroy_plan(plan)

Important: Notice that it is the first dimension of the complex output array that is cut in half in Fortran, rather than the last dimension as in C. This is a consequence of the wrapper routines reversing the order of the array dimensions passed to FFTW so that the Fortran program can use its ordinary column-major order.


Go to the first, previous, next, last section, table of contents.