Sophie

Sophie

distrib > Fedora > 18 > i386 > by-pkgid > f949c940bce372b03cb072dd0ee090a1 > files > 28

fftw-doc-3.3.3-5.fc18.noarch.rpm

<html lang="en">
<head>
<title>Complex Multi-Dimensional DFTs - FFTW 3.3.3</title>
<meta http-equiv="Content-Type" content="text/html">
<meta name="description" content="FFTW 3.3.3">
<meta name="generator" content="makeinfo 4.13">
<link title="Top" rel="start" href="index.html#Top">
<link rel="up" href="Tutorial.html#Tutorial" title="Tutorial">
<link rel="prev" href="Complex-One_002dDimensional-DFTs.html#Complex-One_002dDimensional-DFTs" title="Complex One-Dimensional DFTs">
<link rel="next" href="One_002dDimensional-DFTs-of-Real-Data.html#One_002dDimensional-DFTs-of-Real-Data" title="One-Dimensional DFTs of Real Data">
<link href="http://www.gnu.org/software/texinfo/" rel="generator-home" title="Texinfo Homepage">
<!--
This manual is for FFTW
(version 3.3.3, 25 November 2012).

Copyright (C) 2003 Matteo Frigo.

Copyright (C) 2003 Massachusetts Institute of Technology.

     Permission is granted to make and distribute verbatim copies of
     this manual provided the copyright notice and this permission
     notice are preserved on all copies.

     Permission is granted to copy and distribute modified versions of
     this manual under the conditions for verbatim copying, provided
     that the entire resulting derived work is distributed under the
     terms of a permission notice identical to this one.

     Permission is granted to copy and distribute translations of this
     manual into another language, under the above conditions for
     modified versions, except that this permission notice may be
     stated in a translation approved by the Free Software Foundation.
   -->
<meta http-equiv="Content-Style-Type" content="text/css">
<style type="text/css"><!--
  pre.display { font-family:inherit }
  pre.format  { font-family:inherit }
  pre.smalldisplay { font-family:inherit; font-size:smaller }
  pre.smallformat  { font-family:inherit; font-size:smaller }
  pre.smallexample { font-size:smaller }
  pre.smalllisp    { font-size:smaller }
  span.sc    { font-variant:small-caps }
  span.roman { font-family:serif; font-weight:normal; } 
  span.sansserif { font-family:sans-serif; font-weight:normal; } 
--></style>
</head>
<body>
<div class="node">
<a name="Complex-Multi-Dimensional-DFTs"></a>
<a name="Complex-Multi_002dDimensional-DFTs"></a>
<p>
Next:&nbsp;<a rel="next" accesskey="n" href="One_002dDimensional-DFTs-of-Real-Data.html#One_002dDimensional-DFTs-of-Real-Data">One-Dimensional DFTs of Real Data</a>,
Previous:&nbsp;<a rel="previous" accesskey="p" href="Complex-One_002dDimensional-DFTs.html#Complex-One_002dDimensional-DFTs">Complex One-Dimensional DFTs</a>,
Up:&nbsp;<a rel="up" accesskey="u" href="Tutorial.html#Tutorial">Tutorial</a>
<hr>
</div>

<h3 class="section">2.2 Complex Multi-Dimensional DFTs</h3>

<p>Multi-dimensional transforms work much the same way as one-dimensional
transforms: you allocate arrays of <code>fftw_complex</code> (preferably
using <code>fftw_malloc</code>), create an <code>fftw_plan</code>, execute it as
many times as you want with <code>fftw_execute(plan)</code>, and clean up
with <code>fftw_destroy_plan(plan)</code> (and <code>fftw_free</code>).

   <p>FFTW provides two routines for creating plans for 2d and 3d transforms,
