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, quiet=False, min_wt=0.0, allow_huge=False, fftn=<function fftn>, ifftn=<function ifftn>, complex_dtype=<class 'complex'>)[source]¶ Convolve an ndarray with an nd-kernel. Returns a convolved image with
shape = array.shape. Assumes kernel is centered.convolve_fftis very similar toconvolvein that it replacesNaNvalues 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_fftdiffers fromscipy.signal.fftconvolvein a few ways:- It can treat
NaNvalues as zeros or interpolate over them. infvalues are treated asNaN- (optionally) It pads to the nearest 2^n size to improve FFT speed.
- Its only valid
modeis ‘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
fftnandifftnkeywords to override the default, which isnumpy.fft.fftandnumpy.fft.ifft.
Parameters: - array :
numpy.ndarray Array to be convolved with
kernel. It can be of any dimensionality, though only 1, 2, and 3d arrays have been tested.- kernel :
numpy.ndarrayorastropy.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
convolvefunction. The kernel will be normalized ifnormalize_kernelis 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
Noneand ‘extend’ parameters are not supported for FFT-based convolution- fill_value : float, optional
The value to use outside the array when using boundary=’fill’
- nan_treatment : ‘interpolate’, ‘fill’
interpolatewill result in renormalization of the kernel at each position ignoring (pixels that are NaN in the image) in both the image and the kernel.fillwill replace the NaN pixels with a fixed numerical value (default zero, seefill_value) prior to convolution. Note that if the kernel has a sum equal to zero, NaN interpolation is not possible and will raise an exception.- normalize_kernel : function or boolean, optional
If specified, this is the function to divide kernel by to normalize it. e.g.,
normalize_kernel=np.summeans that kernel will be modified to be:kernel = kernel / np.sum(kernel). If True, defaults tonormalize_kernel = np.sum.- normalization_zero_tol: float, 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_nan : bool
After performing convolution, should pixels that were originally NaN again become NaN?
- mask :
Noneornumpy.ndarray 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. IfNone, no masking will be performed unlessarrayis a masked array. Ifmaskis notNoneandarrayis a masked array, a pixel is masked of it is masked in eithermaskorarray.mask.
Returns: - default : ndarray
arrayconvolved withkernel. Ifreturn_fftis set, returnsfft(array) * fft(kernel). If crop is not set, returns the image, but with the fft-padded size instead of the input size
Other Parameters: - min_wt : float, optional
If ignoring
NaN/ zeros, force all grid points with a weight less than this value toNaN(the weight of a grid point with no ignored neighbors is 1.0). Ifmin_wtis zero, then all zero-weight points will be set to zero instead ofNaN(which they would be otherwise, because 1/0 = nan). See the examples below- fft_pad : bool, optional
Default on. Zero-pad image to the nearest 2^n. With
boundary='wrap', this will be disabled.- psf_pad : bool, 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'andpsf_pad=Trueare not compatible.- crop : bool, 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_fft : bool, optional
Return the
fft(image)*fft(kernel)instead of the convolution (which isifft(fft(image)*fft(kernel))). Useful for making PSDs.- fftn, ifftn : functions, optional
The fft and inverse fft functions. Can be overridden to use your own ffts, e.g. an fftw3 wrapper or scipy’s fftn,
fft=scipy.fftpack.fftn- complex_dtype : numpy.complex, optional
Which complex dtype to use.
numpyhas a range of options, from 64 to 256.- quiet : bool, optional
Silence warning message about NaN interpolation
- allow_huge : bool, optional
Allow huge arrays in the FFT? If False, will raise an exception if the array or kernel size is >1 GB
Raises: - ValueError:
If the array is bigger than 1 GB after padding, will raise this exception unless
allow_hugeis 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=Trueand a large PSF, the resulting data can become very large and consume a lot of memory. See Issue https://github.com/astropy/astropy/pull/4366 for further detail.Examples
>>> convolve_fft([1, 0, 3], [1, 1, 1]) array([ 1., 4., 3.])
>>> convolve_fft([1, np.nan, 3], [1, 1, 1]) array([ 1., 4., 3.])
>>> convolve_fft([1, 0, 3], [0, 1, 0]) array([ 1., 0., 3.])
>>> 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([ 1., 4., 3.])
>>> convolve_fft([1, np.nan, 3], [1, 1, 1], nan_treatment='interpolate', ... normalize_kernel=True) array([ 1., 2., 3.])
>>> import scipy.fftpack # optional - requires scipy >>> convolve_fft([1, np.nan, 3], [1, 1, 1], nan_treatment='interpolate', ... normalize_kernel=True, ... fftn=scipy.fftpack.fft, ifftn=scipy.fftpack.ifft) array([ 1., 2., 3.])
- It can treat