convolve_fft#

astropy.convolution.convolve_fft(array, kernel, boundary='fill', fill_value=0.0, nan_treatment='interpolate', normalize_kernel=True, normalization_zero_tol=1e-08, preserve_nan=False, mask=None, crop=True, return_fft=False, fft_pad=None, psf_pad=None, min_wt=0.0, allow_huge=False, fftn=<function fftn>, ifftn=<function ifftn>, complex_dtype=<class 'complex'>, dealias=False)[source]#

Convolve an ndarray with an nd-kernel. Returns a convolved image with shape = array.shape. Assumes kernel is centered.

convolve_fft is very similar to convolve in that it replaces NaN values in the original image with interpolated values using the kernel as an interpolation function. However, it also includes many additional options specific to the implementation.

convolve_fft differs from scipy.signal.fftconvolve in a few ways:

  • It can treat NaN values as zeros or interpolate over them.

  • inf values are treated as NaN

  • It optionally pads to the nearest faster sizes to improve FFT speed. These sizes are optimized for the numpy and scipy implementations, and fftconvolve uses them by default as well; when using other external functions (see below), results may vary.

  • Its only valid mode is ‘same’ (i.e., the same shape array is returned)

  • It lets you use your own fft, e.g., pyFFTW or pyFFTW3 , which can lead to performance improvements, depending on your system configuration. pyFFTW3 is threaded, and therefore may yield significant performance benefits on multi-core machines at the cost of greater memory requirements. Specify the fftn and ifftn keywords to override the default, which is numpy.fft.fftn and numpy.fft.ifftn. The scipy.fft functions also offer somewhat better performance and a multi-threaded option.

Parameters:
arraynumpy.ndarray

Array to be convolved with kernel. It can be of any dimensionality, though only 1, 2, and 3d arrays have been tested.

kernelnumpy.ndarray or astropy.convolution.Kernel

The convolution kernel. The number of dimensions should match those for the array. The dimensions do not have to be odd in all directions, unlike in the non-fft convolve function. The kernel will be normalized if normalize_kernel is set. It is assumed to be centered (i.e., shifts may result if your kernel is asymmetric)

boundary{‘fill’, ‘wrap’}, optional

A flag indicating how to handle boundaries:

  • ‘fill’: set values outside the array boundary to fill_value (default)

  • ‘wrap’: periodic boundary

The None and ‘extend’ parameters are not supported for FFT-based convolution.

fill_valuefloat, optional

The value to use outside the array when using boundary=’fill’.

nan_treatment{‘interpolate’, ‘fill’}, optional
The method used to handle NaNs in the input array:
  • 'interpolate': NaN values are replaced with interpolated values using the kernel as an interpolation function. Note that if the kernel has a sum equal to zero, NaN interpolation is not possible and will raise an exception.

  • 'fill': NaN values are replaced by fill_value prior to convolution.

normalize_kernelcallable() or bool, optional

If specified, this is the function to divide kernel by to normalize it. e.g., normalize_kernel=np.sum means that kernel will be modified to be: kernel = kernel / np.sum(kernel). If True, defaults to normalize_kernel = np.sum.

normalization_zero_tolfloat, optional

The absolute tolerance on whether the kernel is different than zero. If the kernel sums to zero to within this precision, it cannot be normalized. Default is “1e-8”.

preserve_nanbool, optional

After performing convolution, should pixels that were originally NaN again become NaN?

maskNone or ndarray, optional

A “mask” array. Shape must match array, and anything that is masked (i.e., not 0/False) will be set to NaN for the convolution. If None, no masking will be performed unless array is a masked array. If mask is not None and array is a masked array, a pixel is masked if it is masked in either mask or array.mask.

cropbool, optional

Default on. Return an image of the size of the larger of the input image and the kernel. If the image and kernel are asymmetric in opposite directions, will return the largest image in both directions. For example, if an input image has shape [100,3] but a kernel with shape [6,6] is used, the output will be [100,6].

return_fftbool, optional

Return the fft(image)*fft(kernel) instead of the convolution (which is ifft(fft(image)*fft(kernel))). Useful for making PSDs.

fft_padbool, optional

Default on. Zero-pad image to the nearest size supporting more efficient execution of the FFT, generally values factorizable into the first 3-5 prime numbers. With boundary='wrap', this will be disabled.

psf_padbool, optional

Zero-pad image to be at least the sum of the image sizes to avoid edge-wrapping when smoothing. This is enabled by default with boundary='fill', but it can be overridden with a boolean option. boundary='wrap' and psf_pad=True are not compatible.

