convolve_fft

astropy.nddata.convolution.convolve.convolve_fft(array, kernel, boundary='fill', fill_value=0, crop=True, return_fft=False, fft_pad=True, psf_pad=False, interpolate_nan=False, quiet=False, ignore_edge_zeros=False, min_wt=0.0, normalize_kernel=False, fft_type=None, nthreads=None) [edit on github][source]

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

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

  • can treat NaN’s as zeros or interpolate over them
  • defaults to using the faster FFTW algorithm if installed
  • (optionally) pads to the nearest 2^n size to improve FFT speed
  • only operates in mode=’same’ (i.e., the same shape array is returned) mode
Parameters :

array : numpy.ndarray

Array to be convolved with kernel

kernel : numpy.ndarray

Will be normalized if normalize_kernel is set. Assumed to be centered (i.e., shifts may result if your kernel is asymmetric)

boundary : {‘fill’, ‘wrap’}

A flag indicating how to handle boundaries:

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

interpolate_nan : bool

The convolution will be re-weighted assuming NAN values are meant to be ignored, not treated as zero. If this is off, all NaN values will be treated as zero.

ignore_edge_zeros : bool

Ignore the zero-pad-created zeros. This will effectively decrease the kernel area on the edges but will not re-normalize the kernel. This parameter may result in ‘edge-brightening’ effects if you’re using a normalized kernel

min_wt : float

If ignoring NANs/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 == 0.0, 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

normalize_kernel : function or boolean

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

Returns :

default : ndarray

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

Other Parameters:
 

fft_pad : bool

Default on. Zero-pad image to the nearest 2^n

psf_pad : bool

Default off. Zero-pad image to be at least the sum of the image sizes (in order to avoid edge-wrapping when smoothing)

crop : bool

Default on. Return an image of the size of the largest input image. If the images 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

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

nthreads : int

if fftw3 is installed, can specify the number of threads to allow FFTs to use. Probably only helpful for large arrays

fft_type : [None, ‘fftw’, ‘scipy’, ‘numpy’]

Which FFT implementation to use. If not specified, defaults to the type specified in the FFT_TYPE ConfigurationItem

See also

convolve
Convolve is a non-fft version of this code.

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], interpolate_nan=True)
array([ 1.,  0.,  3.])
>>> convolve_fft([1,np.nan,3],[0,1,0], interpolate_nan=True, min_wt=1e-8)
array([ 1.,  nan,  3.])
>>> convolve_fft([1,np.nan,3],[1,1,1], interpolate_nan=True)
array([ 1.,  4.,  3.])
>>> convolve_fft([1,np.nan,3],[1,1,1], interpolate_nan=True, normalize_kernel=True)
array([ 1.,  2.,  3.])

Page Contents