find_peaks#
- photutils.detection.find_peaks(data, threshold, *, box_size=3, footprint=None, mask=None, border_width=None, npeaks=inf, centroid_func=None, error=None, wcs=None)[source]#
Find local peaks in an image that are above a specified threshold value.
Peaks are the maxima above the
threshold
within a local region. The local regions are defined by either thebox_size
orfootprint
parameters.box_size
defines the local region around each pixel as a square box.footprint
is a boolean array whereTrue
values specify the region shape.If multiple pixels within a local region have identical intensities, then the coordinates of all such pixels are returned. Otherwise, there will be only one peak pixel per local region. Thus, the defined region effectively imposes a minimum separation between peaks unless there are identical peaks within the region.
If
centroid_func
is input, then it will be used to calculate a centroid within the defined local region centered on each detected peak pixel. In this case, the centroid will also be returned in the output table.- Parameters:
- dataarray_like
The 2D array of the image.
- thresholdfloat, scalar
Quantity
or array_like The data value or pixel-wise data values to be used for the detection threshold. If
data
is aQuantity
array, thenthreshold
must have the same units asdata
. A 2Dthreshold
must have the same shape asdata
. Seedetect_threshold
for one way to create athreshold
image.- box_sizescalar or tuple, optional
The size of the local region to search for peaks at every point in
data
. Ifbox_size
is a scalar, then the region shape will be(box_size, box_size)
. Eitherbox_size
orfootprint
must be defined. If they are both defined, thenfootprint
overridesbox_size
.- footprint
ndarray
of bools, optional A boolean array where
True
values describe the local footprint region within which to search for peaks at every point indata
.box_size=(n, m)
is equivalent tofootprint=np.ones((n, m))
. Eitherbox_size
orfootprint
must be defined. If they are both defined, thenfootprint
overridesbox_size
.- maskarray_like, bool, optional
A boolean mask with the same shape as
data
, where aTrue
value indicates the corresponding element ofdata
is masked.- border_widthint, array_like of int, or None, optional
The width in pixels to exclude around the border of the
data
. Ifborder_width
is a scalar thenborder_width
will be applied to all sides. Ifborder_width
has two elements, they must be in(ny, nx)
order. IfNone
, then no border is excluded. The border width values must be non-negative integers.- npeaksint, optional
The maximum number of peaks to return. When the number of detected peaks exceeds
npeaks
, the peaks with the highest peak intensities will be returned.- centroid_funccallable, optional
A callable object (e.g., function or class) that is used to calculate the centroid of a 2D array. The
centroid_func
must accept a 2Dndarray
, have amask
keyword, and optionally anerror
keyword. The callable object must return a tuple of two 1Dndarray
objects, representing the x and y centroids, respectively.- errorarray_like, optional
The 2D array of the 1-sigma errors of the input
data
.error
is used only ifcentroid_func
is input (theerror
array is passed directly to thecentroid_func
). Ifdata
is aQuantity
array, thenerror
must have the same units asdata
.- wcs
None
or WCS object, optional A world coordinate system (WCS) transformation that supports the astropy shared interface for WCS (e.g.,
astropy.wcs.WCS
,gwcs.wcs.WCS
). IfNone
, then the sky coordinates will not be returned in the outputTable
.
- Returns:
Notes
By default, the returned pixel coordinates are the integer indices of the maximum pixel value within the input
box_size
orfootprint
(i.e., only the peak pixel is identified). However, a centroiding function can be input via thecentroid_func
keyword to compute centroid coordinates with subpixel precision within the inputbox_size
orfootprint
.