bundles / scipy 1.17.1 / scipy / signal / _peak_finding / find_peaks
function
scipy.signal._peak_finding:find_peaks
Signature
def find_peaks ( x , height = None , threshold = None , distance = None , prominence = None , width = None , wlen = None , rel_height = 0.5 , plateau_size = None ) Summary
Find peaks inside a signal based on peak properties.
Extended Summary
This function takes a 1-D array and finds all local maxima by simple comparison of neighboring values. Optionally, a subset of these peaks can be selected by specifying conditions for a peak's properties.
Parameters
x: sequenceA signal with peaks.
height: number or ndarray or sequence, optionalRequired height of peaks. Either a number,
None, an array matchingxor a 2-element sequence of the former. The first element is always interpreted as the minimal and the second, if supplied, as the maximal required height.threshold: number or ndarray or sequence, optionalRequired threshold of peaks, the vertical distance to its neighboring samples. Either a number,
None, an array matchingxor a 2-element sequence of the former. The first element is always interpreted as the minimal and the second, if supplied, as the maximal required threshold.distance: number, optionalRequired minimal horizontal distance (>= 1) in samples between neighbouring peaks. Smaller peaks are removed first until the condition is fulfilled for all remaining peaks.
prominence: number or ndarray or sequence, optionalRequired prominence of peaks. Either a number,
None, an array matchingxor a 2-element sequence of the former. The first element is always interpreted as the minimal and the second, if supplied, as the maximal required prominence.width: number or ndarray or sequence, optionalRequired width of peaks in samples. Either a number,
None, an array matchingxor a 2-element sequence of the former. The first element is always interpreted as the minimal and the second, if supplied, as the maximal required width.wlen: int, optionalUsed for calculation of the peaks prominences, thus it is only used if one of the arguments
prominenceorwidthis given. See argumentwlenin peak_prominences for a full description of its effects.rel_height: float, optionalUsed for calculation of the peaks width, thus it is only used if
widthis given. See argumentrel_heightin peak_widths for a full description of its effects.plateau_size: number or ndarray or sequence, optionalRequired size of the flat top of peaks in samples. Either a number,
None, an array matchingxor a 2-element sequence of the former. The first element is always interpreted as the minimal and the second, if supplied as the maximal required plateau size.
Returns
peaks: ndarrayIndices of peaks in
xthat satisfy all given conditions.properties: dictA dictionary containing properties of the returned peaks which were calculated as intermediate results during evaluation of the specified conditions:
'peak_heights'
If
heightis given, the height of each peak inx.
'left_thresholds', 'right_thresholds'
If
thresholdis given, these keys contain a peaks vertical distance to its neighbouring samples.
'prominences', 'right_bases', 'left_bases'
If
prominenceis given, these keys are accessible. See peak_prominences for a description of their content.
'widths', 'width_heights', 'left_ips', 'right_ips'
If
widthis given, these keys are accessible. See peak_widths for a description of their content.
'plateau_sizes', left_edges', 'right_edges'
If
plateau_sizeis given, these keys are accessible and contain the indices of a peak's edges (edges are still part of the plateau) and the calculated plateau sizes.
To calculate and return properties without excluding peaks, provide the open interval
(None, None)as a value to the appropriate argument (excludingdistance).
Warns
: PeakPropertyWarningRaised if a peak's properties have unexpected values (see peak_prominences and peak_widths).
Warnings
This function may return unexpected results for data containing NaNs. To avoid this, NaNs should either be removed or replaced.
Notes
In the context of this function, a peak or local maximum is defined as any sample whose two direct neighbours have a smaller amplitude. For flat peaks (more than one sample of equal amplitude wide) the index of the middle sample is returned (rounded down in case the number of samples is even). For noisy signals the peak locations can be off because the noise might change the position of local maxima. In those cases consider smoothing the signal before searching for peaks or use other peak finding and fitting methods (like find_peaks_cwt).
Some additional comments on specifying conditions:
Almost all conditions (excluding
distance) can be given as half-open or closed intervals, e.g.,1or(1, None)defines the half-open interval while(None, 1)defines the interval . The open interval(None, None)can be specified as well, which returns the matching properties without exclusion of peaks.The border is always included in the interval used to select valid peaks.
For several conditions the interval borders can be specified with arrays matching
xin shape which enables dynamic constrains based on the sample position.The conditions are evaluated in the following order:
plateau_size,height,threshold,distance,prominence,width. In most cases this order is the fastest one because faster operations are applied first to reduce the number of peaks that need to be evaluated later.While indices in peaks are guaranteed to be at least
distancesamples apart, edges of flat peaks may be closer than the alloweddistance.Use
wlento reduce the time it takes to evaluate the conditions forprominenceorwidthifxis large or has many local maxima (see peak_prominences).
Array API Standard Support
find_peaks has experimental support for Python Array API Standard compatible backends in addition to NumPy. Please consider testing these features by setting an environment variable SCIPY_ARRAY_API=1 and providing CuPy, PyTorch, JAX, or Dask arrays as array arguments. The following combinations of backend and device (or other capability) are supported.
==================== ==================== ==================== Library CPU GPU ==================== ==================== ==================== NumPy ✅ n/a CuPy n/a ⛔ PyTorch ⛔ ⛔ JAX ⛔ ⛔ Dask ⛔ n/a ==================== ==================== ====================
See
dev-arrayapifor more information.
Examples
To demonstrate this function's usage we use a signal `x` supplied with SciPy (see `scipy.datasets.electrocardiogram`). Let's find all peaks (local maxima) in `x` whose amplitude lies above 0.import numpy as np import matplotlib.pyplot as plt from scipy.datasets import electrocardiogram from scipy.signal import find_peaks✓
x = electrocardiogram()[2000:4000] peaks, _ = find_peaks(x, height=0) plt.plot(x) plt.plot(peaks, x[peaks], "x") plt.plot(np.zeros_like(x), "--", color="gray")⚠
plt.show()
✓border = np.sin(np.linspace(0, 3 * np.pi, x.size)) peaks, _ = find_peaks(x, height=(-border, border)) plt.plot(x) plt.plot(-border, "--", color="gray") plt.plot(border, ":", color="gray") plt.plot(peaks, x[peaks], "x")⚠
plt.show()
✓peaks, _ = find_peaks(x, distance=150) np.diff(peaks) plt.plot(x) plt.plot(peaks, x[peaks], "x")⚠
plt.show()
✓peaks, properties = find_peaks(x, prominence=(None, 0.6)) properties["prominences"].max() plt.plot(x) plt.plot(peaks, x[peaks], "x")⚠
plt.show()
✓x = electrocardiogram()[17000:18000] peaks, properties = find_peaks(x, prominence=1, width=20) properties["prominences"], properties["widths"] plt.plot(x) plt.plot(peaks, x[peaks], "x") plt.vlines(x=peaks, ymin=x[peaks] - properties["prominences"], ymax = x[peaks], color = "C1") plt.hlines(y=properties["width_heights"], xmin=properties["left_ips"], xmax=properties["right_ips"], color = "C1")⚠
plt.show()
✓See also
- find_peaks_cwt
Find peaks using the wavelet transformation.
- peak_prominences
Directly calculate the prominence of peaks.
- peak_widths
Directly calculate the width of peaks.
Aliases
-
scipy.signal.find_peaks