pyfftw_call¶
- odl.trafos.backends.pyfftw_bindings.pyfftw_call(array_in, array_out, direction='forward', axes=None, halfcomplex=False, **kwargs)[source]¶
Calculate the DFT with pyfftw.
The discrete Fourier (forward) transform calcuates the sum:
f_hat[k] = sum_j( f[j] * exp(-2*pi*1j * j*k/N) )
where the summation is taken over all indices
j = (j[0], ..., j[d-1])in the range0 <= j < N(component-wise), withNbeing the shape of the input array.The output indices
klie in the same range, except for half-complex transforms, where the last axisiinaxesis shortened to0 <= k[i] < floor(N[i]/2) + 1.In the backward transform, sign of the the exponential argument is flipped.
- Parameters:
- array_in
numpy.ndarray Array to be transformed
- array_out
numpy.ndarray Output array storing the transformed values, may be aliased with
array_in.- direction{'forward', 'backward'}, optional
Direction of the transform
- axesint or sequence of ints, optional
Dimensions along which to take the transform.
Nonemeans using all axes and is equivalent tonp.arange(ndim).- halfcomplexbool, optional
If
True, calculate only the negative frequency part along the last axis. IfFalse, calculate the full complex FFT. This option can only be used with real input data.
- array_in
- Returns:
- fftw_plan
pyfftw.FFTW The plan object created from the input arguments. It can be reused for transforms of the same size with the same data types. Note that reuse only gives a speedup if the initial plan used a planner flag other than
'estimate'. Iffftw_planwas specified, the returned object is a reference to it.
- fftw_plan
- Other Parameters:
- fftw_plan
pyfftw.FFTW, optional Use this plan instead of calculating a new one. If specified, the options
planning_effort,planning_timelimitandthreadshave no effect.- planning_effortstr, optional
Flag for the amount of effort put into finding an optimal FFTW plan. See the FFTW doc on planner flags. Available options: {'estimate', 'measure', 'patient', 'exhaustive'} Default: 'estimate'
- planning_timelimitfloat or
None, optional Limit planning time to roughly this many seconds. Default:
None(no limit)- threadsint, optional
Number of threads to use. Default: Number of CPUs if the number of data points is larger than 4096, else 1.
- normalise_idftbool, optional
If
True, the result of the backward transform is divided by1 / N, whereNis the total number of points inarray_in[axes]. This ensures that the IDFT is the true inverse of the forward DFT. Default:False- import_wisdomfilename or file handle, optional
File to load FFTW wisdom from. If the file does not exist, it is ignored.
- export_wisdomfilename or file handle, optional
File to append the accumulated FFTW wisdom to
- fftw_plan
Notes
The planning and direction flags can also be specified as capitalized and prepended by
'FFTW_', i.e. in the original FFTW form.For a
halfcomplexforward transform, the arrays must fulfillarray_out.shape[axes[-1]] == array_in.shape[axes[-1]] // 2 + 1, and vice versa for backward transforms.All planning schemes except
'estimate'require an internal copy of the input array but are often several times faster after the first call (measuring results are cached). Typically, 'measure' is a good compromise. If you cannot afford the copy, use'estimate'.If a plan is provided via the
fftw_planparameter, no copy is needed internally.