Copyright (C) 2000-2012 |
GNU Info (fftw.info)fftwComputing the One-dimensional Transform --------------------------------------- #include <fftw.h> void fftw(fftw_plan plan, int howmany, fftw_complex *in, int istride, int idist, fftw_complex *out, int ostride, int odist); void fftw_one(fftw_plan plan, fftw_complex *in, fftw_complex *out); The function `fftw' computes the one-dimensional Fourier transform, using a plan created by `fftw_create_plan' (Note: Plan Creation for One-dimensional Transforms.) The function `fftw_one' provides a simplified interface for the common case of single input array of stride 1. Arguments ......... * `plan' is the plan created by `fftw_create_plan' (Note: Plan Creation for One-dimensional Transforms.). * `howmany' is the number of transforms `fftw' will compute. It is faster to tell FFTW to compute many transforms, instead of simply calling `fftw' many times. * `in', `istride' and `idist' describe the input array(s). There are `howmany' input arrays; the first one is pointed to by `in', the second one is pointed to by `in + idist', and so on, up to `in + (howmany - 1) * idist'. Each input array consists of complex numbers (Note: Data Types.), which are not necessarily contiguous in memory. Specifically, `in[0]' is the first element of the first array, `in[istride]' is the second element of the first array, and so on. In general, the `i'-th element of the `j'-th input array will be in position `in[i * istride + j * idist]'. * `out', `ostride' and `odist' describe the output array(s). The format is the same as for the input array. - *In-place transforms*: If the `plan' specifies an in-place transform, `ostride' and `odist' are always ignored. If `out' is `NULL', `out' is ignored, too. Otherwise, `out' is interpreted as a pointer to an array of `n' complex numbers, that FFTW will use as temporary space to perform the in-place computation. `out' is used as scratch space and its contents destroyed. In this case, `out' must be an ordinary array whose elements are contiguous in memory (no striding). The function `fftw_one' transforms a single, contiguous input array to a contiguous output array. By definition, the call fftw_one(plan, in, out) is equivalent to fftw(plan, 1, in, 1, 1, out, 1, 1) automatically generated by info2www version 1.2.2.9 |