min_wtfloat, optional

If ignoring NaN / zeros, force all grid points with a weight less than this value to NaN (the weight of a grid point with no ignored neighbors is 1.0). If min_wt is zero, then all zero-weight points will be set to zero instead of NaN (which they would be otherwise, because 1/0 = nan). See the examples below.

allow_hugebool, optional

Allow huge arrays in the FFT? If False, will raise an exception if the array or kernel size is >1 GB.

fftncallable(), optional

The fft function. Can be overridden to use your own ffts, e.g. an fftw3 wrapper or scipy’s fftn, fft=scipy.fftpack.fftn.

ifftncallable(), optional

The inverse fft function. Can be overridden the same way fttn.

complex_dtypecomplex type, optional

Which complex dtype to use. numpy has a range of options, from 64 to 256.

dealias: bool, optional

Default off. Zero-pad image to enable explicit dealiasing of convolution. With boundary='wrap', this will be disabled. Note that for an input of nd dimensions this will increase the size of the temporary arrays by at least 1.5**nd. This may result in significantly more memory usage.

Returns:
defaultndarray

array convolved with kernel. If return_fft is set, returns fft(array) * fft(kernel). If crop is not set, returns the image, but with the fft-padded size instead of the input size.

Raises:
ValueError

If the array is bigger than 1 GB after padding, will raise this exception unless allow_huge is True.

See also

convolve

Convolve is a non-fft version of this code. It is more memory efficient and for small kernels can be faster.

Notes

With psf_pad=True and a large PSF, the resulting data can become large and consume a lot of memory. See Issue astropy/astropy#4366 and the update in astropy/astropy#11533 for further details.

Dealiasing of pseudospectral convolutions is necessary for numerical stability of the underlying algorithms. A common method for handling this is to zero pad the image by at least 1/2 to eliminate the wavenumbers which have been aliased by convolution. This is so that the aliased 1/3 of the results of the convolution computation can be thrown out. See https://doi.org/10.1175/1520-0469(1971)028%3C1074:OTEOAI%3E2.0.CO;2 https://iopscience.iop.org/article/10.1088/1742-6596/318/7/072037

Note that if dealiasing is necessary to your application, but your process is memory constrained, you may want to consider using FFTW++: dealias/fftwpp. It includes python wrappers for a pseudospectral convolution which will implicitly dealias your convolution without the need for additional padding. Note that one cannot use FFTW++’s convlution directly in this method as in handles the entire convolution process internally. Additionally, FFTW++ includes other useful pseudospectral methods to consider.

Examples

>>> convolve_fft([1, 0, 3], [1, 1, 1])
array([0.33333333, 1.33333333, 1.        ])
>>> convolve_fft([1, np.nan, 3], [1, 1, 1])
array([0.5, 2. , 1.5])
>>> convolve_fft([1, 0, 3], [0, 1, 0])  
array([ 1.00000000e+00, -3.70074342e-17,  3.00000000e+00])
>>> convolve_fft([1, 2, 3], [1])
array([1., 2., 3.])
>>> convolve_fft([1, np.nan, 3], [0, 1, 0], nan_treatment='interpolate')
array([1., 0., 3.])
>>> convolve_fft([1, np.nan, 3], [0, 1, 0], nan_treatment='interpolate',
...              min_wt=1e-8)
array([ 1., nan,  3.])
>>> convolve_fft([1, np.nan, 3], [1, 1, 1], nan_treatment='interpolate')
array([0.5, 2. , 1.5])
>>> convolve_fft([1, np.nan, 3], [1, 1, 1], nan_treatment='interpolate',
...               normalize_kernel=True)
array([0.5, 2. , 1.5])
>>> import scipy.fft  # optional - requires scipy
>>> convolve_fft([1, np.nan, 3], [1, 1, 1], nan_treatment='interpolate',
...               normalize_kernel=True,
...               fftn=scipy.fft.fftn, ifftn=scipy.fft.ifftn)
array([0.5, 2. , 1.5])
>>> fft_mp = lambda a: scipy.fft.fftn(a, workers=-1)  # use all available cores
>>> ifft_mp = lambda a: scipy.fft.ifftn(a, workers=-1)
>>> convolve_fft([1, np.nan, 3], [1, 1, 1], nan_treatment='interpolate',
...               normalize_kernel=True, fftn=fft_mp, ifftn=ifft_mp)
array([0.5, 2. , 1.5])