fftw3: FFTW Execution in Fortran

 
 8.3 FFTW Execution in Fortran
 =============================
 
 In C, in order to use a plan, one normally calls 'fftw_execute', which
 executes the plan to perform the transform on the input/output arrays
 passed when the plan was created (SeeUsing Plans).  The
 corresponding subroutine call in legacy Fortran is:
              call dfftw_execute(plan)
 
    However, we have had reports that this causes problems with some
 recent optimizing Fortran compilers.  The problem is, because the
 input/output arrays are not passed as explicit arguments to
 'dfftw_execute', the semantics of Fortran (unlike C) allow the compiler
 to assume that the input/output arrays are not changed by
 'dfftw_execute'.  As a consequence, certain compilers end up optimizing
 out or repositioning the call to 'dfftw_execute', assuming incorrectly
 that it does nothing.
 
    There are various workarounds to this, but the safest and simplest
 thing is to not use 'dfftw_execute' in Fortran.  Instead, use the
 functions described in SeeNew-array Execute Functions, which take
 the input/output arrays as explicit arguments.  For example, if the plan
 is for a complex-data DFT and was created for the arrays 'in' and 'out',
 you would do:
              call dfftw_execute_dft(plan, in, out)
 
    There are a few things to be careful of, however:
 
    * You must use the correct type of execute function, matching the way
      the plan was created.  Complex DFT plans should use
      'dfftw_execute_dft', Real-input (r2c) DFT plans should use use
      'dfftw_execute_dft_r2c', and real-output (c2r) DFT plans should use
      'dfftw_execute_dft_c2r'.  The various r2r plans should use
      'dfftw_execute_r2r'.
 
    * You should normally pass the same input/output arrays that were
      used when creating the plan.  This is always safe.
 
    * _If_ you pass _different_ input/output arrays compared to those
      used when creating the plan, you must abide by all the restrictions
      of the new-array execute functions (SeeNew-array Execute
      Functions).  The most difficult of these, in Fortran, is the
      requirement that the new arrays have the same alignment as the
      original arrays, because there seems to be no way in legacy Fortran
      to obtain guaranteed-aligned arrays (analogous to 'fftw_malloc' in
      C). You can, of course, use the 'FFTW_UNALIGNED' flag when creating
      the plan, in which case the plan does not depend on the alignment,
      but this may sacrifice substantial performance on architectures
      (like x86) with SIMD instructions (SeeSIMD alignment and
      fftw_malloc).