This file is indexed.

/usr/share/doc/libfftw3-doc/html/Advanced-Complex-DFTs.html is in libfftw3-doc 3.3-1ubuntu1.

This file is owned by root:root, with mode 0o644.

The actual contents of the file can be viewed below.

  1
  2
  3
  4
  5
  6
  7
  8
  9
 10
 11
 12
 13
 14
 15
 16
 17
 18
 19
 20
 21
 22
 23
 24
 25
 26
 27
 28
 29
 30
 31
 32
 33
 34
 35
 36
 37
 38
 39
 40
 41
 42
 43
 44
 45
 46
 47
 48
 49
 50
 51
 52
 53
 54
 55
 56
 57
 58
 59
 60
 61
 62
 63
 64
 65
 66
 67
 68
 69
 70
 71
 72
 73
 74
 75
 76
 77
 78
 79
 80
 81
 82
 83
 84
 85
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
<html lang="en">
<head>
<title>Advanced Complex DFTs - FFTW 3.3</title>
<meta http-equiv="Content-Type" content="text/html">
<meta name="description" content="FFTW 3.3">
<meta name="generator" content="makeinfo 4.13">
<link title="Top" rel="start" href="index.html#Top">
<link rel="up" href="Advanced-Interface.html#Advanced-Interface" title="Advanced Interface">
<link rel="prev" href="Advanced-Interface.html#Advanced-Interface" title="Advanced Interface">
<link rel="next" href="Advanced-Real_002ddata-DFTs.html#Advanced-Real_002ddata-DFTs" title="Advanced Real-data DFTs">
<link href="http://www.gnu.org/software/texinfo/" rel="generator-home" title="Texinfo Homepage">
<!--
This manual is for FFTW
(version 3.3, 26 July 2011).

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.
   -->
<meta http-equiv="Content-Style-Type" content="text/css">
<style type="text/css"><!--
  pre.display { font-family:inherit }
  pre.format  { font-family:inherit }
  pre.smalldisplay { font-family:inherit; font-size:smaller }
  pre.smallformat  { font-family:inherit; font-size:smaller }
  pre.smallexample { font-size:smaller }
  pre.smalllisp    { font-size:smaller }
  span.sc    { font-variant:small-caps }
  span.roman { font-family:serif; font-weight:normal; } 
  span.sansserif { font-family:sans-serif; font-weight:normal; } 
--></style>
</head>
<body>
<div class="node">
<a name="Advanced-Complex-DFTs"></a>
<p>
Next:&nbsp;<a rel="next" accesskey="n" href="Advanced-Real_002ddata-DFTs.html#Advanced-Real_002ddata-DFTs">Advanced Real-data DFTs</a>,
Previous:&nbsp;<a rel="previous" accesskey="p" href="Advanced-Interface.html#Advanced-Interface">Advanced Interface</a>,
Up:&nbsp;<a rel="up" accesskey="u" href="Advanced-Interface.html#Advanced-Interface">Advanced Interface</a>
<hr>
</div>

<h4 class="subsection">4.4.1 Advanced Complex DFTs</h4>

<pre class="example">     fftw_plan fftw_plan_many_dft(int rank, const int *n, int howmany,
                                  fftw_complex *in, const int *inembed,
                                  int istride, int idist,
                                  fftw_complex *out, const int *onembed,
                                  int ostride, int odist,
                                  int sign, unsigned flags);
</pre>
   <p><a name="index-fftw_005fplan_005fmany_005fdft-232"></a>
