Chris@82: Chris@82: Chris@82: Chris@82: Chris@82:
Chris@82:Chris@82: Next: Wisdom of Fortran?, Previous: FFTW Execution in Fortran, Up: Calling FFTW from Legacy Fortran [Contents][Index]
Chris@82:In C, you might have something like the following to transform a Chris@82: one-dimensional complex array: Chris@82:
Chris@82:fftw_complex in[N], out[N]; Chris@82: fftw_plan plan; Chris@82: Chris@82: plan = fftw_plan_dft_1d(N,in,out,FFTW_FORWARD,FFTW_ESTIMATE); Chris@82: fftw_execute(plan); Chris@82: fftw_destroy_plan(plan); Chris@82:
In Fortran, you would use the following to accomplish the same thing: Chris@82:
Chris@82:double complex in, out Chris@82: dimension in(N), out(N) Chris@82: integer*8 plan Chris@82: Chris@82: call dfftw_plan_dft_1d(plan,N,in,out,FFTW_FORWARD,FFTW_ESTIMATE) Chris@82: call dfftw_execute_dft(plan, in, out) Chris@82: call dfftw_destroy_plan(plan) Chris@82:
Notice how all routines are called as Fortran subroutines, and the
Chris@82: plan is returned via the first argument to dfftw_plan_dft_1d
.
Chris@82: Notice also that we changed fftw_execute
to
Chris@82: dfftw_execute_dft
(see FFTW Execution in Fortran). To do
Chris@82: the same thing, but using 8 threads in parallel (see Multi-threaded FFTW), you would simply prefix these calls with:
Chris@82:
integer iret Chris@82: call dfftw_init_threads(iret) Chris@82: call dfftw_plan_with_nthreads(8) Chris@82:
(You might want to check the value of iret
: if it is zero, it
Chris@82: indicates an unlikely error during thread initialization.)
Chris@82:
To transform a three-dimensional array in-place with C, you might do: Chris@82:
Chris@82:fftw_complex arr[L][M][N]; Chris@82: fftw_plan plan; Chris@82: Chris@82: plan = fftw_plan_dft_3d(L,M,N, arr,arr, Chris@82: FFTW_FORWARD, FFTW_ESTIMATE); Chris@82: fftw_execute(plan); Chris@82: fftw_destroy_plan(plan); Chris@82:
In Fortran, you would use this instead: Chris@82:
Chris@82:double complex arr Chris@82: dimension arr(L,M,N) Chris@82: integer*8 plan Chris@82: Chris@82: call dfftw_plan_dft_3d(plan, L,M,N, arr,arr, Chris@82: & FFTW_FORWARD, FFTW_ESTIMATE) Chris@82: call dfftw_execute_dft(plan, arr, arr) Chris@82: call dfftw_destroy_plan(plan) Chris@82:
Note that we pass the array dimensions in the “natural” order in both C Chris@82: and Fortran. Chris@82:
Chris@82:To transform a one-dimensional real array in Fortran, you might do: Chris@82:
Chris@82:double precision in Chris@82: dimension in(N) Chris@82: double complex out Chris@82: dimension out(N/2 + 1) Chris@82: integer*8 plan Chris@82: Chris@82: call dfftw_plan_dft_r2c_1d(plan,N,in,out,FFTW_ESTIMATE) Chris@82: call dfftw_execute_dft_r2c(plan, in, out) Chris@82: call dfftw_destroy_plan(plan) Chris@82:
To transform a two-dimensional real array, out of place, you might use Chris@82: the following: Chris@82:
Chris@82:double precision in Chris@82: dimension in(M,N) Chris@82: double complex out Chris@82: dimension out(M/2 + 1, N) Chris@82: integer*8 plan Chris@82: Chris@82: call dfftw_plan_dft_r2c_2d(plan,M,N,in,out,FFTW_ESTIMATE) Chris@82: call dfftw_execute_dft_r2c(plan, in, out) Chris@82: call dfftw_destroy_plan(plan) Chris@82:
Important: Notice that it is the first dimension of the Chris@82: complex output array that is cut in half in Fortran, rather than the Chris@82: last dimension as in C. This is a consequence of the interface routines Chris@82: reversing the order of the array dimensions passed to FFTW so that the Chris@82: Fortran program can use its ordinary column-major order. Chris@82: Chris@82: Chris@82:
Chris@82:Chris@82: Next: Wisdom of Fortran?, Previous: FFTW Execution in Fortran, Up: Calling FFTW from Legacy Fortran [Contents][Index]
Chris@82: