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