Mercurial > hg > sv-dependency-builds
view src/fftw-3.3.5/doc/html/Complex-One_002dDimensional-DFTs.html @ 42:2cd0e3b3e1fd
Current fftw source
author | Chris Cannam |
---|---|
date | Tue, 18 Oct 2016 13:40:26 +0100 |
parents | |
children |
line wrap: on
line source
<!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.5, 30 July 2016). 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 5.2, http://www.gnu.org/software/texinfo/ --> <head> <title>FFTW 3.3.5: Complex One-Dimensional DFTs</title> <meta name="description" content="FFTW 3.3.5: Complex One-Dimensional DFTs"> <meta name="keywords" content="FFTW 3.3.5: Complex One-Dimensional DFTs"> <meta name="resource-type" content="document"> <meta name="distribution" content="global"> <meta name="Generator" content="makeinfo"> <meta http-equiv="Content-Type" content="text/html; charset=utf-8"> <link href="index.html#Top" rel="start" title="Top"> <link href="Concept-Index.html#Concept-Index" rel="index" title="Concept Index"> <link href="index.html#SEC_Contents" rel="contents" title="Table of Contents"> <link href="Tutorial.html#Tutorial" rel="up" title="Tutorial"> <link href="Complex-Multi_002dDimensional-DFTs.html#Complex-Multi_002dDimensional-DFTs" rel="next" title="Complex Multi-Dimensional DFTs"> <link href="Tutorial.html#Tutorial" rel="prev" title="Tutorial"> <style type="text/css"> <!-- a.summary-letter {text-decoration: none} blockquote.smallquotation {font-size: smaller} div.display {margin-left: 3.2em} div.example {margin-left: 3.2em} div.indentedblock {margin-left: 3.2em} div.lisp {margin-left: 3.2em} div.smalldisplay {margin-left: 3.2em} div.smallexample {margin-left: 3.2em} div.smallindentedblock {margin-left: 3.2em; font-size: smaller} div.smalllisp {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} pre.smalldisplay {font-family: inherit; font-size: smaller} pre.smallexample {font-size: smaller} pre.smallformat {font-family: inherit; font-size: smaller} pre.smalllisp {font-size: smaller} span.nocodebreak {white-space:nowrap} span.nolinebreak {white-space:nowrap} span.roman {font-family:serif; font-weight:normal} span.sansserif {font-family:sans-serif; font-weight:normal} ul.no-bullet {list-style: none} --> </style> </head> <body lang="en" bgcolor="#FFFFFF" text="#000000" link="#0000FF" vlink="#800080" alink="#FF0000"> <a name="Complex-One_002dDimensional-DFTs"></a> <div class="header"> <p> Next: <a href="Complex-Multi_002dDimensional-DFTs.html#Complex-Multi_002dDimensional-DFTs" accesskey="n" rel="next">Complex Multi-Dimensional DFTs</a>, Previous: <a href="Tutorial.html#Tutorial" accesskey="p" rel="prev">Tutorial</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> </div> <hr> <a name="Complex-One_002dDimensional-DFTs-1"></a> <h3 class="section">2.1 Complex One-Dimensional DFTs</h3> <blockquote> <p>Plan: To bother about the best method of accomplishing an accidental result. [Ambrose Bierce, <cite>The Enlarged Devil’s Dictionary</cite>.] <a name="index-Devil"></a> </p></blockquote> <p>The basic usage of FFTW to compute a one-dimensional DFT of size <code>N</code> is simple, and it typically looks something like this code: </p> <div class="example"> <pre class="example">#include <fftw3.h> ... { fftw_complex *in, *out; fftw_plan p; ... in = (fftw_complex*) fftw_malloc(sizeof(fftw_complex) * N); out = (fftw_complex*) fftw_malloc(sizeof(fftw_complex) * N); p = fftw_plan_dft_1d(N, in, out, FFTW_FORWARD, FFTW_ESTIMATE); ... fftw_execute(p); /* <span class="roman">repeat as needed</span> */ ... fftw_destroy_plan(p); fftw_free(in); fftw_free(out); } </pre></div> <p>You must link this code with the <code>fftw3</code> library. On Unix systems, link with <code>-lfftw3 -lm</code>. </p> <p>The example code first allocates the input and output arrays. You can allocate them in any way that you like, but we recommend using <code>fftw_malloc</code>, which behaves like <a name="index-fftw_005fmalloc"></a> <code>malloc</code> except that it properly aligns the array when SIMD instructions (such as SSE and Altivec) are available (see <a href="SIMD-alignment-and-fftw_005fmalloc.html#SIMD-alignment-and-fftw_005fmalloc">SIMD alignment and fftw_malloc</a>). [Alternatively, we provide a convenient wrapper function <code>fftw_alloc_complex(N)</code> which has the same effect.] <a name="index-fftw_005falloc_005fcomplex"></a> <a name="index-SIMD"></a> </p> <p>The data is an array of type <code>fftw_complex</code>, which is by default a <code>double[2]</code> composed of the real (<code>in[i][0]</code>) and imaginary (<code>in[i][1]</code>) parts of a complex number. <a name="index-fftw_005fcomplex"></a> </p> <p>The next step is to create a <em>plan</em>, which is an object <a name="index-plan-1"></a> that contains all the data that FFTW needs to compute the FFT. This function creates the plan: </p> <div class="example"> <pre class="example">fftw_plan fftw_plan_dft_1d(int n, fftw_complex *in, fftw_complex *out, int sign, unsigned flags); </pre></div> <a name="index-fftw_005fplan_005fdft_005f1d"></a> <a name="index-fftw_005fplan"></a> <p>The first argument, <code>n</code>, is the size of the transform you are trying to compute. The size <code>n</code> can be any positive integer, but sizes that are products of small factors are transformed most efficiently (although prime sizes still use an <i>O</i>(<i>n</i> log <i>n</i>) algorithm). </p> <p>The next two arguments are pointers to the input and output arrays of the transform. These pointers can be equal, indicating an <em>in-place</em> transform. <a name="index-in_002dplace"></a> </p> <p>The fourth argument, <code>sign</code>, can be either <code>FFTW_FORWARD</code> (<code>-1</code>) or <code>FFTW_BACKWARD</code> (<code>+1</code>), <a name="index-FFTW_005fFORWARD"></a> <a name="index-FFTW_005fBACKWARD"></a> and indicates the direction of the transform you are interested in; technically, it is the sign of the exponent in the transform. </p> <p>The <code>flags</code> argument is usually either <code>FFTW_MEASURE</code> or <a name="index-flags"></a> <code>FFTW_ESTIMATE</code>. <code>FFTW_MEASURE</code> instructs FFTW to run <a name="index-FFTW_005fMEASURE"></a> and measure the execution time of several FFTs in order to find the best way to compute the transform of size <code>n</code>. This process takes some time (usually a few seconds), depending on your machine and on the size of the transform. <code>FFTW_ESTIMATE</code>, on the contrary, does not run any computation and just builds a <a name="index-FFTW_005fESTIMATE"></a> reasonable plan that is probably sub-optimal. In short, if your program performs many transforms of the same size and initialization time is not important, use <code>FFTW_MEASURE</code>; otherwise use the estimate. </p> <p><em>You must create the plan before initializing the input</em>, because <code>FFTW_MEASURE</code> overwrites the <code>in</code>/<code>out</code> arrays. (Technically, <code>FFTW_ESTIMATE</code> does not touch your arrays, but you should always create plans first just to be sure.) </p> <p>Once the plan has been created, you can use it as many times as you like for transforms on the specified <code>in</code>/<code>out</code> arrays, computing the actual transforms via <code>fftw_execute(plan)</code>: </p><div class="example"> <pre class="example">void fftw_execute(const fftw_plan plan); </pre></div> <a name="index-fftw_005fexecute"></a> <p>The DFT results are stored in-order in the array <code>out</code>, with the zero-frequency (DC) component in <code>out[0]</code>. <a name="index-frequency"></a> If <code>in != out</code>, the transform is <em>out-of-place</em> and the input array <code>in</code> is not modified. Otherwise, the input array is overwritten with the transform. </p> <a name="index-execute-1"></a> <p>If you want to transform a <em>different</em> array of the same size, you can create a new plan with <code>fftw_plan_dft_1d</code> and FFTW automatically reuses the information from the previous plan, if possible. Alternatively, with the “guru” interface you can apply a given plan to a different array, if you are careful. See <a href="FFTW-Reference.html#FFTW-Reference">FFTW Reference</a>. </p> <p>When you are done with the plan, you deallocate it by calling <code>fftw_destroy_plan(plan)</code>: </p><div class="example"> <pre class="example">void fftw_destroy_plan(fftw_plan plan); </pre></div> <a name="index-fftw_005fdestroy_005fplan"></a> <p>If you allocate an array with <code>fftw_malloc()</code> you must deallocate it with <code>fftw_free()</code>. Do not use <code>free()</code> or, heaven forbid, <code>delete</code>. <a name="index-fftw_005ffree"></a> </p> <p>FFTW computes an <em>unnormalized</em> DFT. Thus, computing a forward followed by a backward transform (or vice versa) results in the original array scaled by <code>n</code>. For the definition of the DFT, see <a href="What-FFTW-Really-Computes.html#What-FFTW-Really-Computes">What FFTW Really Computes</a>. <a name="index-DFT-1"></a> <a name="index-normalization"></a> </p> <p>If you have a C compiler, such as <code>gcc</code>, that supports the C99 standard, and you <code>#include <complex.h></code> <em>before</em> <code><fftw3.h></code>, then <code>fftw_complex</code> is the native double-precision complex type and you can manipulate it with ordinary arithmetic. Otherwise, FFTW defines its own complex type, which is bit-compatible with the C99 complex type. See <a href="Complex-numbers.html#Complex-numbers">Complex numbers</a>. (The C++ <code><complex></code> template class may also be usable via a typecast.) <a name="index-C_002b_002b"></a> </p> <p>To use single or long-double precision versions of FFTW, replace the <code>fftw_</code> prefix by <code>fftwf_</code> or <code>fftwl_</code> and link with <code>-lfftw3f</code> or <code>-lfftw3l</code>, but use the <em>same</em> <code><fftw3.h></code> header file. <a name="index-precision"></a> </p> <p>Many more flags exist besides <code>FFTW_MEASURE</code> and <code>FFTW_ESTIMATE</code>. For example, use <code>FFTW_PATIENT</code> if you’re willing to wait even longer for a possibly even faster plan (see <a href="FFTW-Reference.html#FFTW-Reference">FFTW Reference</a>). <a name="index-FFTW_005fPATIENT"></a> You can also save plans for future use, as described by <a href="Words-of-Wisdom_002dSaving-Plans.html#Words-of-Wisdom_002dSaving-Plans">Words of Wisdom-Saving Plans</a>. </p> <hr> <div class="header"> <p> Next: <a href="Complex-Multi_002dDimensional-DFTs.html#Complex-Multi_002dDimensional-DFTs" accesskey="n" rel="next">Complex Multi-Dimensional DFTs</a>, Previous: <a href="Tutorial.html#Tutorial" accesskey="p" rel="prev">Tutorial</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> </div> </body> </html>