SciPy

scipy.signal.find_peaks_cwt

scipy.signal.find_peaks_cwt(vector, widths, wavelet=None, max_distances=None, gap_thresh=None, min_length=None, min_snr=1, noise_perc=10)[source]

Attempt to find the peaks in a 1-D array.

The general approach is to smooth vector by convolving it with wavelet(width) for each width in widths. Relative maxima which appear at enough length scales, and with sufficiently high SNR, are accepted.

Parameters:

vector : ndarray

1-D array in which to find the peaks.

widths : sequence

1-D array of widths to use for calculating the CWT matrix. In general, this range should cover the expected width of peaks of interest.

wavelet : callable, optional

Should take a single variable and return a 1-D array to convolve with vector. Should be normalized to unit area. Default is the ricker wavelet.

max_distances : ndarray, optional

At each row, a ridge line is only connected if the relative max at row[n] is within max_distances[n] from the relative max at row[n+1]. Default value is widths/4.

gap_thresh : float, optional

If a relative maximum is not found within max_distances, there will be a gap. A ridge line is discontinued if there are more than gap_thresh points without connecting a new relative maximum. Default is 2.

min_length : int, optional

Minimum length a ridge line needs to be acceptable. Default is cwt.shape[0] / 4, ie 1/4-th the number of widths.

min_snr : float, optional

Minimum SNR ratio. Default 1. The signal is the value of the cwt matrix at the shortest length scale (cwt[0, loc]), the noise is the noise_perc`th percentile of datapoints contained within a window of `window_size around cwt[0, loc].

noise_perc : float, optional

When calculating the noise floor, percentile of data points examined below which to consider noise. Calculated using stats.scoreatpercentile. Default is 10.

Returns:

peaks_indices : list

Indices of the locations in the vector where peaks were found. The list is sorted.

See also

cwt

Notes

This approach was designed for finding sharp peaks among noisy data, however with proper parameter selection it should function well for different peak shapes.

The algorithm is as follows:
  1. Perform a continuous wavelet transform on vector, for the supplied widths. This is a convolution of vector with wavelet(width) for each width in widths. See cwt
  2. Identify “ridge lines” in the cwt matrix. These are relative maxima at each row, connected across adjacent rows. See identify_ridge_lines
  3. Filter the ridge_lines using filter_ridge_lines.

New in version 0.11.0.

References

[R174]Bioinformatics (2006) 22 (17): 2059-2065. doi: 10.1093/bioinformatics/btl355 http://bioinformatics.oxfordjournals.org/content/22/17/2059.long

Examples

>>> from scipy import signal
>>> xs = np.arange(0, np.pi, 0.05)
>>> data = np.sin(xs)
>>> peakind = signal.find_peaks_cwt(data, np.arange(1,10))
>>> peakind, xs[peakind], data[peakind]
([32], array([ 1.6]), array([ 0.9995736]))

Previous topic

scipy.signal.cwt

Next topic

scipy.signal.argrelmin