and one routine for creating plans of arbitrary dimensionality. 
The 2d and 3d routines have the following signature:
<pre class="example">     fftw_plan fftw_plan_dft_2d(int n0, int n1,
                                fftw_complex *in, fftw_complex *out,
                                int sign, unsigned flags);
     fftw_plan fftw_plan_dft_3d(int n0, int n1, int n2,
                                fftw_complex *in, fftw_complex *out,
                                int sign, unsigned flags);
</pre>
   <p><a name="index-fftw_005fplan_005fdft_005f2d-39"></a><a name="index-fftw_005fplan_005fdft_005f3d-40"></a>
These routines create plans for <code>n0</code> by <code>n1</code> two-dimensional
(2d) transforms and <code>n0</code> by <code>n1</code> by <code>n2</code> 3d transforms,
respectively.  All of these transforms operate on contiguous arrays in
the C-standard <dfn>row-major</dfn> order, so that the last dimension has the
fastest-varying index in the array.  This layout is described further in
<a href="Multi_002ddimensional-Array-Format.html#Multi_002ddimensional-Array-Format">Multi-dimensional Array Format</a>.

   <p>FFTW can also compute transforms of higher dimensionality.  In order to
avoid confusion between the various meanings of the the word
&ldquo;dimension&rdquo;, we use the term <em>rank</em>
<a name="index-rank-41"></a>to denote the number of independent indices in an array.<a rel="footnote" href="#fn-1" name="fnd-1"><sup>1</sup></a>  For
example, we say that a 2d transform has rank&nbsp;2, a 3d transform has
rank&nbsp;3, and so on.  You can plan transforms of arbitrary rank by
means of the following function:

<pre class="example">     fftw_plan fftw_plan_dft(int rank, const int *n,
                             fftw_complex *in, fftw_complex *out,
                             int sign, unsigned flags);
</pre>
   <p><a name="index-fftw_005fplan_005fdft-42"></a>
Here, <code>n</code> is a pointer to an array <code>n[rank]</code> denoting an
<code>n[0]</code> by <code>n[1]</code> by <small class="dots">...</small> by <code>n[rank-1]</code> transform. 
Thus, for example, the call
<pre class="example">     fftw_plan_dft_2d(n0, n1, in, out, sign, flags);
</pre>
   <p>is equivalent to the following code fragment:
<pre class="example">     int n[2];
     n[0] = n0;
     n[1] = n1;
     fftw_plan_dft(2, n, in, out, sign, flags);
</pre>
   <p><code>fftw_plan_dft</code> is not restricted to 2d and 3d transforms,
however, but it can plan transforms of arbitrary rank.

   <p>You may have noticed that all the planner routines described so far
have overlapping functionality.  For example, you can plan a 1d or 2d
transform by using <code>fftw_plan_dft</code> with a <code>rank</code> of <code>1</code>
or <code>2</code>, or even by calling <code>fftw_plan_dft_3d</code> with <code>n0</code>
and/or <code>n1</code> equal to <code>1</code> (with no loss in efficiency).  This
pattern continues, and FFTW's planning routines in general form a
&ldquo;partial order,&rdquo; sequences of
<a name="index-partial-order-43"></a>interfaces with strictly increasing generality but correspondingly
greater complexity.

   <p><code>fftw_plan_dft</code> is the most general complex-DFT routine that we
describe in this tutorial, but there are also the advanced and guru interfaces,
<a name="index-advanced-interface-44"></a><a name="index-guru-interface-45"></a>which allow one to efficiently combine multiple/strided transforms
into a single FFTW plan, transform a subset of a larger
multi-dimensional array, and/or to handle more general complex-number
formats.  For more information, see <a href="FFTW-Reference.html#FFTW-Reference">FFTW Reference</a>.

<!--  -->
   <div class="footnote">
<hr>
<h4>Footnotes</h4><p class="footnote"><small>[<a name="fn-1" href="#fnd-1">1</a>]</small> The
term &ldquo;rank&rdquo; is commonly used in the APL, FORTRAN, and Common Lisp
traditions, although it is not so common in the C&nbsp;world.</p>

   <hr></div>

   </body></html>