Chris@42: <!DOCTYPE html PUBLIC "-//W3C//DTD HTML 4.01 Transitional//EN" "http://www.w3.org/TR/html4/loose.dtd">
Chris@42: <html>
Chris@42: <!-- This manual is for FFTW
Chris@42: (version 3.3.5, 30 July 2016).
Chris@42: 
Chris@42: Copyright (C) 2003 Matteo Frigo.
Chris@42: 
Chris@42: Copyright (C) 2003 Massachusetts Institute of Technology.
Chris@42: 
Chris@42: Permission is granted to make and distribute verbatim copies of this
Chris@42: manual provided the copyright notice and this permission notice are
Chris@42: preserved on all copies.
Chris@42: 
Chris@42: Permission is granted to copy and distribute modified versions of this
Chris@42: manual under the conditions for verbatim copying, provided that the
Chris@42: entire resulting derived work is distributed under the terms of a
Chris@42: permission notice identical to this one.
Chris@42: 
Chris@42: Permission is granted to copy and distribute translations of this manual
Chris@42: into another language, under the above conditions for modified versions,
Chris@42: except that this permission notice may be stated in a translation
Chris@42: approved by the Free Software Foundation. -->
Chris@42: <!-- Created by GNU Texinfo 5.2, http://www.gnu.org/software/texinfo/ -->
Chris@42: <head>
Chris@42: <title>FFTW 3.3.5: 2d MPI example</title>
Chris@42: 
Chris@42: <meta name="description" content="FFTW 3.3.5: 2d MPI example">
Chris@42: <meta name="keywords" content="FFTW 3.3.5: 2d MPI example">
Chris@42: <meta name="resource-type" content="document">
Chris@42: <meta name="distribution" content="global">
Chris@42: <meta name="Generator" content="makeinfo">
Chris@42: <meta http-equiv="Content-Type" content="text/html; charset=utf-8">
Chris@42: <link href="index.html#Top" rel="start" title="Top">
Chris@42: <link href="Concept-Index.html#Concept-Index" rel="index" title="Concept Index">
Chris@42: <link href="index.html#SEC_Contents" rel="contents" title="Table of Contents">
Chris@42: <link href="Distributed_002dmemory-FFTW-with-MPI.html#Distributed_002dmemory-FFTW-with-MPI" rel="up" title="Distributed-memory FFTW with MPI">
Chris@42: <link href="MPI-Data-Distribution.html#MPI-Data-Distribution" rel="next" title="MPI Data Distribution">
Chris@42: <link href="Linking-and-Initializing-MPI-FFTW.html#Linking-and-Initializing-MPI-FFTW" rel="prev" title="Linking and Initializing MPI FFTW">
Chris@42: <style type="text/css">
Chris@42: <!--
Chris@42: a.summary-letter {text-decoration: none}
Chris@42: blockquote.smallquotation {font-size: smaller}
Chris@42: div.display {margin-left: 3.2em}
Chris@42: div.example {margin-left: 3.2em}
Chris@42: div.indentedblock {margin-left: 3.2em}
Chris@42: div.lisp {margin-left: 3.2em}
Chris@42: div.smalldisplay {margin-left: 3.2em}
Chris@42: div.smallexample {margin-left: 3.2em}
Chris@42: div.smallindentedblock {margin-left: 3.2em; font-size: smaller}
Chris@42: div.smalllisp {margin-left: 3.2em}
Chris@42: kbd {font-style:oblique}
Chris@42: pre.display {font-family: inherit}
Chris@42: pre.format {font-family: inherit}
Chris@42: pre.menu-comment {font-family: serif}
Chris@42: pre.menu-preformatted {font-family: serif}
Chris@42: pre.smalldisplay {font-family: inherit; font-size: smaller}
Chris@42: pre.smallexample {font-size: smaller}
Chris@42: pre.smallformat {font-family: inherit; font-size: smaller}
Chris@42: pre.smalllisp {font-size: smaller}
Chris@42: span.nocodebreak {white-space:nowrap}
Chris@42: span.nolinebreak {white-space:nowrap}
Chris@42: span.roman {font-family:serif; font-weight:normal}
Chris@42: span.sansserif {font-family:sans-serif; font-weight:normal}
Chris@42: ul.no-bullet {list-style: none}
Chris@42: -->
Chris@42: </style>
Chris@42: 
Chris@42: 
Chris@42: </head>
Chris@42: 
Chris@42: <body lang="en" bgcolor="#FFFFFF" text="#000000" link="#0000FF" vlink="#800080" alink="#FF0000">
Chris@42: <a name="g_t2d-MPI-example"></a>
Chris@42: <div class="header">
Chris@42: <p>
Chris@42: Next: <a href="MPI-Data-Distribution.html#MPI-Data-Distribution" accesskey="n" rel="next">MPI Data Distribution</a>, Previous: <a href="Linking-and-Initializing-MPI-FFTW.html#Linking-and-Initializing-MPI-FFTW" accesskey="p" rel="prev">Linking and Initializing MPI FFTW</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> &nbsp; [<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: </div>
Chris@42: <hr>
Chris@42: <a name="g_t2d-MPI-example-1"></a>
Chris@42: <h3 class="section">6.3 2d MPI example</h3>
Chris@42: 
Chris@42: <p>Before we document the FFTW MPI interface in detail, we begin with a
Chris@42: simple example outlining how one would perform a two-dimensional
Chris@42: <code>N0</code> by <code>N1</code> complex DFT. 
Chris@42: </p>
Chris@42: <div class="example">
Chris@42: <pre class="example">#include &lt;fftw3-mpi.h&gt;
Chris@42: 
Chris@42: int main(int argc, char **argv)
Chris@42: {
Chris@42:     const ptrdiff_t N0 = ..., N1 = ...;
Chris@42:     fftw_plan plan;
Chris@42:     fftw_complex *data;
Chris@42:     ptrdiff_t alloc_local, local_n0, local_0_start, i, j;
Chris@42: 
Chris@42:     MPI_Init(&amp;argc, &amp;argv);
Chris@42:     fftw_mpi_init();
Chris@42: 
Chris@42:     /* <span class="roman">get local data size and allocate</span> */
Chris@42:     alloc_local = fftw_mpi_local_size_2d(N0, N1, MPI_COMM_WORLD,
Chris@42:                                          &amp;local_n0, &amp;local_0_start);
Chris@42:     data = fftw_alloc_complex(alloc_local);
Chris@42: 
Chris@42:     /* <span class="roman">create plan for in-place forward DFT</span> */
Chris@42:     plan = fftw_mpi_plan_dft_2d(N0, N1, data, data, MPI_COMM_WORLD,
Chris@42:                                 FFTW_FORWARD, FFTW_ESTIMATE);    
Chris@42: 
Chris@42:     /* <span class="roman">initialize data to some function</span> my_function(x,y) */
Chris@42:     for (i = 0; i &lt; local_n0; ++i) for (j = 0; j &lt; N1; ++j)
Chris@42:        data[i*N1 + j] = my_function(local_0_start + i, j);
Chris@42: 
Chris@42:     /* <span class="roman">compute transforms, in-place, as many times as desired</span> */
Chris@42:     fftw_execute(plan);
Chris@42: 
Chris@42:     fftw_destroy_plan(plan);
Chris@42: 
Chris@42:     MPI_Finalize();
Chris@42: }
Chris@42: </pre></div>
Chris@42: 
Chris@42: <p>As can be seen above, the MPI interface follows the same basic style
Chris@42: of allocate/plan/execute/destroy as the serial FFTW routines.  All of
Chris@42: the MPI-specific routines are prefixed with &lsquo;<samp>fftw_mpi_</samp>&rsquo; instead
Chris@42: of &lsquo;<samp>fftw_</samp>&rsquo;.  There are a few important differences, however:
Chris@42: </p>
Chris@42: <p>First, we must call <code>fftw_mpi_init()</code> after calling
Chris@42: <code>MPI_Init</code> (required in all MPI programs) and before calling any
Chris@42: other &lsquo;<samp>fftw_mpi_</samp>&rsquo; routine.
Chris@42: <a name="index-MPI_005fInit"></a>
Chris@42: <a name="index-fftw_005fmpi_005finit-1"></a>
Chris@42: </p>
Chris@42: 
Chris@42: <p>Second, when we create the plan with <code>fftw_mpi_plan_dft_2d</code>,
Chris@42: analogous to <code>fftw_plan_dft_2d</code>, we pass an additional argument:
Chris@42: the communicator, indicating which processes will participate in the
Chris@42: transform (here <code>MPI_COMM_WORLD</code>, indicating all processes).
Chris@42: Whenever you create, execute, or destroy a plan for an MPI transform,
Chris@42: you must call the corresponding FFTW routine on <em>all</em> processes
Chris@42: in the communicator for that transform.  (That is, these are
Chris@42: <em>collective</em> calls.)  Note that the plan for the MPI transform
Chris@42: uses the standard <code>fftw_execute</code> and <code>fftw_destroy</code> routines
Chris@42: (on the other hand, there are MPI-specific new-array execute functions
Chris@42: documented below).
Chris@42: <a name="index-collective-function"></a>
Chris@42: <a name="index-fftw_005fmpi_005fplan_005fdft_005f2d"></a>
Chris@42: <a name="index-MPI_005fCOMM_005fWORLD-1"></a>
Chris@42: </p>
Chris@42: 
Chris@42: <p>Third, all of the FFTW MPI routines take <code>ptrdiff_t</code> arguments
Chris@42: instead of <code>int</code> as for the serial FFTW.  <code>ptrdiff_t</code> is a
Chris@42: standard C integer type which is (at least) 32 bits wide on a 32-bit
Chris@42: machine and 64 bits wide on a 64-bit machine.  This is to make it easy
Chris@42: to specify very large parallel transforms on a 64-bit machine.  (You
Chris@42: can specify 64-bit transform sizes in the serial FFTW, too, but only
Chris@42: by using the &lsquo;<samp>guru64</samp>&rsquo; planner interface.  See <a href="64_002dbit-Guru-Interface.html#g_t64_002dbit-Guru-Interface">64-bit Guru Interface</a>.)
Chris@42: <a name="index-ptrdiff_005ft-1"></a>
Chris@42: <a name="index-64_002dbit-architecture-1"></a>
Chris@42: </p>
Chris@42: 
Chris@42: <p>Fourth, and most importantly, you don&rsquo;t allocate the entire
Chris@42: two-dimensional array on each process.  Instead, you call
Chris@42: <code>fftw_mpi_local_size_2d</code> to find out what <em>portion</em> of the
Chris@42: array resides on each processor, and how much space to allocate.
Chris@42: Here, the portion of the array on each process is a <code>local_n0</code> by
Chris@42: <code>N1</code> slice of the total array, starting at index
Chris@42: <code>local_0_start</code>.  The total number of <code>fftw_complex</code> numbers
Chris@42: to allocate is given by the <code>alloc_local</code> return value, which
Chris@42: <em>may</em> be greater than <code>local_n0 * N1</code> (in case some
Chris@42: intermediate calculations require additional storage).  The data
Chris@42: distribution in FFTW&rsquo;s MPI interface is described in more detail by
Chris@42: the next section.
Chris@42: <a name="index-fftw_005fmpi_005flocal_005fsize_005f2d"></a>
Chris@42: <a name="index-data-distribution-1"></a>
Chris@42: </p>
Chris@42: 
Chris@42: <p>Given the portion of the array that resides on the local process, it
Chris@42: is straightforward to initialize the data (here to a function
Chris@42: <code>myfunction</code>) and otherwise manipulate it.  Of course, at the end
Chris@42: of the program you may want to output the data somehow, but
Chris@42: synchronizing this output is up to you and is beyond the scope of this
Chris@42: manual.  (One good way to output a large multi-dimensional distributed
Chris@42: array in MPI to a portable binary file is to use the free HDF5
Chris@42: library; see the <a href="http://www.hdfgroup.org/">HDF home page</a>.)
Chris@42: <a name="index-HDF5"></a>
Chris@42: <a name="index-MPI-I_002fO"></a>
Chris@42: </p>
Chris@42: <hr>
Chris@42: <div class="header">
Chris@42: <p>
Chris@42: Next: <a href="MPI-Data-Distribution.html#MPI-Data-Distribution" accesskey="n" rel="next">MPI Data Distribution</a>, Previous: <a href="Linking-and-Initializing-MPI-FFTW.html#Linking-and-Initializing-MPI-FFTW" accesskey="p" rel="prev">Linking and Initializing MPI FFTW</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> &nbsp; [<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: </div>
Chris@42: 
Chris@42: 
Chris@42: 
Chris@42: </body>
Chris@42: </html>