3
<title>One-Dimensional DFTs of Real Data - FFTW 3.1.1</title>
4
<meta http-equiv="Content-Type" content="text/html">
5
<meta name="description" content="FFTW 3.1.1">
6
<meta name="generator" content="makeinfo 4.8">
7
<link title="Top" rel="start" href="index.html#Top">
8
<link rel="up" href="Tutorial.html#Tutorial" title="Tutorial">
9
<link rel="prev" href="Complex-Multi_002dDimensional-DFTs.html#Complex-Multi_002dDimensional-DFTs" title="Complex Multi-Dimensional DFTs">
10
<link rel="next" href="Multi_002dDimensional-DFTs-of-Real-Data.html#Multi_002dDimensional-DFTs-of-Real-Data" title="Multi-Dimensional DFTs of Real Data">
11
<link href="http://www.gnu.org/software/texinfo/" rel="generator-home" title="Texinfo Homepage">
13
This manual is for FFTW
14
(version 3.1.1, 4 March 2006).
16
Copyright (C) 2003 Matteo Frigo.
18
Copyright (C) 2003 Massachusetts Institute of Technology.
20
Permission is granted to make and distribute verbatim copies of
21
this manual provided the copyright notice and this permission
22
notice are preserved on all copies.
24
Permission is granted to copy and distribute modified versions of
25
this manual under the conditions for verbatim copying, provided
26
that the entire resulting derived work is distributed under the
27
terms of a permission notice identical to this one.
29
Permission is granted to copy and distribute translations of this
30
manual into another language, under the above conditions for
31
modified versions, except that this permission notice may be
32
stated in a translation approved by the Free Software Foundation.
34
<meta http-equiv="Content-Style-Type" content="text/css">
35
<style type="text/css"><!--
36
pre.display { font-family:inherit }
37
pre.format { font-family:inherit }
38
pre.smalldisplay { font-family:inherit; font-size:smaller }
39
pre.smallformat { font-family:inherit; font-size:smaller }
40
pre.smallexample { font-size:smaller }
41
pre.smalllisp { font-size:smaller }
42
span.sc { font-variant:small-caps }
43
span.roman { font-family:serif; font-weight:normal; }
44
span.sansserif { font-family:sans-serif; font-weight:normal; }
50
<a name="One-Dimensional-DFTs-of-Real-Data"></a>
51
<a name="One_002dDimensional-DFTs-of-Real-Data"></a>
52
Next: <a rel="next" accesskey="n" href="Multi_002dDimensional-DFTs-of-Real-Data.html#Multi_002dDimensional-DFTs-of-Real-Data">Multi-Dimensional DFTs of Real Data</a>,
53
Previous: <a rel="previous" accesskey="p" href="Complex-Multi_002dDimensional-DFTs.html#Complex-Multi_002dDimensional-DFTs">Complex Multi-Dimensional DFTs</a>,
54
Up: <a rel="up" accesskey="u" href="Tutorial.html#Tutorial">Tutorial</a>
58
<h3 class="section">2.3 One-Dimensional DFTs of Real Data</h3>
60
<p>In many practical applications, the input data <code>in[i]</code> are purely
61
real numbers, in which case the DFT output satisfies the “Hermitian”
62
<a name="index-Hermitian-45"></a>redundancy: <code>out[i]</code> is the conjugate of <code>out[n-i]</code>. It is
63
possible to take advantage of these circumstances in order to achieve
64
roughly a factor of two improvement in both speed and memory usage.
66
<p>In exchange for these speed and space advantages, the user sacrifices
67
some of the simplicity of FFTW's complex transforms. First of all, the
68
input and output arrays are of <em>different sizes and types</em>: the
69
input is <code>n</code> real numbers, while the output is <code>n/2+1</code>
70
complex numbers (the non-redundant outputs); this also requires slight
71
“padding” of the input array for
72
<a name="index-padding-46"></a>in-place transforms. Second, the inverse transform (complex to real)
73
has the side-effect of <em>destroying its input array</em>, by default.
74
Neither of these inconveniences should pose a serious problem for
75
users, but it is important to be aware of them.
77
<p>The routines to perform real-data transforms are almost the same as
78
those for complex transforms: you allocate arrays of <code>double</code>
79
and/or <code>fftw_complex</code> (preferably using <code>fftw_malloc</code>),
80
create an <code>fftw_plan</code>, execute it as many times as you want with
81
<code>fftw_execute(plan)</code>, and clean up with
82
<code>fftw_destroy_plan(plan)</code> (and <code>fftw_free</code>). The only
83
differences are that the input (or output) is of type <code>double</code>
84
and there are new routines to create the plan. In one dimension:
86
<pre class="example"> fftw_plan fftw_plan_dft_r2c_1d(int n, double *in, fftw_complex *out,
88
fftw_plan fftw_plan_dft_c2r_1d(int n, fftw_complex *in, double *out,
91
<p><a name="index-fftw_005fplan_005fdft_005fr2c_005f1d-47"></a><a name="index-fftw_005fplan_005fdft_005fc2r_005f1d-48"></a>
92
for the real input to complex-Hermitian output (<dfn>r2c</dfn>) and
93
complex-Hermitian input to real output (<dfn>c2r</dfn>) transforms.
94
<a name="index-r2c-49"></a><a name="index-c2r-50"></a>Unlike the complex DFT planner, there is no <code>sign</code> argument.
95
Instead, r2c DFTs are always <code>FFTW_FORWARD</code> and c2r DFTs are
96
always <code>FFTW_BACKWARD</code>.
97
<a name="index-FFTW_005fFORWARD-51"></a><a name="index-FFTW_005fBACKWARD-52"></a>(For single/long-double precision
98
<code>fftwf</code> and <code>fftwl</code>, <code>double</code> should be replaced by
99
<code>float</code> and <code>long double</code>, respectively.)
100
<a name="index-precision-53"></a>
101
Here, <code>n</code> is the “logical” size of the DFT, not necessarily the
102
physical size of the array. In particular, the real (<code>double</code>)
103
array has <code>n</code> elements, while the complex (<code>fftw_complex</code>)
104
array has <code>n/2+1</code> elements (where the division is rounded down).
105
For an in-place transform,
106
<a name="index-in_002dplace-54"></a><code>in</code> and <code>out</code> are aliased to the same array, which must be
107
big enough to hold both; so, the real array would actually have
108
<code>2*(n/2+1)</code> elements, where the elements beyond the first <code>n</code>
109
are unused padding. The kth element of the complex array is
110
exactly the same as the kth element of the corresponding complex
111
DFT. All positive <code>n</code> are supported; products of small factors are
112
most efficient, but an <i>O</i>(<i>n</i> log <i>n</i>) algorithm is used even for prime
115
<p>As noted above, the c2r transform destroys its input array even for
116
out-of-place transforms. This can be prevented, if necessary, by
117
including <code>FFTW_PRESERVE_INPUT</code> in the <code>flags</code>, with
118
unfortunately some sacrifice in performance.
119
<a name="index-flags-55"></a><a name="index-FFTW_005fPRESERVE_005fINPUT-56"></a>This flag is also not currently supported for multi-dimensional real
122
<p>Readers familiar with DFTs of real data will recall that the 0th (the
123
“DC”) and <code>n/2</code>-th (the “Nyquist” frequency, when <code>n</code> is
124
even) elements of the complex output are purely real. Some
125
implementations therefore store the Nyquist element where the DC
126
imaginary part would go, in order to make the input and output arrays
127
the same size. Such packing, however, does not generalize well to
128
multi-dimensional transforms, and the space savings are miniscule in
129
any case; FFTW does not support it.
131
<p>An alternative interface for one-dimensional r2c and c2r DFTs can be
132
found in the `<samp><span class="samp">r2r</span></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
133
(and type) as the input array.
134
<a name="index-halfcomplex-format-57"></a>That interface, although it is not very useful for multi-dimensional
135
transforms, may sometimes yield better performance.