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