biosppy.signals
This sub-package provides methods to process common physiological signals (biosignals).
Modules
biosppy.signals.abp
This module provides methods to process Arterial Blood Pressure (ABP) signals.
- copyright:
2015-2018 by Instituto de Telecomunicacoes
- license:
BSD 3-clause, see LICENSE for more details.
- biosppy.signals.abp.abp(signal=None, sampling_rate=1000.0, show=True)[source]
Process a raw ABP signal and extract relevant signal features using default parameters.
- Parameters:
signal (array) – Raw ABP signal.
sampling_rate (int, float, optional) – Sampling frequency (Hz).
show (bool, optional) – If True, show a summary plot.
- Returns:
ts (array) – Signal time axis reference (seconds).
filtered (array) – Filtered ABP signal.
onsets (array) – Indices of ABP pulse onsets.
heart_rate_ts (array) – Heart rate time axis reference (seconds).
heart_rate (array) – Instantaneous heart rate (bpm).
- biosppy.signals.abp.find_onsets_zong2003(signal=None, sampling_rate=1000.0, sm_size=None, size=None, alpha=2.0, wrange=None, d1_th=0, d2_th=None)[source]
Determine onsets of ABP pulses. Skips corrupted signal parts. Based on the approach by Zong et al. [Zong03].
- Parameters:
signal (array) – Input filtered ABP signal.
sampling_rate (int, float, optional) – Sampling frequency (Hz).
sm_size (int, optional) – Size of smoother kernel (seconds). Defaults to 0.25
size (int, optional) – Window to search for maxima (seconds). Defaults to 5
alpha (float, optional) – Normalization parameter. Defaults to 2.0
wrange (int, optional) – The window in which to search for a peak (seconds). Defaults to 0.1
d1_th (int, optional) – Smallest allowed difference between maxima and minima. Defaults to 0
d2_th (int, optional) – Smallest allowed time between maxima and minima (seconds), Defaults to 0.15
- Returns:
onsets (array) – Indices of ABP pulse onsets.
References
[Zong03]W Zong, T Heldt, GB Moody and RG Mark, “An Open-source Algorithm to Detect Onset of Arterial Blood Pressure Pulses”, IEEE Comp. in Cardiology, vol. 30, pp. 259-262, 2003
biosppy.signals.acc
This module provides methods to process Acceleration (ACC) signals. Implemented code assumes ACC acquisition from a 3 orthogonal axis reference system.
- copyright:
2015-2023 by Instituto de Telecomunicacoes
- license:
BSD 3-clause, see LICENSE for more details.
Authors: Afonso Ferreira Diogo Vieira
- biosppy.signals.acc.acc(signal=None, sampling_rate=100.0, units=None, path=None, show=True, interactive=False)[source]
Process a raw ACC signal and extract relevant signal features using default parameters.
- Parameters:
signal (array) – Raw ACC signal.
sampling_rate (int, float, optional) – Sampling frequency (Hz).
units (str, optional) – The units of the input signal. If specified, the plot will have the y-axis labeled with the corresponding units.
path (str, optional) – If provided, the plot will be saved to the specified file.
show (bool, optional) – If True, show a summary plot.
interactive (bool, optional) – If True, shows an interactive plot.
- Returns:
ts (array) – Signal time axis reference (seconds).
signal (array) – Raw (unfiltered) ACC signal.
vm (array) – Vector Magnitude feature of the signal.
sm (array) – Signal Magnitude feature of the signal.
freq_features (dict) – Positive Frequency domains (Hz) of the signal.
amp_features (dict) – Normalized Absolute Amplitudes of the signal.
- biosppy.signals.acc.activity_index(signal=None, sampling_rate=100.0, window_1=5, window_2=60)[source]
Compute the activity index of an ACC signal. Follows the method described in [Lin18], the activity index is computed as follows: 1) Calculate the ACC magnitude if the signal is triaxial 2) Calculate the standard deviation of the ACC magnitude for each ‘window_1’ seconds 3) The activity index will be the mean standard deviation for each ‘window_2’ seconds
- Parameters:
signal (array) – Raw ACC signal.
sampling_rate (int, float, optional) – Sampling frequency (Hz).
window_1 (int, float, optional) – Window length (seconds) for the first moving average filter. Default: 5
window_2 (int, float, optional) – Window length (seconds) for the second moving average filter. Default: 60
- Returns:
ts (array) – Time axis reference (seconds).
activity_index (array) – Activity index of the signal.
References
[Lin18]W.-Y. Lin, V. Verma, M.-Y. Lee, C.-S. Lai, Activity
Monitoring with a Wrist-Worn, Accelerometer-Based Device, Micromachines 9 (2018) 450
- biosppy.signals.acc.frequency_domain_feature_extractor(signal=None, sampling_rate=100.0)[source]
Extracts the FFT from each ACC sub-signal (x, y, z), given the signal itself.
- Parameters:
signal (array) – Input ACC signal.
sampling_rate (int, float, optional) – Sampling frequency (Hz).
- Returns:
freq_features (dict) – Dictionary of positive frequencies (Hz) for all sub-signals.
amp_features (dict) – Dictionary of Normalized Absolute Amplitudes for all sub-signals.
- biosppy.signals.acc.time_domain_feature_extractor(signal=None)[source]
Extracts the vector magnitude and signal magnitude features from an input ACC signal, given the signal itself.
- Parameters:
signal (array) – Input ACC signal.
- Returns:
vm_features (array) – Extracted Vector Magnitude (VM) feature.
sm_features (array) – Extracted Signal Magnitude (SM) feature.
biosppy.signals.bvp
This module provides methods to process Blood Volume Pulse (BVP) signals.
——– DEPRECATED ——– PLEASE, USE THE PPG MODULE This module was left for compatibility —————————-
- copyright:
2015-2018 by Instituto de Telecomunicacoes
- license:
BSD 3-clause, see LICENSE for more details.
- biosppy.signals.bvp.bvp(signal=None, sampling_rate=1000.0, path=None, show=True)[source]
Process a raw BVP signal and extract relevant signal features using default parameters. :param signal: Raw BVP signal. :type signal: array :param sampling_rate: Sampling frequency (Hz). :type sampling_rate: int, float, optional :param path: If provided, the plot will be saved to the specified file. :type path: str, optional :param show: If True, show a summary plot. :type show: bool, optional
- Returns:
ts (array) – Signal time axis reference (seconds).
filtered (array) – Filtered BVP signal.
onsets (array) – Indices of BVP pulse onsets.
heart_rate_ts (array) – Heart rate time axis reference (seconds).
heart_rate (array) – Instantaneous heart rate (bpm).
biosppy.signals.ppg
This module provides methods to process Photoplethysmogram (PPG) signals.
- copyright:
2015-2018 by Instituto de Telecomunicacoes
- license:
BSD 3-clause, see LICENSE for more details.
- biosppy.signals.ppg.find_onsets_elgendi2013(signal=None, sampling_rate=1000.0, peakwindow=0.111, beatwindow=0.667, beatoffset=0.02, mindelay=0.3)[source]
Determines onsets of PPG pulses.
- Parameters:
signal (array) – Input filtered PPG signal.
sampling_rate (int, float, optional) – Sampling frequency (Hz).
peakwindow (float) – Parameter W1 on referenced article Optimized at 0.111
beatwindow (float) – Parameter W2 on referenced article Optimized at 0.667
beatoffset (float) – Parameter beta on referenced article Optimized at 0.2
mindelay (float) – Minimum delay between peaks. Avoids false positives
- Returns:
onsets (array) – Indices of PPG pulse onsets.
params (dict) – Input parameters of the function
References
Elgendi M, Norton I, Brearley M, Abbott D, Schuurmans D (2013) Systolic Peak Detection in
Acceleration Photoplethysmograms Measured from Emergency Responders in Tropical Conditions. PLoS ONE 8(10): e76585. doi:10.1371/journal.pone.0076585.
Notes
Optimal ranges for signal filtering (from Elgendi et al. 2013): “Optimization of the beat detector’s spectral window for the lower frequency resulted in a value within 0.5– 1 Hz with the higher frequency within 7–15 Hz”
All the number references below between curly brackets {…} by the code refer to the line numbers of code in “Table 2 Algorithm IV: DETECTOR (PPG signal, F1, F2, W1, W2, b)” from Elgendi et al. 2013 for a better comparison of the algorithm
- biosppy.signals.ppg.find_onsets_kavsaoglu2016(signal=None, sampling_rate=1000.0, alpha=0.2, k=4, init_bpm=90, min_delay=0.6, max_BPM=150)[source]
Determines onsets of PPG pulses.
- Parameters:
signal (array) – Input filtered PPG signal.
sampling_rate (int, float, optional) – Sampling frequency (Hz).
alpha (float, optional) – Low-pass filter factor. Avoids abrupt changes of BPM.
k (int, float, optional) – Number of segments by pulse. Width of each segment = Period of pulse according to current BPM / k
init_bpm (int, float, optional) – Initial BPM. Higher value results in a smaller segment width.
min_delay (float) – Minimum delay between peaks as percentage of current BPM pulse period. Avoids false positives
max_bpm (int, float, optional) – Maximum BPM. Maximum value accepted as valid BPM.
- Returns:
onsets (array) – Indices of PPG pulse onsets.
window_marks (array) – Indices of segments window boundaries.
params (dict) – Input parameters of the function
References
Kavsaoğlu, Ahmet & Polat, Kemal & Bozkurt, Mehmet. (2016). An innovative peak detection algorithm for
photoplethysmography signals: An adaptive segmentation method. TURKISH JOURNAL OF ELECTRICAL ENGINEERING & COMPUTER SCIENCES. 24. 1782-1796. 10.3906/elk-1310-177.
Notes
This algorithm is an adaption of the one described on Kavsaoğlu et al. (2016). This version takes into account a minimum delay between peaks and builds upon the adaptive segmentation by using a low-pass filter for BPM changes. This way, even if the algorithm wrongly detects a peak, the BPM value will stay relatively constant so the next pulse can be correctly segmented.
- biosppy.signals.ppg.ppg(signal=None, sampling_rate=1000.0, units=None, show=True)[source]
Process a raw PPG signal and extract relevant signal features using default parameters.
- Parameters:
signal (array) – Raw PPG signal.
sampling_rate (int, float, optional) – Sampling frequency (Hz).
units (str, optional) – The units of the input signal. If specified, the plot will have the y-axis labeled with the corresponding units.
show (bool, optional) – If True, show a summary plot.
- Returns:
ts (array) – Signal time axis reference (seconds).
filtered (array) – Filtered PPG signal.
peaks (array) – Indices of PPG pulse peaks.
templates_ts (array) – Templates time axis reference (seconds).
templates (array) – Extracted heartbeat templates.
heart_rate_ts (array) – Heart rate time axis reference (seconds).
heart_rate (array) – Instantaneous heart rate (bpm).
- biosppy.signals.ppg.ppg_segmentation(signal=None, sampling_rate=1000.0, peaks=None, selection=False, peak_threshold=None)[source]
Segments a filtered PPG signal. Segmentation filtering is achieved by taking into account segments selected by peak height and pulse morphology.
- Parameters:
signal (array) – Filtered PPG signal.
sampling_rate (int, float, optional) – Sampling frequency (Hz).
peaks (array) – List of PPG systolic peaks.
selection (bool, optional) – If True, performs selection with peak height and pulse morphology.
peak_threshold (int, float, optional) – If selection is True, selects peaks with height greater than defined threshold.
- Returns:
onsets (array) – Indices of PPG pulse onsets (i.e., start of beats) of the selected segments.
peaks (array) – List of PPG systolic peaks of the selected segments.
segments_loc (array) – Start and end indices for each selected pulse segment.
biosppy.signals.tools
This module provides various signal analysis methods in the time and frequency domains.
- copyright:
2015-2023 by Instituto de Telecomunicacoes
- license:
BSD 3-clause, see LICENSE for more details.
- class biosppy.signals.tools.OnlineFilter(b=None, a=None)[source]
Bases:
objectOnline filtering.
- Parameters:
b (array) – Numerator coefficients.
a (array) – Denominator coefficients.
- biosppy.signals.tools.analytic_signal(signal=None, N=None)[source]
Compute analytic signal, using the Hilbert Transform.
- Parameters:
signal (array) – Input signal.
N (int, optional) – Number of Fourier components; default is len(signal).
- Returns:
amplitude (array) – Amplitude envelope of the analytic signal.
phase (array) – Instantaneous phase component of the analystic signal.
- biosppy.signals.tools.band_power(freqs=None, power=None, frequency=None, decibel=True)[source]
Compute the avearge power in a frequency band.
- Parameters:
freqs (array) – Array of frequencies (Hz) at which the power was computed.
power (array) – Input power spectrum.
frequency (list, array) – Pair of frequencies defining the band.
decibel (bool, optional) – If True, input power is in decibels.
- Returns:
avg_power (float) – The average power in the band.
- biosppy.signals.tools.detrend_smoothness_priors(signal, smoothing_factor=10)[source]
Detrending method based on smoothness priors applied to HRV signal analysis.
Follows the approach by Tarvainen et al. [Tarivainen2002].
- Parameters:
signal (array) – Signal to be detrended.
smoothing_factor (int, float, optional) – Smoothing parameter lambda. Default: 10.
- Returns:
detrended (array) – Detrended signal.
trend (array) – The trend component of the signal.
References
[Tarivainen2002]M. P. Tarvainen, P. O. Ranta-aho and P. A. Karjalainen, “An advanced detrending method with application to HRV analysis,” in IEEE Transactions on Biomedical Engineering, vol. 49, no. 2, pp. 172-175, Feb. 2002, doi: 10.1109/10.979357.
- biosppy.signals.tools.distance_profile(query=None, signal=None, metric='euclidean')[source]
Compute the distance profile of a query sequence against a signal.
Implements the algorithm described in [Mueen2014].
- Parameters:
query (array) – Input query signal sequence.
signal (array) – Input target time series signal.
metric (str, optional) – The distance metric to use; one of ‘euclidean’ or ‘pearson’; default is ‘euclidean’.
- Returns:
dist (array) – Distance of the query sequence to every sub-sequnce in the signal.
Notes
Computes distances on z-normalized data.
References
[Mueen2014]Abdullah Mueen, Hossein Hamooni, “Trilce Estrada: Time Series Join on Subsequence Correlation”, ICDM 2014: 450-459
- biosppy.signals.tools.filter_signal(signal=None, ftype='FIR', band='lowpass', order=None, frequency=None, sampling_rate=1000.0, **kwargs)[source]
Filter a signal according to the given parameters.
- Parameters:
signal (array) – Signal to filter.
ftype (str) –
- Filter type:
Finite Impulse Response filter (‘FIR’);
Butterworth filter (‘butter’);
Chebyshev filters (‘cheby1’, ‘cheby2’);
Elliptic filter (‘ellip’);
Bessel filter (‘bessel’).
Notch filter (‘notch’).
band (str) –
- Band type:
Low-pass filter (‘lowpass’);
High-pass filter (‘highpass’);
Band-pass filter (‘bandpass’);
Band-stop filter (‘bandstop’).
order (int) – Order of the filter.
frequency (int, float, list, array) –
- Cutoff frequencies; format depends on type of band:
’lowpass’ or ‘bandpass’: single frequency;
’bandpass’ or ‘bandstop’: pair of frequencies.
sampling_rate (int, float, optional) – Sampling frequency (Hz).
**kwargs (dict, optional) –
Additional keyword arguments are passed to the underlying scipy.signal function. - Q : float
Quality factor (only for ‘notch’ filter). Default: 30.
- Returns:
signal (array) – Filtered signal.
sampling_rate (float) – Sampling frequency (Hz).
params (dict) – Filter parameters.
Notes
Uses a forward-backward filter implementation. Therefore, the combined filter has linear phase.
- biosppy.signals.tools.find_extrema(signal=None, mode='both')[source]
Locate local extrema points in a signal.
Based on Fermat’s Theorem [Ferm].
- Parameters:
signal (array) – Input signal.
mode (str, optional) – Whether to find maxima (‘max’), minima (‘min’), or both (‘both’).
- Returns:
extrema (array) – Indices of the extrema points.
values (array) – Signal values at the extrema points.
References
[Ferm]Wikipedia, “Fermat’s theorem (stationary points)”, https://en.wikipedia.org/wiki/Fermat%27s_theorem_(stationary_points)
- biosppy.signals.tools.find_intersection(x1=None, y1=None, x2=None, y2=None, alpha=1.5, xtol=1e-06, ytol=1e-06)[source]
Find the intersection points between two lines using piecewise polynomial interpolation.
- Parameters:
x1 (array) – Array of x-coordinates of the first line.
y1 (array) – Array of y-coordinates of the first line.
x2 (array) – Array of x-coordinates of the second line.
y2 (array) – Array of y-coordinates of the second line.
alpha (float, optional) – Resolution factor for the x-axis; fraction of total number of x-coordinates.
xtol (float, optional) – Tolerance for the x-axis.
ytol (float, optional) – Tolerance for the y-axis.
- Returns:
roots (array) – Array of x-coordinates of found intersection points.
values (array) – Array of y-coordinates of found intersection points.
Notes
If no intersection is found, returns the closest point.
- biosppy.signals.tools.finite_difference(signal=None, weights=None)[source]
Apply the Finite Difference method to compute derivatives.
- Parameters:
signal (array) – Signal to differentiate.
weights (list, array) – Finite difference weight coefficients.
- Returns:
index (array) – Indices from signal for which the derivative was computed.
derivative (array) – Computed derivative.
Notes
The method assumes central differences weights.
The method accounts for the delay introduced by the algorithm.
- Raises:
ValueError – If the number of weights is not odd.
- biosppy.signals.tools.get_filter(ftype='FIR', band='lowpass', order=None, frequency=None, sampling_rate=1000.0, **kwargs)[source]
Compute digital (FIR or IIR) filter coefficients with the given parameters.
- Parameters:
ftype (str) –
- Filter type:
Finite Impulse Response filter (‘FIR’);
Butterworth filter (‘butter’);
Chebyshev filters (‘cheby1’, ‘cheby2’);
Elliptic filter (‘ellip’);
Bessel filter (‘bessel’).
Notch filter (‘notch’).
band (str) –
- Band type:
Low-pass filter (‘lowpass’);
High-pass filter (‘highpass’);
Band-pass filter (‘bandpass’);
Band-stop filter (‘bandstop’).
order (int) – Order of the filter.
frequency (int, float, list, array) –
- Cutoff frequencies; format depends on type of band:
’lowpass’ or ‘highpass’: single frequency;
’bandpass’ or ‘bandstop’: pair of frequencies.
sampling_rate (int, float, optional) – Sampling frequency (Hz).
**kwargs (dict, optional) –
Additional keyword arguments are passed to the underlying scipy.signal function. - Q : float
Quality factor (only for ‘notch’ filter). Default: 30.
- Returns:
b (array) – Numerator coefficients.
a (array) – Denominator coefficients.
See Also – scipy.signal
- biosppy.signals.tools.get_heart_rate(beats=None, sampling_rate=1000.0, smooth=False, size=3)[source]
Compute instantaneous heart rate from an array of beat indices.
- Parameters:
beats (array) – Beat location indices.
sampling_rate (int, float, optional) – Sampling frequency (Hz).
smooth (bool, optional) – If True, perform smoothing on the resulting heart rate.
size (int, optional) – Size of smoothing window; ignored if smooth is False.
- Returns:
index (array) – Heart rate location indices.
heart_rate (array) – Instantaneous heart rate (bpm).
Notes
Assumes normal human heart rate to be between 40 and 200 bpm.
- biosppy.signals.tools.mean_waves(data=None, size=None, step=None)[source]
Extract mean samples from a data set.
- Parameters:
data (array) – An m by n array of m data samples in an n-dimensional space.
size (int) – Number of samples to use for each mean sample.
step (int, optional) – Number of samples to jump, controlling overlap; default is equal to size (no overlap).
- Returns:
waves (array) – An k by n array of mean samples.
Notes
Discards trailing samples if they are not enough to satify the size parameter.
- Raises:
ValueError – If step is an invalid value.
ValueError – If there are not enough samples for the given size.
- biosppy.signals.tools.median_waves(data=None, size=None, step=None)[source]
Extract median samples from a data set.
- Parameters:
data (array) – An m by n array of m data samples in an n-dimensional space.
size (int) – Number of samples to use for each median sample.
step (int, optional) – Number of samples to jump, controlling overlap; default is equal to size (no overlap).
- Returns:
waves (array) – An k by n array of median samples.
Notes
Discards trailing samples if they are not enough to satify the size parameter.
- Raises:
ValueError – If step is an invalid value.
ValueError – If there are not enough samples for the given size.
- biosppy.signals.tools.normalize(signal=None, ddof=1)[source]
Normalize a signal to zero mean and unitary standard deviation.
- Parameters:
signal (array) – Input signal.
ddof (int, optional) – Delta degrees of freedom for standard deviation computation; the divisor is N - ddof, where N is the number of elements; default is one.
- Returns:
signal (array) – Normalized signal.
- biosppy.signals.tools.pearson_correlation(x=None, y=None)[source]
Compute the Pearson Correlation Coefficient bertween two signals.
The coefficient is given by:
![r_{xy} = \frac{E[(X - \mu_X) (Y - \mu_Y)]}{\sigma_X \sigma_Y}](_images/math/25984ac0b99af0b57d7d197c63928096548adcb5.png)
- Parameters:
x (array) – First input signal.
y (array) – Second input signal.
- Returns:
rxy (float) – Pearson correlation coefficient, ranging between -1 and +1.
- Raises:
ValueError – If the input signals do not have the same length.
- biosppy.signals.tools.phase_locking(signal1=None, signal2=None, N=None)[source]
Compute the Phase-Locking Factor (PLF) between two signals.
- Parameters:
signal1 (array) – First input signal.
signal2 (array) – Second input signal.
N (int, optional) – Number of Fourier components.
- Returns:
plf (float) – The PLF between the two signals.
- biosppy.signals.tools.power_spectrum(signal=None, sampling_rate=1000.0, pad=None, pow2=False, decibel=True)[source]
Compute the power spectrum of a signal (one-sided).
- Parameters:
signal (array) – Input signal.
sampling_rate (int, float, optional) – Sampling frequency (Hz).
pad (int, optional) – Padding for the Fourier Transform (number of zeros added); defaults to no padding..
pow2 (bool, optional) – If True, rounds the number of points N = len(signal) + pad to the nearest power of 2 greater than N.
decibel (bool, optional) – If True, returns the power in decibels.
- Returns:
freqs (array) – Array of frequencies (Hz) at which the power was computed.
power (array) – Power spectrum.
- biosppy.signals.tools.rms_error(x=None, y=None)[source]
Compute the Root-Mean-Square Error between two signals.
The error is given by:
![rmse = \sqrt{E[(X - Y)^2]}](_images/math/54caf68afcbb74e741a359c4a7afb94b23e61956.png)
- Parameters:
x (array) – First input signal.
y (array) – Second input signal.
- Returns:
rmse (float) – Root-mean-square error.
- Raises:
ValueError – If the input signals do not have the same length.
- biosppy.signals.tools.signal_cross_join(signal1=None, signal2=None, size=None, index=None, limit=None)[source]
Compute the matrix profile for a similarity join of two time series.
Computes the nearest sub-sequence in signal2 for each sub-sequence in signal1. Implements the algorithm described in [Yeh2016_c].
- Parameters:
signal1 (array) – Fisrt input time series signal.
signal2 (array) – Second input time series signal.
size (int) – Size of the query sub-sequences.
index (list, array, optional) – Starting indices for query sub-sequences; the default is to search all sub-sequences.
limit (int, optional) – Upper limit for the number of query sub-sequences; the default is to search all sub-sequences.
- Returns:
matrix_index (array) – Matric profile index.
matrix_profile (array) – Computed matrix profile (distances).
Notes
Computes euclidean distances on z-normalized data.
References
[Yeh2016_c]Chin-Chia Michael Yeh, Yan Zhu, Liudmila Ulanova, Nurjahan Begum, Yifei Ding, Hoang Anh Dau, Diego Furtado Silva, Abdullah Mueen, Eamonn Keogh, “Matrix Profile I: All Pairs Similarity Joins for Time Series: A Unifying View that Includes Motifs, Discords and Shapelets”, IEEE ICDM 2016
- biosppy.signals.tools.signal_self_join(signal=None, size=None, index=None, limit=None)[source]
Compute the matrix profile for a self-similarity join of a time series.
Implements the algorithm described in [Yeh2016_b].
- Parameters:
signal (array) – Input target time series signal.
size (int) – Size of the query sub-sequences.
index (list, array, optional) – Starting indices for query sub-sequences; the default is to search all sub-sequences.
limit (int, optional) – Upper limit for the number of query sub-sequences; the default is to search all sub-sequences.
- Returns:
matrix_index (array) – Matric profile index.
matrix_profile (array) – Computed matrix profile (distances).
Notes
Computes euclidean distances on z-normalized data.
References
[Yeh2016_b]Chin-Chia Michael Yeh, Yan Zhu, Liudmila Ulanova, Nurjahan Begum, Yifei Ding, Hoang Anh Dau, Diego Furtado Silva, Abdullah Mueen, Eamonn Keogh, “Matrix Profile I: All Pairs Similarity Joins for Time Series: A Unifying View that Includes Motifs, Discords and Shapelets”, IEEE ICDM 2016
- biosppy.signals.tools.signal_stats(signal=None)[source]
Compute various metrics describing the signal.
- Parameters:
signal (array) – Input signal.
- Returns:
mean (float) – Mean of the signal.
median (float) – Median of the signal.
min (float) – Minimum signal value.
max (float) – Maximum signal value.
max_amp (float) – Maximum absolute signal amplitude, in relation to the mean.
range (float) – Signal range (max - min).
q1 (float) – First quartile of the signal.
q3 (float) – Third quartile of the signal.
var (float) – Signal variance (unbiased).
std_dev (float) – Standard signal deviation (unbiased).
abs_dev (float) – Mean absolute signal deviation around the median.
rms (float) – Root-mean-square of the signal.
kurtosis (float) – Signal kurtosis (unbiased).
skew (float) – Signal skewness (unbiased).
- biosppy.signals.tools.smoother(signal=None, kernel='boxzen', size=10, mirror=True, **kwargs)[source]
Smooth a signal using an N-point moving average [MAvg] filter.
This implementation uses the convolution of a filter kernel with the input signal to compute the smoothed signal [Smit97].
Availabel kernels: median, boxzen, boxcar, triang, blackman, hamming, hann, bartlett, flattop, parzen, bohman, blackmanharris, nuttall, barthann, kaiser (needs beta), gaussian (needs std), general_gaussian (needs power, width), slepian (needs width), chebwin (needs attenuation).
- Parameters:
signal (array) – Signal to smooth.
kernel (str, array, optional) – Type of kernel to use; if array, use directly as the kernel.
size (int, optional) – Size of the kernel; ignored if kernel is an array.
mirror (bool, optional) – If True, signal edges are extended to avoid boundary effects.
**kwargs (dict, optional) – Additional keyword arguments are passed to the underlying scipy.signal.windows function.
- Returns:
signal (array) – Smoothed signal.
params (dict) – Smoother parameters.
Notes
When the kernel is ‘median’, mirror is ignored.
References
[MAvg]Wikipedia, “Moving Average”, http://en.wikipedia.org/wiki/Moving_average
[Smit97]S. W. Smith, “Moving Average Filters - Implementation by Convolution”, http://www.dspguide.com/ch15/1.htm, 1997
- biosppy.signals.tools.synchronize(x=None, y=None, detrend=True)[source]
Align two signals based on cross-correlation.
- Parameters:
x (array) – First input signal.
y (array) – Second input signal.
detrend (bool, optional) – If True, remove signal means before computation.
- Returns:
delay (int) – Delay (number of samples) of ‘x’ in relation to ‘y’; if ‘delay’ < 0 , ‘x’ is ahead in relation to ‘y’; if ‘delay’ > 0 , ‘x’ is delayed in relation to ‘y’.
corr (float) – Value of maximum correlation.
synch_x (array) – Biggest possible portion of ‘x’ in synchronization.
synch_y (array) – Biggest possible portion of ‘y’ in synchronization.
- biosppy.signals.tools.welch_spectrum(signal=None, sampling_rate=1000.0, size=None, overlap=None, window='hanning', window_kwargs=None, pad=None, decibel=True)[source]
Compute the power spectrum of a signal using Welch’s method (one-sided).
- Parameters:
signal (array) – Input signal.
sampling_rate (int, float, optional) – Sampling frequency (Hz).
size (int, optional) – Number of points in each Welch segment; defaults to the equivalent of 1 second; ignored when ‘window’ is an array.
overlap (int, optional) – Number of points to overlap between segments; defaults to size / 2.
window (str, array, optional) – Type of window to use.
window_kwargs (dict, optional) – Additional keyword arguments to pass on window creation; ignored if ‘window’ is an array.
pad (int, optional) – Padding for the Fourier Transform (number of zeros added); defaults to no padding.
decibel (bool, optional) – If True, returns the power in decibels.
- Returns:
freqs (array) – Array of frequencies (Hz) at which the power was computed.
power (array) – Power spectrum.
Notes
Detrends each Welch segment by removing the mean.
- biosppy.signals.tools.windower(signal=None, size=None, step=None, fcn=None, fcn_kwargs=None, kernel='boxcar', kernel_kwargs=None)[source]
Apply a function to a signal in sequential windows, with optional overlap.
Availabel window kernels: boxcar, triang, blackman, hamming, hann, bartlett, flattop, parzen, bohman, blackmanharris, nuttall, barthann, kaiser (needs beta), gaussian (needs std), general_gaussian (needs power, width), slepian (needs width), chebwin (needs attenuation).
- Parameters:
signal (array) – Input signal.
size (int) – Size of the signal window.
step (int, optional) – Size of window shift; if None, there is no overlap.
fcn (callable) – Function to apply to each window.
fcn_kwargs (dict, optional) – Additional keyword arguments to pass to ‘fcn’.
kernel (str, array, optional) – Type of kernel to use; if array, use directly as the kernel.
kernel_kwargs (dict, optional) – Additional keyword arguments to pass on window creation; ignored if ‘kernel’ is an array.
- Returns:
index (array) – Indices characterizing window locations (start of the window).
values (array) – Concatenated output of calling ‘fcn’ on each window.
- biosppy.signals.tools.zero_cross(signal=None, detrend=False)[source]
Locate the indices where the signal crosses zero.
- Parameters:
signal (array) – Input signal.
detrend (bool, optional) – If True, remove signal mean before computation.
- Returns:
zeros (array) – Indices of zero crossings.
Notes
When the signal crosses zero between samples, the first index is returned.