This routine plans multiple multidimensional complex DFTs, and it
extends the <code>fftw_plan_dft</code> routine (see <a href="Complex-DFTs.html#Complex-DFTs">Complex DFTs</a>) to
compute <code>howmany</code> transforms, each having rank <code>rank</code> and size
<code>n</code>.  In addition, the transform data need not be contiguous, but
it may be laid out in memory with an arbitrary stride.  To account for
these possibilities, <code>fftw_plan_many_dft</code> adds the new parameters
<code>howmany</code>, {<code>i</code>,<code>o</code>}<code>nembed</code>,
{<code>i</code>,<code>o</code>}<code>stride</code>, and
{<code>i</code>,<code>o</code>}<code>dist</code>.  The FFTW basic interface
(see <a href="Complex-DFTs.html#Complex-DFTs">Complex DFTs</a>) provides routines specialized for ranks 1, 2,
and&nbsp;3, but the advanced interface handles only the general-rank
case.

   <p><code>howmany</code> is the number of transforms to compute.  The resulting
plan computes <code>howmany</code> transforms, where the input of the
<code>k</code>-th transform is at location <code>in+k*idist</code> (in C pointer
arithmetic), and its output is at location <code>out+k*odist</code>.  Plans
obtained in this way can often be faster than calling FFTW multiple
times for the individual transforms.  The basic <code>fftw_plan_dft</code>
interface corresponds to <code>howmany=1</code> (in which case the <code>dist</code>
parameters are ignored). 
<a name="index-howmany-parameter-233"></a><a name="index-dist-234"></a>

   <p>Each of the <code>howmany</code> transforms has rank <code>rank</code> and size
<code>n</code>, as in the basic interface.  In addition, the advanced
interface allows the input and output arrays of each transform to be
row-major subarrays of larger rank-<code>rank</code> arrays, described by
<code>inembed</code> and <code>onembed</code> parameters, respectively. 
{<code>i</code>,<code>o</code>}<code>nembed</code> must be arrays of length <code>rank</code>,
and <code>n</code> should be elementwise less than or equal to
{<code>i</code>,<code>o</code>}<code>nembed</code>.  Passing <code>NULL</code> for an
<code>nembed</code> parameter is equivalent to passing <code>n</code> (i.e. same
physical and logical dimensions, as in the basic interface.)

   <p>The <code>stride</code> parameters indicate that the <code>j</code>-th element of
the input or output arrays is located at <code>j*istride</code> or
<code>j*ostride</code>, respectively.  (For a multi-dimensional array,
<code>j</code> is the ordinary row-major index.)  When combined with the
<code>k</code>-th transform in a <code>howmany</code> loop, from above, this means
that the (<code>j</code>,<code>k</code>)-th element is at <code>j*stride+k*dist</code>. 
(The basic <code>fftw_plan_dft</code> interface corresponds to a stride of 1.) 
<a name="index-stride-235"></a>

   <p>For in-place transforms, the input and output <code>stride</code> and
<code>dist</code> parameters should be the same; otherwise, the planner may
return <code>NULL</code>.

   <p>Arrays <code>n</code>, <code>inembed</code>, and <code>onembed</code> are not used after
this function returns.  You can safely free or reuse them.

   <p><strong>Examples</strong>:
One transform of one 5 by 6 array contiguous in memory:
<pre class="example">        int rank = 2;
        int n[] = {5, 6};
        int howmany = 1;
        int idist = odist = 0; /* unused because howmany = 1 */
        int istride = ostride = 1; /* array is contiguous in memory */
        int *inembed = n, *onembed = n;
</pre>
   <p>Transform of three 5 by 6 arrays, each contiguous in memory,
stored in memory one after another:
<pre class="example">        int rank = 2;
        int n[] = {5, 6};
        int howmany = 3;
        int idist = odist = n[0]*n[1]; /* = 30, the distance in memory
                                          between the first element
                                          of the first array and the
                                          first element of the second array */
        int istride = ostride = 1; /* array is contiguous in memory */
        int *inembed = n, *onembed = n;
</pre>
   <p>Transform each column of a 2d array with 10 rows and 3 columns:
<pre class="example">        int rank = 1; /* not 2: we are computing 1d transforms */
        int n[] = {10}; /* 1d transforms of length 10 */
        int howmany = 3;
        int idist = odist = 1;
        int istride = ostride = 3; /* distance between two elements in
                                      the same column */
        int *inembed = n, *onembed = n;
</pre>
   <!-- =========> -->
   </body></html>