astropy:docs

convolve_fft

astropy.convolution.convolve_fft(array, kernel, boundary=u'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, allow_huge=False, fftn=<function fftn at 0x7f32c3302398>, ifftn=<function ifftn at 0x7f32c3302410>, complex_dtype=<type 'complex'>)[source] [edit on github]

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:

  • It can treat NaN values as zeros or interpolate over them.
  • inf values are treated as NaN
  • (optionally) It pads to the nearest 2^n size to improve FFT speed.
  • 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.fft and numpy.fft.ifft.
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’}, optional

A flag indicating how to handle boundaries:

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

interpolate_nan : bool, optional

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, optional

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, 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

normalize_kernel : function or boolean, 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.

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, optional

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

psf_pad : bool, optional

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, optional

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, optional

Return the fft(image)*fft(kernel) instead of the convolution (which is ifft(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, e.g. fftn=scipy.fftpack.fftn

complex_dtype : np.complex, optional

Which complex dtype to use. numpy has 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_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.

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, ignore_edge_zeros=True)
array([ 1.,  2.,  3.])
>>> import scipy.fftpack  # optional - requires scipy
>>> convolve_fft([1, np.nan, 3], [1, 1, 1], interpolate_nan=True,
...               normalize_kernel=True, ignore_edge_zeros=True,
...               fftn=scipy.fftpack.fft, ifftn=scipy.fftpack.ifft)
array([ 1.,  2.,  3.])

Page Contents