| Chris@42 | 1 <!DOCTYPE html PUBLIC "-//W3C//DTD HTML 4.01 Transitional//EN" "http://www.w3.org/TR/html4/loose.dtd"> | 
| Chris@42 | 2 <html> | 
| Chris@42 | 3 <!-- This manual is for FFTW | 
| Chris@42 | 4 (version 3.3.5, 30 July 2016). | 
| Chris@42 | 5 | 
| Chris@42 | 6 Copyright (C) 2003 Matteo Frigo. | 
| Chris@42 | 7 | 
| Chris@42 | 8 Copyright (C) 2003 Massachusetts Institute of Technology. | 
| Chris@42 | 9 | 
| Chris@42 | 10 Permission is granted to make and distribute verbatim copies of this | 
| Chris@42 | 11 manual provided the copyright notice and this permission notice are | 
| Chris@42 | 12 preserved on all copies. | 
| Chris@42 | 13 | 
| Chris@42 | 14 Permission is granted to copy and distribute modified versions of this | 
| Chris@42 | 15 manual under the conditions for verbatim copying, provided that the | 
| Chris@42 | 16 entire resulting derived work is distributed under the terms of a | 
| Chris@42 | 17 permission notice identical to this one. | 
| Chris@42 | 18 | 
| Chris@42 | 19 Permission is granted to copy and distribute translations of this manual | 
| Chris@42 | 20 into another language, under the above conditions for modified versions, | 
| Chris@42 | 21 except that this permission notice may be stated in a translation | 
| Chris@42 | 22 approved by the Free Software Foundation. --> | 
| Chris@42 | 23 <!-- Created by GNU Texinfo 5.2, http://www.gnu.org/software/texinfo/ --> | 
| Chris@42 | 24 <head> | 
| Chris@42 | 25 <title>FFTW 3.3.5: FFTW MPI Fortran Interface</title> | 
| Chris@42 | 26 | 
| Chris@42 | 27 <meta name="description" content="FFTW 3.3.5: FFTW MPI Fortran Interface"> | 
| Chris@42 | 28 <meta name="keywords" content="FFTW 3.3.5: FFTW MPI Fortran Interface"> | 
| Chris@42 | 29 <meta name="resource-type" content="document"> | 
| Chris@42 | 30 <meta name="distribution" content="global"> | 
| Chris@42 | 31 <meta name="Generator" content="makeinfo"> | 
| Chris@42 | 32 <meta http-equiv="Content-Type" content="text/html; charset=utf-8"> | 
| Chris@42 | 33 <link href="index.html#Top" rel="start" title="Top"> | 
| Chris@42 | 34 <link href="Concept-Index.html#Concept-Index" rel="index" title="Concept Index"> | 
| Chris@42 | 35 <link href="index.html#SEC_Contents" rel="contents" title="Table of Contents"> | 
| Chris@42 | 36 <link href="Distributed_002dmemory-FFTW-with-MPI.html#Distributed_002dmemory-FFTW-with-MPI" rel="up" title="Distributed-memory FFTW with MPI"> | 
| Chris@42 | 37 <link href="Calling-FFTW-from-Modern-Fortran.html#Calling-FFTW-from-Modern-Fortran" rel="next" title="Calling FFTW from Modern Fortran"> | 
| Chris@42 | 38 <link href="MPI-Wisdom-Communication.html#MPI-Wisdom-Communication" rel="prev" title="MPI Wisdom Communication"> | 
| Chris@42 | 39 <style type="text/css"> | 
| Chris@42 | 40 <!-- | 
| Chris@42 | 41 a.summary-letter {text-decoration: none} | 
| Chris@42 | 42 blockquote.smallquotation {font-size: smaller} | 
| Chris@42 | 43 div.display {margin-left: 3.2em} | 
| Chris@42 | 44 div.example {margin-left: 3.2em} | 
| Chris@42 | 45 div.indentedblock {margin-left: 3.2em} | 
| Chris@42 | 46 div.lisp {margin-left: 3.2em} | 
| Chris@42 | 47 div.smalldisplay {margin-left: 3.2em} | 
| Chris@42 | 48 div.smallexample {margin-left: 3.2em} | 
| Chris@42 | 49 div.smallindentedblock {margin-left: 3.2em; font-size: smaller} | 
| Chris@42 | 50 div.smalllisp {margin-left: 3.2em} | 
| Chris@42 | 51 kbd {font-style:oblique} | 
| Chris@42 | 52 pre.display {font-family: inherit} | 
| Chris@42 | 53 pre.format {font-family: inherit} | 
| Chris@42 | 54 pre.menu-comment {font-family: serif} | 
| Chris@42 | 55 pre.menu-preformatted {font-family: serif} | 
| Chris@42 | 56 pre.smalldisplay {font-family: inherit; font-size: smaller} | 
| Chris@42 | 57 pre.smallexample {font-size: smaller} | 
| Chris@42 | 58 pre.smallformat {font-family: inherit; font-size: smaller} | 
| Chris@42 | 59 pre.smalllisp {font-size: smaller} | 
| Chris@42 | 60 span.nocodebreak {white-space:nowrap} | 
| Chris@42 | 61 span.nolinebreak {white-space:nowrap} | 
| Chris@42 | 62 span.roman {font-family:serif; font-weight:normal} | 
| Chris@42 | 63 span.sansserif {font-family:sans-serif; font-weight:normal} | 
| Chris@42 | 64 ul.no-bullet {list-style: none} | 
| Chris@42 | 65 --> | 
| Chris@42 | 66 </style> | 
| Chris@42 | 67 | 
| Chris@42 | 68 | 
| Chris@42 | 69 </head> | 
| Chris@42 | 70 | 
| Chris@42 | 71 <body lang="en" bgcolor="#FFFFFF" text="#000000" link="#0000FF" vlink="#800080" alink="#FF0000"> | 
| Chris@42 | 72 <a name="FFTW-MPI-Fortran-Interface"></a> | 
| Chris@42 | 73 <div class="header"> | 
| Chris@42 | 74 <p> | 
| Chris@42 | 75 Previous: <a href="FFTW-MPI-Reference.html#FFTW-MPI-Reference" accesskey="p" rel="prev">FFTW MPI Reference</a>, Up: <a href="Distributed_002dmemory-FFTW-with-MPI.html#Distributed_002dmemory-FFTW-with-MPI" accesskey="u" rel="up">Distributed-memory FFTW with MPI</a>   [<a href="index.html#SEC_Contents" title="Table of contents" rel="contents">Contents</a>][<a href="Concept-Index.html#Concept-Index" title="Index" rel="index">Index</a>]</p> | 
| Chris@42 | 76 </div> | 
| Chris@42 | 77 <hr> | 
| Chris@42 | 78 <a name="FFTW-MPI-Fortran-Interface-1"></a> | 
| Chris@42 | 79 <h3 class="section">6.13 FFTW MPI Fortran Interface</h3> | 
| Chris@42 | 80 <a name="index-Fortran-interface-1"></a> | 
| Chris@42 | 81 | 
| Chris@42 | 82 <a name="index-iso_005fc_005fbinding"></a> | 
| Chris@42 | 83 <p>The FFTW MPI interface is callable from modern Fortran compilers | 
| Chris@42 | 84 supporting the Fortran 2003 <code>iso_c_binding</code> standard for calling | 
| Chris@42 | 85 C functions.  As described in <a href="Calling-FFTW-from-Modern-Fortran.html#Calling-FFTW-from-Modern-Fortran">Calling FFTW from Modern Fortran</a>, | 
| Chris@42 | 86 this means that you can directly call FFTW’s C interface from Fortran | 
| Chris@42 | 87 with only minor changes in syntax.  There are, however, a few things | 
| Chris@42 | 88 specific to the MPI interface to keep in mind: | 
| Chris@42 | 89 </p> | 
| Chris@42 | 90 <ul> | 
| Chris@42 | 91 <li> Instead of including <code>fftw3.f03</code> as in <a href="Overview-of-Fortran-interface.html#Overview-of-Fortran-interface">Overview of Fortran interface</a>, you should <code>include 'fftw3-mpi.f03'</code> (after | 
| Chris@42 | 92 <code>use, intrinsic :: iso_c_binding</code> as before).  The | 
| Chris@42 | 93 <code>fftw3-mpi.f03</code> file includes <code>fftw3.f03</code>, so you should | 
| Chris@42 | 94 <em>not</em> <code>include</code> them both yourself.  (You will also want to | 
| Chris@42 | 95 include the MPI header file, usually via <code>include 'mpif.h'</code> or | 
| Chris@42 | 96 similar, although though this is not needed by <code>fftw3-mpi.f03</code> | 
| Chris@42 | 97 <i>per se</i>.)  (To use the ‘<samp>fftwl_</samp>’ <code>long double</code> extended-precision routines in supporting compilers, you should include <code>fftw3f-mpi.f03</code> in <em>addition</em> to <code>fftw3-mpi.f03</code>. See <a href="Extended-and-quadruple-precision-in-Fortran.html#Extended-and-quadruple-precision-in-Fortran">Extended and quadruple precision in Fortran</a>.) | 
| Chris@42 | 98 | 
| Chris@42 | 99 </li><li> Because of the different storage conventions between C and Fortran, | 
| Chris@42 | 100 you reverse the order of your array dimensions when passing them to | 
| Chris@42 | 101 FFTW (see <a href="Reversing-array-dimensions.html#Reversing-array-dimensions">Reversing array dimensions</a>).  This is merely a | 
| Chris@42 | 102 difference in notation and incurs no performance overhead.  However, | 
| Chris@42 | 103 it means that, whereas in C the <em>first</em> dimension is distributed, | 
| Chris@42 | 104 in Fortran the <em>last</em> dimension of your array is distributed. | 
| Chris@42 | 105 | 
| Chris@42 | 106 </li><li> <a name="index-MPI-communicator-3"></a> | 
| Chris@42 | 107 In Fortran, communicators are stored as <code>integer</code> types; there is | 
| Chris@42 | 108 no <code>MPI_Comm</code> type, nor is there any way to access a C | 
| Chris@42 | 109 <code>MPI_Comm</code>.  Fortunately, this is taken care of for you by the | 
| Chris@42 | 110 FFTW Fortran interface: whenever the C interface expects an | 
| Chris@42 | 111 <code>MPI_Comm</code> type, you should pass the Fortran communicator as an | 
| Chris@42 | 112 <code>integer</code>.<a name="DOCF8" href="#FOOT8"><sup>8</sup></a> | 
| Chris@42 | 113 | 
| Chris@42 | 114 </li><li> Because you need to call the ‘<samp>local_size</samp>’ function to find out | 
| Chris@42 | 115 how much space to allocate, and this may be <em>larger</em> than the | 
| Chris@42 | 116 local portion of the array (see <a href="MPI-Data-Distribution.html#MPI-Data-Distribution">MPI Data Distribution</a>), you should | 
| Chris@42 | 117 <em>always</em> allocate your arrays dynamically using FFTW’s allocation | 
| Chris@42 | 118 routines as described in <a href="Allocating-aligned-memory-in-Fortran.html#Allocating-aligned-memory-in-Fortran">Allocating aligned memory in Fortran</a>. | 
| Chris@42 | 119 (Coincidentally, this also provides the best performance by | 
| Chris@42 | 120 guaranteeding proper data alignment.) | 
| Chris@42 | 121 | 
| Chris@42 | 122 </li><li> Because all sizes in the MPI FFTW interface are declared as | 
| Chris@42 | 123 <code>ptrdiff_t</code> in C, you should use <code>integer(C_INTPTR_T)</code> in | 
| Chris@42 | 124 Fortran (see <a href="FFTW-Fortran-type-reference.html#FFTW-Fortran-type-reference">FFTW Fortran type reference</a>). | 
| Chris@42 | 125 | 
| Chris@42 | 126 </li><li> <a name="index-fftw_005fexecute_005fdft-1"></a> | 
| Chris@42 | 127 <a name="index-fftw_005fmpi_005fexecute_005fdft-1"></a> | 
| Chris@42 | 128 <a name="index-new_002darray-execution-3"></a> | 
| Chris@42 | 129 In Fortran, because of the language semantics, we generally recommend | 
| Chris@42 | 130 using the new-array execute functions for all plans, even in the | 
| Chris@42 | 131 common case where you are executing the plan on the same arrays for | 
| Chris@42 | 132 which the plan was created (see <a href="Plan-execution-in-Fortran.html#Plan-execution-in-Fortran">Plan execution in Fortran</a>). | 
| Chris@42 | 133 However, note that in the MPI interface these functions are changed: | 
| Chris@42 | 134 <code>fftw_execute_dft</code> becomes <code>fftw_mpi_execute_dft</code>, | 
| Chris@42 | 135 etcetera. See <a href="Using-MPI-Plans.html#Using-MPI-Plans">Using MPI Plans</a>. | 
| Chris@42 | 136 | 
| Chris@42 | 137 </li></ul> | 
| Chris@42 | 138 | 
| Chris@42 | 139 <p>For example, here is a Fortran code snippet to perform a distributed | 
| Chris@42 | 140 L × M complex DFT in-place.  (This assumes you have already | 
| Chris@42 | 141 initialized MPI with <code>MPI_init</code> and have also performed | 
| Chris@42 | 142 <code>call fftw_mpi_init</code>.) | 
| Chris@42 | 143 </p> | 
| Chris@42 | 144 <div class="example"> | 
| Chris@42 | 145 <pre class="example">  use, intrinsic :: iso_c_binding | 
| Chris@42 | 146   include 'fftw3-mpi.f03' | 
| Chris@42 | 147   integer(C_INTPTR_T), parameter :: L = ... | 
| Chris@42 | 148   integer(C_INTPTR_T), parameter :: M = ... | 
| Chris@42 | 149   type(C_PTR) :: plan, cdata | 
| Chris@42 | 150   complex(C_DOUBLE_COMPLEX), pointer :: data(:,:) | 
| Chris@42 | 151   integer(C_INTPTR_T) :: i, j, alloc_local, local_M, local_j_offset | 
| Chris@42 | 152 | 
| Chris@42 | 153 !   <span class="roman">get local data size and allocate (note dimension reversal)</span> | 
| Chris@42 | 154   alloc_local = fftw_mpi_local_size_2d(M, L, MPI_COMM_WORLD, & | 
| Chris@42 | 155                                        local_M, local_j_offset) | 
| Chris@42 | 156   cdata = fftw_alloc_complex(alloc_local) | 
| Chris@42 | 157   call c_f_pointer(cdata, data, [L,local_M]) | 
| Chris@42 | 158 | 
| Chris@42 | 159 !   <span class="roman">create MPI plan for in-place forward DFT (note dimension reversal)</span> | 
| Chris@42 | 160   plan = fftw_mpi_plan_dft_2d(M, L, data, data, MPI_COMM_WORLD, & | 
| Chris@42 | 161                               FFTW_FORWARD, FFTW_MEASURE) | 
| Chris@42 | 162 | 
| Chris@42 | 163 ! <span class="roman">initialize data to some function</span> my_function(i,j) | 
| Chris@42 | 164   do j = 1, local_M | 
| Chris@42 | 165     do i = 1, L | 
| Chris@42 | 166       data(i, j) = my_function(i, j + local_j_offset) | 
| Chris@42 | 167     end do | 
| Chris@42 | 168   end do | 
| Chris@42 | 169 | 
| Chris@42 | 170 ! <span class="roman">compute transform (as many times as desired)</span> | 
| Chris@42 | 171   call fftw_mpi_execute_dft(plan, data, data) | 
| Chris@42 | 172 | 
| Chris@42 | 173   call fftw_destroy_plan(plan) | 
| Chris@42 | 174   call fftw_free(cdata) | 
| Chris@42 | 175 </pre></div> | 
| Chris@42 | 176 | 
| Chris@42 | 177 <p>Note that when we called <code>fftw_mpi_local_size_2d</code> and | 
| Chris@42 | 178 <code>fftw_mpi_plan_dft_2d</code> with the dimensions in reversed order, | 
| Chris@42 | 179 since a L × M Fortran array is viewed by FFTW in C as a | 
| Chris@42 | 180 M × L array.  This means that the array was distributed over | 
| Chris@42 | 181 the <code>M</code> dimension, the local portion of which is a | 
| Chris@42 | 182 L × local_M array in Fortran.  (You must <em>not</em> use an | 
| Chris@42 | 183 <code>allocate</code> statement to allocate an L × local_M array, | 
| Chris@42 | 184 however; you must allocate <code>alloc_local</code> complex numbers, which | 
| Chris@42 | 185 may be greater than <code>L * local_M</code>, in order to reserve space for | 
| Chris@42 | 186 intermediate steps of the transform.)  Finally, we mention that | 
| Chris@42 | 187 because C’s array indices are zero-based, the <code>local_j_offset</code> | 
| Chris@42 | 188 argument can conveniently be interpreted as an offset in the 1-based | 
| Chris@42 | 189 <code>j</code> index (rather than as a starting index as in C). | 
| Chris@42 | 190 </p> | 
| Chris@42 | 191 <p>If instead you had used the <code>ior(FFTW_MEASURE, | 
| Chris@42 | 192 FFTW_MPI_TRANSPOSED_OUT)</code> flag, the output of the transform would be a | 
| Chris@42 | 193 transposed M × local_L array, associated with the <em>same</em> | 
| Chris@42 | 194 <code>cdata</code> allocation (since the transform is in-place), and which | 
| Chris@42 | 195 you could declare with: | 
| Chris@42 | 196 </p> | 
| Chris@42 | 197 <div class="example"> | 
| Chris@42 | 198 <pre class="example">  complex(C_DOUBLE_COMPLEX), pointer :: tdata(:,:) | 
| Chris@42 | 199   ... | 
| Chris@42 | 200   call c_f_pointer(cdata, tdata, [M,local_L]) | 
| Chris@42 | 201 </pre></div> | 
| Chris@42 | 202 | 
| Chris@42 | 203 <p>where <code>local_L</code> would have been obtained by changing the | 
| Chris@42 | 204 <code>fftw_mpi_local_size_2d</code> call to: | 
| Chris@42 | 205 </p> | 
| Chris@42 | 206 <div class="example"> | 
| Chris@42 | 207 <pre class="example">  alloc_local = fftw_mpi_local_size_2d_transposed(M, L, MPI_COMM_WORLD, & | 
| Chris@42 | 208                            local_M, local_j_offset, local_L, local_i_offset) | 
| Chris@42 | 209 </pre></div> | 
| Chris@42 | 210 <div class="footnote"> | 
| Chris@42 | 211 <hr> | 
| Chris@42 | 212 <h4 class="footnotes-heading">Footnotes</h4> | 
| Chris@42 | 213 | 
| Chris@42 | 214 <h3><a name="FOOT8" href="#DOCF8">(8)</a></h3> | 
| Chris@42 | 215 <p>Technically, this is because you aren’t | 
| Chris@42 | 216 actually calling the C functions directly. You are calling wrapper | 
| Chris@42 | 217 functions that translate the communicator with <code>MPI_Comm_f2c</code> | 
| Chris@42 | 218 before calling the ordinary C interface.  This is all done | 
| Chris@42 | 219 transparently, however, since the <code>fftw3-mpi.f03</code> interface file | 
| Chris@42 | 220 renames the wrappers so that they are called in Fortran with the same | 
| Chris@42 | 221 names as the C interface functions.</p> | 
| Chris@42 | 222 </div> | 
| Chris@42 | 223 <hr> | 
| Chris@42 | 224 <div class="header"> | 
| Chris@42 | 225 <p> | 
| Chris@42 | 226 Previous: <a href="FFTW-MPI-Reference.html#FFTW-MPI-Reference" accesskey="p" rel="prev">FFTW MPI Reference</a>, Up: <a href="Distributed_002dmemory-FFTW-with-MPI.html#Distributed_002dmemory-FFTW-with-MPI" accesskey="u" rel="up">Distributed-memory FFTW with MPI</a>   [<a href="index.html#SEC_Contents" title="Table of contents" rel="contents">Contents</a>][<a href="Concept-Index.html#Concept-Index" title="Index" rel="index">Index</a>]</p> | 
| Chris@42 | 227 </div> | 
| Chris@42 | 228 | 
| Chris@42 | 229 | 
| Chris@42 | 230 | 
| Chris@42 | 231 </body> | 
| Chris@42 | 232 </html> |