Copyright (C) 2000-2012 |
GNU Info (fftw.info)rfftw_create_planPlan Creation for Real One-dimensional Transforms ------------------------------------------------- #include <rfftw.h> rfftw_plan rfftw_create_plan(int n, fftw_direction dir, int flags); rfftw_plan rfftw_create_plan_specific(int n, fftw_direction dir, int flags, fftw_real *in, int istride, fftw_real *out, int ostride); The function `rfftw_create_plan' creates a plan, which is a data structure containing all the information that `rfftw' needs in order to compute the 1D real Fourier transform. You can create as many plans as you need, but only one plan for a given array size is required (a plan can be reused many times). `rfftw_create_plan' returns a valid plan, or `NULL' if, for some reason, the plan can't be created. In the default installation, this cannot happen, but it is possible to configure RFFTW in such a way that some input sizes are forbidden, and RFFTW cannot create a plan. The `rfftw_create_plan_specific' variant takes as additional arguments specific input/output arrays and their strides. For the last four arguments, you should pass the arrays and strides that you will eventually be passing to `rfftw'. The resulting plans will be optimized for those arrays and strides, although they may be used on other arrays as well. Note: the contents of the in and out arrays are *destroyed* by the specific planner (the initial contents are ignored, so the arrays need not have been initialized). Note: Discussion on Specific Plans, for a discussion on specific plans. Arguments ......... * `n' is the size of the transform. It can be any positive integer. - RFFTW is best at handling sizes of the form 2^a 3^b 5^c 7^d 11^e 13^f, where e+f is either 0 or 1, and the other exponents are arbitrary. Other sizes are computed by means of a slow, general-purpose routine (reducing to O(n^2) performance for prime sizes). (It is possible to customize RFFTW for different array sizes. Note: Installation and Customization, for more information.) Transforms whose sizes are powers of 2 are especially fast. * `dir' is the direction of the desired transform, either `FFTW_REAL_TO_COMPLEX' or `FFTW_COMPLEX_TO_REAL', corresponding to `FFTW_FORWARD' or `FFTW_BACKWARD', respectively. * `flags' is a boolean OR (`|') of zero or more of the following: - `FFTW_MEASURE': this flag tells RFFTW to find the optimal plan by actually *computing* several FFTs and measuring their execution time. Depending on the installation, this can take some time. - `FFTW_ESTIMATE': do not run any FFT and provide a "reasonable" plan (for a RISC processor with many registers). If neither `FFTW_ESTIMATE' nor `FFTW_MEASURE' is provided, the default is `FFTW_ESTIMATE'. - `FFTW_OUT_OF_PLACE': produce a plan assuming that the input and output arrays will be distinct (this is the default). - `FFTW_IN_PLACE': produce a plan assuming that you want the output in the input array. The algorithm used is not necessarily in place: RFFTW is able to compute true in-place transforms only for small values of `n'. If RFFTW is not able to compute the transform in-place, it will allocate a temporary array (unless you provide one yourself), compute the transform out of place, and copy the result back. *Warning: This option changes the meaning of some parameters of `rfftw'* (Note: Computing the Real One-dimensional Transform.). The default mode of operation is `FFTW_OUT_OF_PLACE'. - `FFTW_USE_WISDOM': use any `wisdom' that is available to help in the creation of the plan. (Note: Words of Wisdom.) This can greatly speed the creation of plans, especially with the `FFTW_MEASURE' option. `FFTW_ESTIMATE' plans can also take advantage of `wisdom' to produce a more optimal plan (based on past measurements) than the estimation heuristic would normally generate. When the `FFTW_MEASURE' option is used, new `wisdom' will also be generated if the current transform size is not completely understood by existing `wisdom'. * `in', `out', `istride', `ostride' (only for `rfftw_create_plan_specific'): see corresponding arguments in the description of `rfftw'. (Note: Computing the Real One-dimensional Transform.) In particular, the `out' and `ostride' parameters have the same special meaning for `FFTW_IN_PLACE' transforms as they have for `rfftw'. automatically generated by info2www version 1.2.2.9 |