furnace/extern/fftw/doc/html/The-Halfcomplex_002dformat-DFT.html

137 lines
6.9 KiB
HTML
Raw Normal View History

<!DOCTYPE html PUBLIC "-//W3C//DTD HTML 4.01 Transitional//EN" "http://www.w3.org/TR/html4/loose.dtd">
<html>
<!-- This manual is for FFTW
(version 3.3.10, 10 December 2020).
Copyright (C) 2003 Matteo Frigo.
Copyright (C) 2003 Massachusetts Institute of Technology.
Permission is granted to make and distribute verbatim copies of this
manual provided the copyright notice and this permission notice are
preserved on all copies.
Permission is granted to copy and distribute modified versions of this
manual under the conditions for verbatim copying, provided that the
entire resulting derived work is distributed under the terms of a
permission notice identical to this one.
Permission is granted to copy and distribute translations of this manual
into another language, under the above conditions for modified versions,
except that this permission notice may be stated in a translation
approved by the Free Software Foundation. -->
<!-- Created by GNU Texinfo 6.7, http://www.gnu.org/software/texinfo/ -->
<head>
<meta http-equiv="Content-Type" content="text/html; charset=utf-8">
<title>The Halfcomplex-format DFT (FFTW 3.3.10)</title>
<meta name="description" content="The Halfcomplex-format DFT (FFTW 3.3.10)">
<meta name="keywords" content="The Halfcomplex-format DFT (FFTW 3.3.10)">
<meta name="resource-type" content="document">
<meta name="distribution" content="global">
<meta name="Generator" content="makeinfo">
<link href="index.html" rel="start" title="Top">
<link href="Concept-Index.html" rel="index" title="Concept Index">
<link href="index.html#SEC_Contents" rel="contents" title="Table of Contents">
<link href="More-DFTs-of-Real-Data.html" rel="up" title="More DFTs of Real Data">
<link href="Real-even_002fodd-DFTs-_0028cosine_002fsine-transforms_0029.html" rel="next" title="Real even/odd DFTs (cosine/sine transforms)">
<link href="More-DFTs-of-Real-Data.html" rel="prev" title="More DFTs of Real Data">
<style type="text/css">
<!--
a.summary-letter {text-decoration: none}
blockquote.indentedblock {margin-right: 0em}
div.display {margin-left: 3.2em}
div.example {margin-left: 3.2em}
div.lisp {margin-left: 3.2em}
kbd {font-style: oblique}
pre.display {font-family: inherit}
pre.format {font-family: inherit}
pre.menu-comment {font-family: serif}
pre.menu-preformatted {font-family: serif}
span.nolinebreak {white-space: nowrap}
span.roman {font-family: initial; font-weight: normal}
span.sansserif {font-family: sans-serif; font-weight: normal}
ul.no-bullet {list-style: none}
-->
</style>
</head>
<body lang="en">
<span id="The-Halfcomplex_002dformat-DFT"></span><div class="header">
<p>
Next: <a href="Real-even_002fodd-DFTs-_0028cosine_002fsine-transforms_0029.html" accesskey="n" rel="next">Real even/odd DFTs (cosine/sine transforms)</a>, Previous: <a href="More-DFTs-of-Real-Data.html" accesskey="p" rel="prev">More DFTs of Real Data</a>, Up: <a href="More-DFTs-of-Real-Data.html" accesskey="u" rel="up">More DFTs of Real Data</a> &nbsp; [<a href="index.html#SEC_Contents" title="Table of contents" rel="contents">Contents</a>][<a href="Concept-Index.html" title="Index" rel="index">Index</a>]</p>
</div>
<hr>
<span id="The-Halfcomplex_002dformat-DFT-1"></span><h4 class="subsection">2.5.1 The Halfcomplex-format DFT</h4>
<p>An r2r kind of <code>FFTW_R2HC</code> (<em>r2hc</em>) corresponds to an r2c DFT
<span id="index-FFTW_005fR2HC"></span>
<span id="index-r2c-1"></span>
<span id="index-r2hc"></span>
(see <a href="One_002dDimensional-DFTs-of-Real-Data.html">One-Dimensional DFTs of Real Data</a>) but with &ldquo;halfcomplex&rdquo;
format output, and may sometimes be faster and/or more convenient than
the latter.
<span id="index-halfcomplex-format-1"></span>
The inverse <em>hc2r</em> transform is of kind <code>FFTW_HC2R</code>.
<span id="index-FFTW_005fHC2R"></span>
<span id="index-hc2r"></span>
This consists of the non-redundant half of the complex output for a 1d
real-input DFT of size <code>n</code>, stored as a sequence of <code>n</code> real
numbers (<code>double</code>) in the format:
</p>
<p align=center>
r<sub>0</sub>, r<sub>1</sub>, r<sub>2</sub>, ..., r<sub>n/2</sub>, i<sub>(n+1)/2-1</sub>, ..., i<sub>2</sub>, i<sub>1</sub>
</p>
<p>Here,
r<sub>k</sub>
is the real part of the <em>k</em>th output, and
i<sub>k</sub>
is the imaginary part. (Division by 2 is rounded down.) For a
halfcomplex array <code>hc[n]</code>, the <em>k</em>th component thus has its
real part in <code>hc[k]</code> and its imaginary part in <code>hc[n-k]</code>, with
the exception of <code>k</code> <code>==</code> <code>0</code> or <code>n/2</code> (the latter
only if <code>n</code> is even)&mdash;in these two cases, the imaginary part is
zero due to symmetries of the real-input DFT, and is not stored.
Thus, the r2hc transform of <code>n</code> real values is a halfcomplex array of
length <code>n</code>, and vice versa for hc2r.
<span id="index-normalization-2"></span>
</p>
<p>Aside from the differing format, the output of
<code>FFTW_R2HC</code>/<code>FFTW_HC2R</code> is otherwise exactly the same as for
the corresponding 1d r2c/c2r transform
(i.e. <code>FFTW_FORWARD</code>/<code>FFTW_BACKWARD</code> transforms, respectively).
Recall that these transforms are unnormalized, so r2hc followed by hc2r
will result in the original data multiplied by <code>n</code>. Furthermore,
like the c2r transform, an out-of-place hc2r transform will
<em>destroy its input</em> array.
</p>
<p>Although these halfcomplex transforms can be used with the
multi-dimensional r2r interface, the interpretation of such a separable
product of transforms along each dimension is problematic. For example,
consider a two-dimensional <code>n0</code> by <code>n1</code>, r2hc by r2hc
transform planned by <code>fftw_plan_r2r_2d(n0, n1, in, out, FFTW_R2HC,
FFTW_R2HC, FFTW_MEASURE)</code>. Conceptually, FFTW first transforms the rows
(of size <code>n1</code>) to produce halfcomplex rows, and then transforms the
columns (of size <code>n0</code>). Half of these column transforms, however,
are of imaginary parts, and should therefore be multiplied by <em>i</em>
and combined with the r2hc transforms of the real columns to produce the
2d DFT amplitudes; FFTW&rsquo;s r2r transform does <em>not</em> perform this
combination for you. Thus, if a multi-dimensional real-input/output DFT
is required, we recommend using the ordinary r2c/c2r
interface (see <a href="Multi_002dDimensional-DFTs-of-Real-Data.html">Multi-Dimensional DFTs of Real Data</a>).
</p>
<hr>
<div class="header">
<p>
Next: <a href="Real-even_002fodd-DFTs-_0028cosine_002fsine-transforms_0029.html" accesskey="n" rel="next">Real even/odd DFTs (cosine/sine transforms)</a>, Previous: <a href="More-DFTs-of-Real-Data.html" accesskey="p" rel="prev">More DFTs of Real Data</a>, Up: <a href="More-DFTs-of-Real-Data.html" accesskey="u" rel="up">More DFTs of Real Data</a> &nbsp; [<a href="index.html#SEC_Contents" title="Table of contents" rel="contents">Contents</a>][<a href="Concept-Index.html" title="Index" rel="index">Index</a>]</p>
</div>
</body>
</html>