annotate fft/fftw/fftw-3.3.4/doc/html/New_002darray-Execute-Functions.html @ 40:223f770b5341 kissfft-double tip

Try a double-precision kissfft
author Chris Cannam
date Wed, 07 Sep 2016 10:40:32 +0100
parents 26056e866c29
children
rev   line source
Chris@19 1 <html lang="en">
Chris@19 2 <head>
Chris@19 3 <title>New-array Execute Functions - FFTW 3.3.4</title>
Chris@19 4 <meta http-equiv="Content-Type" content="text/html">
Chris@19 5 <meta name="description" content="FFTW 3.3.4">
Chris@19 6 <meta name="generator" content="makeinfo 4.13">
Chris@19 7 <link title="Top" rel="start" href="index.html#Top">
Chris@19 8 <link rel="up" href="FFTW-Reference.html#FFTW-Reference" title="FFTW Reference">
Chris@19 9 <link rel="prev" href="Guru-Interface.html#Guru-Interface" title="Guru Interface">
Chris@19 10 <link rel="next" href="Wisdom.html#Wisdom" title="Wisdom">
Chris@19 11 <link href="http://www.gnu.org/software/texinfo/" rel="generator-home" title="Texinfo Homepage">
Chris@19 12 <!--
Chris@19 13 This manual is for FFTW
Chris@19 14 (version 3.3.4, 20 September 2013).
Chris@19 15
Chris@19 16 Copyright (C) 2003 Matteo Frigo.
Chris@19 17
Chris@19 18 Copyright (C) 2003 Massachusetts Institute of Technology.
Chris@19 19
Chris@19 20 Permission is granted to make and distribute verbatim copies of
Chris@19 21 this manual provided the copyright notice and this permission
Chris@19 22 notice are preserved on all copies.
Chris@19 23
Chris@19 24 Permission is granted to copy and distribute modified versions of
Chris@19 25 this manual under the conditions for verbatim copying, provided
Chris@19 26 that the entire resulting derived work is distributed under the
Chris@19 27 terms of a permission notice identical to this one.
Chris@19 28
Chris@19 29 Permission is granted to copy and distribute translations of this
Chris@19 30 manual into another language, under the above conditions for
Chris@19 31 modified versions, except that this permission notice may be
Chris@19 32 stated in a translation approved by the Free Software Foundation.
Chris@19 33 -->
Chris@19 34 <meta http-equiv="Content-Style-Type" content="text/css">
Chris@19 35 <style type="text/css"><!--
Chris@19 36 pre.display { font-family:inherit }
Chris@19 37 pre.format { font-family:inherit }
Chris@19 38 pre.smalldisplay { font-family:inherit; font-size:smaller }
Chris@19 39 pre.smallformat { font-family:inherit; font-size:smaller }
Chris@19 40 pre.smallexample { font-size:smaller }
Chris@19 41 pre.smalllisp { font-size:smaller }
Chris@19 42 span.sc { font-variant:small-caps }
Chris@19 43 span.roman { font-family:serif; font-weight:normal; }
Chris@19 44 span.sansserif { font-family:sans-serif; font-weight:normal; }
Chris@19 45 --></style>
Chris@19 46 </head>
Chris@19 47 <body>
Chris@19 48 <div class="node">
Chris@19 49 <a name="New-array-Execute-Functions"></a>
Chris@19 50 <a name="New_002darray-Execute-Functions"></a>
Chris@19 51 <p>
Chris@19 52 Next:&nbsp;<a rel="next" accesskey="n" href="Wisdom.html#Wisdom">Wisdom</a>,
Chris@19 53 Previous:&nbsp;<a rel="previous" accesskey="p" href="Guru-Interface.html#Guru-Interface">Guru Interface</a>,
Chris@19 54 Up:&nbsp;<a rel="up" accesskey="u" href="FFTW-Reference.html#FFTW-Reference">FFTW Reference</a>
Chris@19 55 <hr>
Chris@19 56 </div>
Chris@19 57
Chris@19 58 <h3 class="section">4.6 New-array Execute Functions</h3>
Chris@19 59
Chris@19 60 <p><a name="index-execute-266"></a><a name="index-new_002darray-execution-267"></a>
Chris@19 61 Normally, one executes a plan for the arrays with which the plan was
Chris@19 62 created, by calling <code>fftw_execute(plan)</code> as described in <a href="Using-Plans.html#Using-Plans">Using Plans</a>.
Chris@19 63 <a name="index-fftw_005fexecute-268"></a>However, it is possible for sophisticated users to apply a given plan
Chris@19 64 to a <em>different</em> array using the &ldquo;new-array execute&rdquo; functions
Chris@19 65 detailed below, provided that the following conditions are met:
Chris@19 66
Chris@19 67 <ul>
Chris@19 68 <li>The array size, strides, etcetera are the same (since those are set by
Chris@19 69 the plan).
Chris@19 70
Chris@19 71 <li>The input and output arrays are the same (in-place) or different
Chris@19 72 (out-of-place) if the plan was originally created to be in-place or
Chris@19 73 out-of-place, respectively.
Chris@19 74
Chris@19 75 <li>For split arrays, the separations between the real and imaginary
Chris@19 76 parts, <code>ii-ri</code> and <code>io-ro</code>, are the same as they were for
Chris@19 77 the input and output arrays when the plan was created. (This
Chris@19 78 condition is automatically satisfied for interleaved arrays.)
Chris@19 79
Chris@19 80 <li>The <dfn>alignment</dfn> of the new input/output arrays is the same as that
Chris@19 81 of the input/output arrays when the plan was created, unless the plan
Chris@19 82 was created with the <code>FFTW_UNALIGNED</code> flag.
Chris@19 83 <a name="index-FFTW_005fUNALIGNED-269"></a>Here, the alignment is a platform-dependent quantity (for example, it is
Chris@19 84 the address modulo 16 if SSE SIMD instructions are used, but the address
Chris@19 85 modulo 4 for non-SIMD single-precision FFTW on the same machine). In
Chris@19 86 general, only arrays allocated with <code>fftw_malloc</code> are guaranteed to
Chris@19 87 be equally aligned (see <a href="SIMD-alignment-and-fftw_005fmalloc.html#SIMD-alignment-and-fftw_005fmalloc">SIMD alignment and fftw_malloc</a>).
Chris@19 88
Chris@19 89 </ul>
Chris@19 90
Chris@19 91 <p><a name="index-alignment-270"></a>The alignment issue is especially critical, because if you don't use
Chris@19 92 <code>fftw_malloc</code> then you may have little control over the alignment
Chris@19 93 of arrays in memory. For example, neither the C++ <code>new</code> function
Chris@19 94 nor the Fortran <code>allocate</code> statement provide strong enough
Chris@19 95 guarantees about data alignment. If you don't use <code>fftw_malloc</code>,
Chris@19 96 therefore, you probably have to use <code>FFTW_UNALIGNED</code> (which
Chris@19 97 disables most SIMD support). If possible, it is probably better for
Chris@19 98 you to simply create multiple plans (creating a new plan is quick once
Chris@19 99 one exists for a given size), or better yet re-use the same array for
Chris@19 100 your transforms.
Chris@19 101
Chris@19 102 <p><a name="index-fftw_005falignment_005fof-271"></a>For rare circumstances in which you cannot control the alignment of
Chris@19 103 allocated memory, but wish to determine where a given array is
Chris@19 104 aligned like the original array for which a plan was created, you can
Chris@19 105 use the <code>fftw_alignment_of</code> function:
Chris@19 106 <pre class="example"> int fftw_alignment_of(double *p);
Chris@19 107 </pre>
Chris@19 108 <p>Two arrays have equivalent alignment (for the purposes of applying a
Chris@19 109 plan) if and only if <code>fftw_alignment_of</code> returns the same value
Chris@19 110 for the corresponding pointers to their data (typecast to <code>double*</code>
Chris@19 111 if necessary).
Chris@19 112
Chris@19 113 <p>If you are tempted to use the new-array execute interface because you
Chris@19 114 want to transform a known bunch of arrays of the same size, you should
Chris@19 115 probably go use the advanced interface instead (see <a href="Advanced-Interface.html#Advanced-Interface">Advanced Interface</a>)).
Chris@19 116
Chris@19 117 <p>The new-array execute functions are:
Chris@19 118
Chris@19 119 <pre class="example"> void fftw_execute_dft(
Chris@19 120 const fftw_plan p,
Chris@19 121 fftw_complex *in, fftw_complex *out);
Chris@19 122
Chris@19 123 void fftw_execute_split_dft(
Chris@19 124 const fftw_plan p,
Chris@19 125 double *ri, double *ii, double *ro, double *io);
Chris@19 126
Chris@19 127 void fftw_execute_dft_r2c(
Chris@19 128 const fftw_plan p,
Chris@19 129 double *in, fftw_complex *out);
Chris@19 130
Chris@19 131 void fftw_execute_split_dft_r2c(
Chris@19 132 const fftw_plan p,
Chris@19 133 double *in, double *ro, double *io);
Chris@19 134
Chris@19 135 void fftw_execute_dft_c2r(
Chris@19 136 const fftw_plan p,
Chris@19 137 fftw_complex *in, double *out);
Chris@19 138
Chris@19 139 void fftw_execute_split_dft_c2r(
Chris@19 140 const fftw_plan p,
Chris@19 141 double *ri, double *ii, double *out);
Chris@19 142
Chris@19 143 void fftw_execute_r2r(
Chris@19 144 const fftw_plan p,
Chris@19 145 double *in, double *out);
Chris@19 146 </pre>
Chris@19 147 <p><a name="index-fftw_005fexecute_005fdft-272"></a><a name="index-fftw_005fexecute_005fsplit_005fdft-273"></a><a name="index-fftw_005fexecute_005fdft_005fr2c-274"></a><a name="index-fftw_005fexecute_005fsplit_005fdft_005fr2c-275"></a><a name="index-fftw_005fexecute_005fdft_005fc2r-276"></a><a name="index-fftw_005fexecute_005fsplit_005fdft_005fc2r-277"></a><a name="index-fftw_005fexecute_005fr2r-278"></a>
Chris@19 148 These execute the <code>plan</code> to compute the corresponding transform on
Chris@19 149 the input/output arrays specified by the subsequent arguments. The
Chris@19 150 input/output array arguments have the same meanings as the ones passed
Chris@19 151 to the guru planner routines in the preceding sections. The <code>plan</code>
Chris@19 152 is not modified, and these routines can be called as many times as
Chris@19 153 desired, or intermixed with calls to the ordinary <code>fftw_execute</code>.
Chris@19 154
Chris@19 155 <p>The <code>plan</code> <em>must</em> have been created for the transform type
Chris@19 156 corresponding to the execute function, e.g. it must be a complex-DFT
Chris@19 157 plan for <code>fftw_execute_dft</code>. Any of the planner routines for that
Chris@19 158 transform type, from the basic to the guru interface, could have been
Chris@19 159 used to create the plan, however.
Chris@19 160
Chris@19 161 <!-- -->
Chris@19 162 </body></html>
Chris@19 163