fftw3: Guru vector and transform sizes

 
 4.5.2 Guru vector and transform sizes
 -------------------------------------
 
 The guru interface introduces one basic new data structure,
 'fftw_iodim', that is used to specify sizes and strides for
 multi-dimensional transforms and vectors:
 
      typedef struct {
           int n;
           int is;
           int os;
      } fftw_iodim;
 
    Here, 'n' is the size of the dimension, and 'is' and 'os' are the
 strides of that dimension for the input and output arrays.  (The stride
 is the separation of consecutive elements along this dimension.)
 
    The meaning of the stride parameter depends on the type of the array
 that the stride refers to.  _If the array is interleaved complex,
 strides are expressed in units of complex numbers ('fftw_complex').  If
 the array is split complex or real, strides are expressed in units of
 real numbers ('double')._  This convention is consistent with the usual
 pointer arithmetic in the C language.  An interleaved array is denoted
 by a pointer 'p' to 'fftw_complex', so that 'p+1' points to the next
 complex number.  Split arrays are denoted by pointers to 'double', in
 which case pointer arithmetic operates in units of 'sizeof(double)'.
 
    The guru planner interfaces all take a ('rank', 'dims[rank]') pair
 describing the transform size, and a ('howmany_rank',
 'howmany_dims[howmany_rank]') pair describing the "vector" size (a
 multi-dimensional loop of transforms to perform), where 'dims' and
 'howmany_dims' are arrays of 'fftw_iodim'.
 
    For example, the 'howmany' parameter in the advanced complex-DFT
 interface corresponds to 'howmany_rank' = 1, 'howmany_dims[0].n' =
 'howmany', 'howmany_dims[0].is' = 'idist', and 'howmany_dims[0].os' =
 'odist'.  (To compute a single transform, you can just use
 'howmany_rank' = 0.)
 
    A row-major multidimensional array with dimensions 'n[rank]' (See
 Row-major Format) corresponds to 'dims[i].n' = 'n[i]' and the
 recurrence 'dims[i].is' = 'n[i+1] * dims[i+1].is' (similarly for 'os').
 The stride of the last ('i=rank-1') dimension is the overall stride of
 the array.  e.g.  to be equivalent to the advanced complex-DFT
 interface, you would have 'dims[rank-1].is' = 'istride' and
 'dims[rank-1].os' = 'ostride'.
 
    In general, we only guarantee FFTW to return a non-'NULL' plan if the
 vector and transform dimensions correspond to a set of distinct indices,
 and for in-place transforms the input/output strides should be the same.