Copyright (C) 2000-2012 |
GNU Info (fftw.info)rfftwComputing the Real One-dimensional Transform -------------------------------------------- #include <rfftw.h> void rfftw(rfftw_plan plan, int howmany, fftw_real *in, int istride, int idist, fftw_real *out, int ostride, int odist); void rfftw_one(rfftw_plan plan, fftw_real *in, fftw_real *out); The function `rfftw' computes the Real One-dimensional Fourier Transform, using a plan created by `rfftw_create_plan' (Note: Plan Creation for Real One-dimensional Transforms.). The function `rfftw_one' provides a simplified interface for the common case of single input array of stride 1. *Important:* When invoked for an out-of-place, `FFTW_COMPLEX_TO_REAL' transform, the input array is overwritten with scratch values by these routines. The input array is not modified for `FFTW_REAL_TO_COMPLEX' transforms. Arguments ......... * `plan' is the plan created by `rfftw_create_plan' (Note: Plan Creation for Real One-dimensional Transforms.). * `howmany' is the number of transforms `rfftw' will compute. It is faster to tell RFFTW to compute many transforms, instead of simply calling `rfftw' many times. * `in', `istride' and `idist' describe the input array(s). There are two cases. If the `plan' defines a `FFTW_REAL_TO_COMPLEX' transform, `in' is a real array. Otherwise, for `FFTW_COMPLEX_TO_REAL' transforms, `in' is a halfcomplex array *whose contents will be destroyed*. * `out', `ostride' and `odist' describe the output array(s), and have the same meaning as the corresponding parameters 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 `rfftw_one' transforms a single, contiguous input array to a contiguous output array. By definition, the call rfftw_one(plan, in, out) is equivalent to rfftw(plan, 1, in, 1, 1, out, 1, 1) automatically generated by info2www version 1.2.2.9 |