pyFFTW is a pythonic wrapper around FFTW, the speedy FFT library. The ultimate aim is to present a unified interface for all the possible transforms that FFTW can perform.
Both the complex DFT and the real DFT are supported, as well as on arbitrary axes of abitrary shaped and strided arrays, which makes it almost feature equivalent to standard and real FFT functions of numpy.fft (indeed, it supports the clongdouble dtype which numpy.fft does not). See the numpy.fft documentation for more information on that module.
The source can be found in github and it’s page in the python package index is here.
Operating FFTW in multithreaded mode is supported.
A comprehensive unittest suite is included with the source on the repository.
A quick (and very much non-comprehensive) usage example:
>>> import pyfftw, numpy
>>> # Create 3 16-byte aligned arays using the aligned array creation functions
>>> # They are cleared during the object creation, so there is no point filling them.
>>> a = pyfftw.n_byte_align_empty((1,4), 16, dtype=numpy.complex128)
>>> b = pyfftw.n_byte_align_empty(a.shape, 16, dtype=a.dtype)
>>> c = pyfftw.n_byte_align_empty(a.shape, 16, dtype=a.dtype)
>>> # Create the DFT and inverse DFT objects
>>> fft = pyfftw.FFTW(a, b)
>>> ifft = pyfftw.FFTW(b, c, direction='FFTW_BACKWARD')
>>> # Fill a with data
>>> a[:] = [1, 2, 3, 4]
>>> print a
[[ 1.+0.j 2.+0.j 3.+0.j 4.+0.j]]
>>> # perform the fft
>>> fft.execute()
>>> print b
[[ 10.+0.j -2.+2.j -2.+0.j -2.-2.j]]
>>> # perform the inverse fft
>>> ifft.execute()
>>> print c/a.size
[[ 1.+0.j 2.+0.j 3.+0.j 4.+0.j]]
Note that what was previously the ComplexFFTW class is now just called FFTW. Simply renaming the class should be sufficient to migrate
FFTW is a class for computing the complex N-Dimensional DFT or inverse DFT of an array using the FFTW library. The interface is designed to be somewhat pythonic, with the correct transform being inferred from the dtypes of the passed arrays.
On instantiation, the dtypes of the input arrays are compared to the set of valid (and implemented) FFTW schemes. If a match is found, the plan that corresponds to that scheme is created, operating on the arrays that are passed in. If no scheme can be created, then ValueError is raised.
The actual FFT or iFFT is performed by calling the execute() method.
The arrays can be updated by calling the update_arrays() method.
input_array and output_array should be numpy arrays. The contents of these arrays will be destroyed by the planning process during initialisation. Information on supported dtypes for the arrays is given below.
axes describes along which axes the DFT should be taken. This should be a valid list of axes. Repeated axes are only transformed once. Invalid axes will raise an exception. This argument is equivalent to the same argument in numpy.fft.fftn, except for the fact that the behaviour of repeated axes is different (numpy.fft will happily take the fft of the same axis if it is repeated in the axes argument). Rudimentary testing has suggested this is down to FFTW and so unlikely to be fixed in these wrappers.
direction should be a string and one of FFTW_FORWARD or FFTW_BACKWARD, which dictate whether to take the DFT (forwards) or the inverse DFT (backwards) respectively (specifically, it dictates the sign of the exponent in the DFT formulation).
Note that only the Complex schemes allow a free choice for direction. The direction must agree with the the table below if a Real scheme is used, otherwise a ValueError is raised.
flags is a list of strings and is a subset of the flags that FFTW allows for the planners. Specifically, FFTW_ESTIMATE, FFTW_MEASURE, FFTW_PATIENT and FFTW_EXHAUSTIVE are supported. These describe the increasing amount of effort spent during the planning stage to create the fastest possible transform. Usually, FFTW_MEASURE is a good compromise and is the default.
In addition the FFTW_UNALIGNED flag is supported. This tells FFTW not to assume anything about the alignment of the data and disabling any SIMD capability (see below).
threads tells the wrapper how many threads to use when invoking FFTW, with a default of 1. If the number of threads is greater than 1, then the GIL is released by necessity.
The currently supported schemes are as follows:
Type | input_array.dtype | output_array.dtype | Direction |
---|---|---|---|
Complex | complex64 | complex64 | Both |
Complex | complex128 | complex128 | Both |
Complex | clongdouble | clongdouble | Both |
Real | float32 | complex64 | Forwards |
Real | float64 | complex128 | Forwards |
Real | longdouble | clongdouble | Forwards |
Real1 | complex64 | float32 | Backwards |
Real1 | complex128 | float64 | Backwards |
Real1 | clongdouble | longdouble | Backwards |
1 Note that the Backwards Real transform will destroy the input array. This is inherent to FFTW and the only general work-around for this is to copy the array prior to performing the transform.
clongdouble typically maps directly to complex256 or complex192, and longdouble to float128 or float96, dependent on platform.
The relative shapes of the arrays should be as follows:
In the above expressions for the Real transform, the axes arguments denotes the the unique set of axes on which we are taking the FFT, in the order passed. It is the last of these axes that is subject to the special case shown.
The shapes for the real transforms corresponds to those stipulated by the FFTW library. Further information can be found in the FFTW documentation on the real DFT.
The actual arrangement in memory is arbitrary and the scheme can be planned for any set of strides on either the input or the output. The user should not have to worry about this and any valid numpy array should work just fine.
What is calculated is exactly what FFTW calculates. Notably, this is an unnormalized transform so should be scaled as necessary (fft followed by ifft will scale the input by N, the product of the dimensions along which the DFT is taken). For further information, see the FFTW documentation.
The FFTW library benefits greatly from the beginning of each DFT axes being aligned on a 16-byte boundary, which enables SIMD instructions. By default, if the data begins on a 16-byte boundary, then FFTW will be allowed to try and enable SIMD instructions. This means that all future changes to the data arrays will be checked for similar alignment. SIMD instructions can be explicitly disabled by setting the FFTW_UNALIGNED flags, to allow for updates with unaligned data.
n_byte_align() and n_byte_align_empty() are two methods included with this module for producing aligned arrays.
Update the arrays upon which the DFT is taken.
The new arrays should be of the same dtypes as the originals, the same shapes as the originals and should have the same strides between axes. If the original data was aligned so as to allow SIMD instructions (by being aligned on a 16-byte boundary), then the new array must also be aligned in the same way.
If all these conditions are not met, a ValueError will be raised and the data will not be updated (though the object will still be in a sane state).
Note that if the original array was not aligned on a 16-byte boundary, then SIMD is disabled and the alignment of the new array can be arbitrary.
Execute the planned operation.
Functions for dealing with FFTW’s ability to export and restore plans, referred to as wisdom. For further information, refer to the FFTW wisdom documentation.
Return the FFTW wisdom as a tuple of strings.
The first string in the tuple is the string for the double precision wisdom. The second string in the tuple is the string for the single precision wisdom. The third string in the tuple is the string for the long double precision wisdom.
The tuple that is returned from this function can be used as the argument to import_wisdom().
Function that imports wisdom from the passed tuple of strings.
The first string in the tuple is the string for the double precision wisdom. The second string in the tuple is the string for the single precision wisdom. The third string in the tuple is the string for the long double precision wisdom.
The tuple that is returned from export_wisdom() can be used as the argument to this function.
This function returns a tuple of boolean values indicating the success of loading each of the wisdom types (double, float and long double, in that order).
Forget all the accumulated wisdom.
Function that takes a numpy array and checks it is aligned on an n-byte boundary, where n is a passed parameter. If it is, the array is returned without further ado. If it is not, a new array is created and the data copied in, but aligned on the n-byte boundary.
Function that returns an empty numpy array that is n-byte aligned.
The alignment is given by the second argument, n. The rest of the arguments are as per numpy.empty.