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: One-Dimensional DFTs of Real Data</title> Chris@82: Chris@82: <meta name="description" content="FFTW 3.3.8: One-Dimensional DFTs of Real Data"> Chris@82: <meta name="keywords" content="FFTW 3.3.8: One-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="Multi_002dDimensional-DFTs-of-Real-Data.html#Multi_002dDimensional-DFTs-of-Real-Data" rel="next" title="Multi-Dimensional DFTs of Real Data"> Chris@82: <link href="Complex-Multi_002dDimensional-DFTs.html#Complex-Multi_002dDimensional-DFTs" rel="prev" title="Complex Multi-Dimensional DFTs"> 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="One_002dDimensional-DFTs-of-Real-Data"></a> Chris@82: <div class="header"> Chris@82: <p> Chris@82: Next: <a href="Multi_002dDimensional-DFTs-of-Real-Data.html#Multi_002dDimensional-DFTs-of-Real-Data" accesskey="n" rel="next">Multi-Dimensional DFTs of Real Data</a>, Previous: <a href="Complex-Multi_002dDimensional-DFTs.html#Complex-Multi_002dDimensional-DFTs" accesskey="p" rel="prev">Complex Multi-Dimensional DFTs</a>, Up: <a href="Tutorial.html#Tutorial" accesskey="u" rel="up">Tutorial</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: </div> Chris@82: <hr> Chris@82: <a name="One_002dDimensional-DFTs-of-Real-Data-1"></a> Chris@82: <h3 class="section">2.3 One-Dimensional DFTs of Real Data</h3> Chris@82: Chris@82: <p>In many practical applications, the input data <code>in[i]</code> are purely Chris@82: real numbers, in which case the DFT output satisfies the “Hermitian” Chris@82: <a name="index-Hermitian"></a> Chris@82: redundancy: <code>out[i]</code> is the conjugate of <code>out[n-i]</code>. It is Chris@82: possible to take advantage of these circumstances in order to achieve Chris@82: roughly a factor of two improvement in both speed and memory usage. Chris@82: </p> Chris@82: <p>In exchange for these speed and space advantages, the user sacrifices Chris@82: some of the simplicity of FFTW’s complex transforms. First of all, the Chris@82: input and output arrays are of <em>different sizes and types</em>: the Chris@82: input is <code>n</code> real numbers, while the output is <code>n/2+1</code> Chris@82: complex numbers (the non-redundant outputs); this also requires slight Chris@82: “padding” of the input array for Chris@82: <a name="index-padding"></a> Chris@82: in-place transforms. Second, the inverse transform (complex to real) Chris@82: has the side-effect of <em>overwriting its input array</em>, by default. Chris@82: Neither of these inconveniences should pose a serious problem for Chris@82: users, but it is important to be aware of them. Chris@82: </p> Chris@82: <p>The routines to perform real-data transforms are almost the same as Chris@82: those for complex transforms: you allocate arrays of <code>double</code> Chris@82: and/or <code>fftw_complex</code> (preferably using <code>fftw_malloc</code> or Chris@82: <code>fftw_alloc_complex</code>), create an <code>fftw_plan</code>, execute it as Chris@82: many times as you want with <code>fftw_execute(plan)</code>, and clean up Chris@82: with <code>fftw_destroy_plan(plan)</code> (and <code>fftw_free</code>). The only Chris@82: differences are that the input (or output) is of type <code>double</code> Chris@82: and there are new routines to create the plan. In one dimension: Chris@82: </p> Chris@82: <div class="example"> Chris@82: <pre class="example">fftw_plan fftw_plan_dft_r2c_1d(int n, double *in, fftw_complex *out, Chris@82: unsigned flags); Chris@82: fftw_plan fftw_plan_dft_c2r_1d(int n, fftw_complex *in, double *out, Chris@82: unsigned flags); Chris@82: </pre></div> Chris@82: <a name="index-fftw_005fplan_005fdft_005fr2c_005f1d"></a> Chris@82: <a name="index-fftw_005fplan_005fdft_005fc2r_005f1d"></a> Chris@82: Chris@82: <p>for the real input to complex-Hermitian output (<em>r2c</em>) and Chris@82: complex-Hermitian input to real output (<em>c2r</em>) transforms. Chris@82: <a name="index-r2c"></a> Chris@82: <a name="index-c2r"></a> Chris@82: Unlike the complex DFT planner, there is no <code>sign</code> argument. Chris@82: Instead, r2c DFTs are always <code>FFTW_FORWARD</code> and c2r DFTs are Chris@82: always <code>FFTW_BACKWARD</code>. Chris@82: <a name="index-FFTW_005fFORWARD-1"></a> Chris@82: <a name="index-FFTW_005fBACKWARD-1"></a> Chris@82: (For single/long-double precision Chris@82: <code>fftwf</code> and <code>fftwl</code>, <code>double</code> should be replaced by Chris@82: <code>float</code> and <code>long double</code>, respectively.) Chris@82: <a name="index-precision-1"></a> Chris@82: </p> Chris@82: Chris@82: <p>Here, <code>n</code> is the “logical” size of the DFT, not necessarily the Chris@82: physical size of the array. In particular, the real (<code>double</code>) Chris@82: array has <code>n</code> elements, while the complex (<code>fftw_complex</code>) Chris@82: array has <code>n/2+1</code> elements (where the division is rounded down). Chris@82: For an in-place transform, Chris@82: <a name="index-in_002dplace-1"></a> Chris@82: <code>in</code> and <code>out</code> are aliased to the same array, which must be Chris@82: big enough to hold both; so, the real array would actually have Chris@82: <code>2*(n/2+1)</code> elements, where the elements beyond the first Chris@82: <code>n</code> are unused padding. (Note that this is very different from Chris@82: the concept of “zero-padding” a transform to a larger length, which Chris@82: changes the logical size of the DFT by actually adding new input Chris@82: data.) The <em>k</em>th element of the complex array is exactly the Chris@82: same as the <em>k</em>th element of the corresponding complex DFT. All Chris@82: positive <code>n</code> are supported; products of small factors are most Chris@82: efficient, but an <i>O</i>(<i>n</i> log <i>n</i>) Chris@82: algorithm is used even for prime sizes. Chris@82: </p> Chris@82: <p>As noted above, the c2r transform destroys its input array even for Chris@82: out-of-place transforms. This can be prevented, if necessary, by Chris@82: including <code>FFTW_PRESERVE_INPUT</code> in the <code>flags</code>, with Chris@82: unfortunately some sacrifice in performance. Chris@82: <a name="index-flags-1"></a> Chris@82: <a name="index-FFTW_005fPRESERVE_005fINPUT"></a> Chris@82: This flag is also not currently supported for multi-dimensional real Chris@82: DFTs (next section). Chris@82: </p> Chris@82: <p>Readers familiar with DFTs of real data will recall that the 0th (the Chris@82: “DC”) and <code>n/2</code>-th (the “Nyquist” frequency, when <code>n</code> is Chris@82: even) elements of the complex output are purely real. Some Chris@82: implementations therefore store the Nyquist element where the DC Chris@82: imaginary part would go, in order to make the input and output arrays Chris@82: the same size. Such packing, however, does not generalize well to Chris@82: multi-dimensional transforms, and the space savings are miniscule in Chris@82: any case; FFTW does not support it. Chris@82: </p> Chris@82: <p>An alternative interface for one-dimensional r2c and c2r DFTs can be Chris@82: found in the ‘<samp>r2r</samp>’ interface (see <a href="The-Halfcomplex_002dformat-DFT.html#The-Halfcomplex_002dformat-DFT">The Halfcomplex-format DFT</a>), with “halfcomplex”-format output that <em>is</em> the same size Chris@82: (and type) as the input array. Chris@82: <a name="index-halfcomplex-format"></a> Chris@82: That interface, although it is not very useful for multi-dimensional Chris@82: transforms, may sometimes yield better performance. Chris@82: </p> Chris@82: <hr> Chris@82: <div class="header"> Chris@82: <p> Chris@82: Next: <a href="Multi_002dDimensional-DFTs-of-Real-Data.html#Multi_002dDimensional-DFTs-of-Real-Data" accesskey="n" rel="next">Multi-Dimensional DFTs of Real Data</a>, Previous: <a href="Complex-Multi_002dDimensional-DFTs.html#Complex-Multi_002dDimensional-DFTs" accesskey="p" rel="prev">Complex Multi-Dimensional DFTs</a>, Up: <a href="Tutorial.html#Tutorial" accesskey="u" rel="up">Tutorial</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: </div> Chris@82: Chris@82: Chris@82: Chris@82: </body> Chris@82: </html>