Chris@82: <!DOCTYPE html PUBLIC "-//W3C//DTD HTML 4.01 Transitional//EN" "http://www.w3.org/TR/html4/loose.dtd">
Chris@82: <html>
Chris@82: <!-- This manual is for FFTW
Chris@82: (version 3.3.8, 24 May 2018).
Chris@82: 
Chris@82: Copyright (C) 2003 Matteo Frigo.
Chris@82: 
Chris@82: Copyright (C) 2003 Massachusetts Institute of Technology.
Chris@82: 
Chris@82: Permission is granted to make and distribute verbatim copies of this
Chris@82: manual provided the copyright notice and this permission notice are
Chris@82: preserved on all copies.
Chris@82: 
Chris@82: Permission is granted to copy and distribute modified versions of this
Chris@82: manual under the conditions for verbatim copying, provided that the
Chris@82: entire resulting derived work is distributed under the terms of a
Chris@82: permission notice identical to this one.
Chris@82: 
Chris@82: Permission is granted to copy and distribute translations of this manual
Chris@82: into another language, under the above conditions for modified versions,
Chris@82: except that this permission notice may be stated in a translation
Chris@82: approved by the Free Software Foundation. -->
Chris@82: <!-- Created by GNU Texinfo 6.3, http://www.gnu.org/software/texinfo/ -->
Chris@82: <head>
Chris@82: <title>FFTW 3.3.8: Multi-Dimensional DFTs of Real Data</title>
Chris@82: 
Chris@82: <meta name="description" content="FFTW 3.3.8: Multi-Dimensional DFTs of Real Data">
Chris@82: <meta name="keywords" content="FFTW 3.3.8: Multi-Dimensional DFTs of Real Data">
Chris@82: <meta name="resource-type" content="document">
Chris@82: <meta name="distribution" content="global">
Chris@82: <meta name="Generator" content="makeinfo">
Chris@82: <meta http-equiv="Content-Type" content="text/html; charset=utf-8">
Chris@82: <link href="index.html#Top" rel="start" title="Top">
Chris@82: <link href="Concept-Index.html#Concept-Index" rel="index" title="Concept Index">
Chris@82: <link href="index.html#SEC_Contents" rel="contents" title="Table of Contents">
Chris@82: <link href="Tutorial.html#Tutorial" rel="up" title="Tutorial">
Chris@82: <link href="More-DFTs-of-Real-Data.html#More-DFTs-of-Real-Data" rel="next" title="More DFTs of Real Data">
Chris@82: <link href="One_002dDimensional-DFTs-of-Real-Data.html#One_002dDimensional-DFTs-of-Real-Data" rel="prev" title="One-Dimensional DFTs of Real Data">
Chris@82: <style type="text/css">
Chris@82: <!--
Chris@82: a.summary-letter {text-decoration: none}
Chris@82: blockquote.indentedblock {margin-right: 0em}
Chris@82: blockquote.smallindentedblock {margin-right: 0em; font-size: smaller}
Chris@82: blockquote.smallquotation {font-size: smaller}
Chris@82: div.display {margin-left: 3.2em}
Chris@82: div.example {margin-left: 3.2em}
Chris@82: div.lisp {margin-left: 3.2em}
Chris@82: div.smalldisplay {margin-left: 3.2em}
Chris@82: div.smallexample {margin-left: 3.2em}
Chris@82: div.smalllisp {margin-left: 3.2em}
Chris@82: kbd {font-style: oblique}
Chris@82: pre.display {font-family: inherit}
Chris@82: pre.format {font-family: inherit}
Chris@82: pre.menu-comment {font-family: serif}
Chris@82: pre.menu-preformatted {font-family: serif}
Chris@82: pre.smalldisplay {font-family: inherit; font-size: smaller}
Chris@82: pre.smallexample {font-size: smaller}
Chris@82: pre.smallformat {font-family: inherit; font-size: smaller}
Chris@82: pre.smalllisp {font-size: smaller}
Chris@82: span.nolinebreak {white-space: nowrap}
Chris@82: span.roman {font-family: initial; font-weight: normal}
Chris@82: span.sansserif {font-family: sans-serif; font-weight: normal}
Chris@82: ul.no-bullet {list-style: none}
Chris@82: -->
Chris@82: </style>
Chris@82: 
Chris@82: 
Chris@82: </head>
Chris@82: 
Chris@82: <body lang="en">
Chris@82: <a name="Multi_002dDimensional-DFTs-of-Real-Data"></a>
Chris@82: <div class="header">
Chris@82: <p>
Chris@82: Next: <a href="More-DFTs-of-Real-Data.html#More-DFTs-of-Real-Data" accesskey="n" rel="next">More DFTs of Real Data</a>, Previous: <a href="One_002dDimensional-DFTs-of-Real-Data.html#One_002dDimensional-DFTs-of-Real-Data" accesskey="p" rel="prev">One-Dimensional DFTs of Real Data</a>, Up: <a href="Tutorial.html#Tutorial" accesskey="u" rel="up">Tutorial</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@82: </div>
Chris@82: <hr>
Chris@82: <a name="Multi_002dDimensional-DFTs-of-Real-Data-1"></a>
Chris@82: <h3 class="section">2.4 Multi-Dimensional DFTs of Real Data</h3>
Chris@82: 
Chris@82: <p>Multi-dimensional DFTs of real data use the following planner routines:
Chris@82: </p>
Chris@82: <div class="example">
Chris@82: <pre class="example">fftw_plan fftw_plan_dft_r2c_2d(int n0, int n1,
Chris@82:                                double *in, fftw_complex *out,
Chris@82:                                unsigned flags);
Chris@82: fftw_plan fftw_plan_dft_r2c_3d(int n0, int n1, int n2,
Chris@82:                                double *in, fftw_complex *out,
Chris@82:                                unsigned flags);
Chris@82: fftw_plan fftw_plan_dft_r2c(int rank, const int *n,
Chris@82:                             double *in, fftw_complex *out,
Chris@82:                             unsigned flags);
Chris@82: </pre></div>
Chris@82: <a name="index-fftw_005fplan_005fdft_005fr2c_005f2d"></a>
Chris@82: <a name="index-fftw_005fplan_005fdft_005fr2c_005f3d"></a>
Chris@82: <a name="index-fftw_005fplan_005fdft_005fr2c"></a>
Chris@82: 
Chris@82: <p>as well as the corresponding <code>c2r</code> routines with the input/output
Chris@82: types swapped.  These routines work similarly to their complex
Chris@82: analogues, except for the fact that here the complex output array is cut
Chris@82: roughly in half and the real array requires padding for in-place
Chris@82: transforms (as in 1d, above).
Chris@82: </p>
Chris@82: <p>As before, <code>n</code> is the logical size of the array, and the
Chris@82: consequences of this on the the format of the complex arrays deserve
Chris@82: careful attention.
Chris@82: <a name="index-r2c_002fc2r-multi_002ddimensional-array-format"></a>
Chris@82: Suppose that the real data has dimensions n<sub>0</sub>&nbsp;&times;&nbsp;n<sub>1</sub>&nbsp;&times;&nbsp;n<sub>2</sub>&nbsp;&times;&nbsp;&hellip;&nbsp;&times;&nbsp;n<sub>d-1</sub>
Chris@82:  (in row-major order).
Chris@82: Then, after an r2c transform, the output is an n<sub>0</sub>&nbsp;&times;&nbsp;n<sub>1</sub>&nbsp;&times;&nbsp;n<sub>2</sub>&nbsp;&times;&nbsp;&hellip;&nbsp;&times;&nbsp;(n<sub>d-1</sub>/2 + 1)
Chris@82:  array of
Chris@82: <code>fftw_complex</code> values in row-major order, corresponding to slightly
Chris@82: over half of the output of the corresponding complex DFT.  (The division
Chris@82: is rounded down.)  The ordering of the data is otherwise exactly the
Chris@82: same as in the complex-DFT case.
Chris@82: </p>
Chris@82: <p>For out-of-place transforms, this is the end of the story: the real
Chris@82: data is stored as a row-major array of size n<sub>0</sub>&nbsp;&times;&nbsp;n<sub>1</sub>&nbsp;&times;&nbsp;n<sub>2</sub>&nbsp;&times;&nbsp;&hellip;&nbsp;&times;&nbsp;n<sub>d-1</sub>
Chris@82:  and the complex
Chris@82: data is stored as a row-major array of size n<sub>0</sub>&nbsp;&times;&nbsp;n<sub>1</sub>&nbsp;&times;&nbsp;n<sub>2</sub>&nbsp;&times;&nbsp;&hellip;&nbsp;&times;&nbsp;(n<sub>d-1</sub>/2 + 1)
Chris@82: .
Chris@82: </p>
Chris@82: <p>For in-place transforms, however, extra padding of the real-data array
Chris@82: is necessary because the complex array is larger than the real array,
Chris@82: and the two arrays share the same memory locations.  Thus, for
Chris@82: in-place transforms, the final dimension of the real-data array must
Chris@82: be padded with extra values to accommodate the size of the complex
Chris@82: data&mdash;two values if the last dimension is even and one if it is odd.
Chris@82: <a name="index-padding-1"></a>
Chris@82: That is, the last dimension of the real data must physically contain
Chris@82: 2 * (n<sub>d-1</sub>/2+1)
Chris@82: <code>double</code> values (exactly enough to hold the complex data).
Chris@82: This physical array size does not, however, change the <em>logical</em>
Chris@82: array size&mdash;only
Chris@82: n<sub>d-1</sub>
Chris@82: values are actually stored in the last dimension, and
Chris@82: n<sub>d-1</sub>
Chris@82: is the last dimension passed to the plan-creation routine.
Chris@82: </p>
Chris@82: <p>For example, consider the transform of a two-dimensional real array of
Chris@82: size <code>n0</code> by <code>n1</code>.  The output of the r2c transform is a
Chris@82: two-dimensional complex array of size <code>n0</code> by <code>n1/2+1</code>, where
Chris@82: the <code>y</code> dimension has been cut nearly in half because of
Chris@82: redundancies in the output.  Because <code>fftw_complex</code> is twice the
Chris@82: size of <code>double</code>, the output array is slightly bigger than the
Chris@82: input array.  Thus, if we want to compute the transform in place, we
Chris@82: must <em>pad</em> the input array so that it is of size <code>n0</code> by
Chris@82: <code>2*(n1/2+1)</code>.  If <code>n1</code> is even, then there are two padding
Chris@82: elements at the end of each row (which need not be initialized, as they
Chris@82: are only used for output).
Chris@82: </p>
Chris@82: <p>The following illustration depicts the input and output arrays just
Chris@82: described, for both the out-of-place and in-place transforms (with the
Chris@82: arrows indicating consecutive memory locations):
Chris@82: <img src="rfftwnd-for-html.png" alt="rfftwnd-for-html">
Chris@82: </p>
Chris@82: <p>These transforms are unnormalized, so an r2c followed by a c2r
Chris@82: transform (or vice versa) will result in the original data scaled by
Chris@82: the number of real data elements&mdash;that is, the product of the
Chris@82: (logical) dimensions of the real data.
Chris@82: <a name="index-normalization-1"></a>
Chris@82: </p>
Chris@82: 
Chris@82: <p>(Because the last dimension is treated specially, if it is equal to
Chris@82: <code>1</code> the transform is <em>not</em> equivalent to a lower-dimensional
Chris@82: r2c/c2r transform.  In that case, the last complex dimension also has
Chris@82: size <code>1</code> (<code>=1/2+1</code>), and no advantage is gained over the
Chris@82: complex transforms.)
Chris@82: </p>
Chris@82: <hr>
Chris@82: <div class="header">
Chris@82: <p>
Chris@82: Next: <a href="More-DFTs-of-Real-Data.html#More-DFTs-of-Real-Data" accesskey="n" rel="next">More DFTs of Real Data</a>, Previous: <a href="One_002dDimensional-DFTs-of-Real-Data.html#One_002dDimensional-DFTs-of-Real-Data" accesskey="p" rel="prev">One-Dimensional DFTs of Real Data</a>, Up: <a href="Tutorial.html#Tutorial" accesskey="u" rel="up">Tutorial</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@82: </div>
Chris@82: 
Chris@82: 
Chris@82: 
Chris@82: </body>
Chris@82: </html>