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