biosppy.signals

This sub-package provides methods to process common physiological signals (biosignals).

Modules

biosppy.signals.bvp

This module provides methods to process Blood Volume Pulse (BVP) signals.

copyright:
  1. 2015 by Instituto de Telecomunicacoes
license:

BSD 3-clause, see LICENSE for more details.

biosppy.signals.bvp.bvp(signal=None, sampling_rate=1000.0, show=True)

Process a raw BVP signal and extract relevant signal features using default parameters.

Parameters:
  • signal (array) – Raw BVP 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 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.bvp.find_onsets(signal=None, sampling_rate=1000.0)

Determine onsets of BVP pulses.

Skips corrupted signal parts.

Parameters:
  • signal (array) – Input filtered BVP signal.
  • sampling_rate (int, float, optional) – Sampling frequency (Hz).
Returns:

onsets (array) – Indices of BVP pulse onsets.

biosppy.signals.ecg

This module provides methods to process Electrocardiographic (ECG) signals. Implemented code assumes a single-channel Lead I like ECG signal.

copyright:
  1. 2015 by Instituto de Telecomunicacoes
license:

BSD 3-clause, see LICENSE for more details.

biosppy.signals.ecg.christov_segmenter(signal=None, sampling_rate=1000.0)

ECG R-peak segmentation algorithm.

Follows the approach by Christov [Chri04].

Parameters:
  • signal (array) – Input filtered ECG signal.
  • sampling_rate (int, float, optional) – Sampling frequency (Hz).
Returns:

rpeaks (array) – R-peak location indices.

References

[Chri04]Ivaylo I. Christov, “Real time electrocardiogram QRS detection using combined adaptive threshold”, BioMedical Engineering OnLine 2004, vol. 3:28, 2004
biosppy.signals.ecg.compare_segmentation(reference=None, test=None, sampling_rate=1000.0, offset=0, minRR=None, tol=0.05)

Compare the segmentation performance of a list of R-peak positions against a reference list.

Parameters:
  • reference (array) – Reference R-peak location indices.
  • test (array) – Test R-peak location indices.
  • sampling_rate (int, float, optional) – Sampling frequency (Hz).
  • offset (int, optional) – Constant a priori offset (number of samples) between reference and test R-peak locations.
  • minRR (float, optional) – Minimum admissible RR interval (seconds).
  • tol (float, optional) – Tolerance between corresponding reference and test R-peak locations (seconds).
Returns:

  • TP (int) – Number of true positive R-peaks.
  • FP (int) – Number of false positive R-peaks.
  • performance (float) – Test performance; TP / len(reference).
  • acc (float) – Accuracy rate; TP / (TP + FP).
  • err (float) – Error rate; FP / (TP + FP).
  • match (list) – Indices of the elements of ‘test’ that match to an R-peak from ‘reference’.
  • deviation (array) – Absolute errors of the matched R-peaks (seconds).
  • mean_deviation (float) – Mean error (seconds).
  • std_deviation (float) – Standard deviation of error (seconds).
  • mean_ref_ibi (float) – Mean of the reference interbeat intervals (seconds).
  • std_ref_ibi (float) – Standard deviation of the reference interbeat intervals (seconds).
  • mean_test_ibi (float) – Mean of the test interbeat intervals (seconds).
  • std_test_ibi (float) – Standard deviation of the test interbeat intervals (seconds).

biosppy.signals.ecg.ecg(signal=None, sampling_rate=1000.0, show=True)

Process a raw ECG signal and extract relevant signal features using default parameters.

Parameters:
  • signal (array) – Raw ECG 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 ECG signal.
  • rpeaks (array) – R-peak location indices.
  • 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.ecg.engzee_segmenter(signal=None, sampling_rate=1000.0, threshold=0.48)

ECG R-peak segmentation algorithm.

Follows the approach by Engelse and Zeelenberg [EnZe79] with the modifications by Lourenco et al. [LSLL12].

Parameters:
  • signal (array) – Input filtered ECG signal.
  • sampling_rate (int, float, optional) – Sampling frequency (Hz).
  • threshold (float, optional) – Detection threshold.
Returns:

rpeaks (array) – R-peak location indices.

References

[EnZe79]W. Engelse and C. Zeelenberg, “A single scan algorithm for QRS detection and feature extraction”, IEEE Comp. in Cardiology, vol. 6, pp. 37-42, 1979
[LSLL12]A. Lourenco, H. Silva, P. Leite, R. Lourenco and A. Fred, “Real Time Electrocardiogram Segmentation for Finger Based ECG Biometrics”, BIOSIGNALS 2012, pp. 49-54, 2012
biosppy.signals.ecg.extract_heartbeats(signal=None, rpeaks=None, sampling_rate=1000.0, before=0.2, after=0.4)

Extract heartbeat templates from an ECG signal, given a list of R-peak locations.

Parameters:
  • signal (array) – Input ECG signal.
  • rpeaks (array) – R-peak location indices.
  • sampling_rate (int, float, optional) – Sampling frequency (Hz).
  • before (float, optional) – Window size to include before the R peak (seconds).
  • after (int, optional) – Window size to include after the R peak (seconds).
Returns:

  • templates (array) – Extracted heartbeat templates.
  • rpeaks (array) – Corresponding R-peak location indices of the extracted heartbeat templates.

biosppy.signals.ecg.gamboa_segmenter(signal=None, sampling_rate=1000.0, tol=0.002)

ECG R-peak segmentation algorithm.

Follows the approach by Gamboa.

Parameters:
  • signal (array) – Input filtered ECG signal.
  • sampling_rate (int, float, optional) – Sampling frequency (Hz).
  • tol (float, optional) – Tolerance parameter.
Returns:

rpeaks (array) – R-peak location indices.

biosppy.signals.ecg.hamilton_segmenter(signal=None, sampling_rate=1000.0)

ECG R-peak segmentation algorithm.

Follows the approach by Hamilton [Hami02].

Parameters:
  • signal (array) – Input filtered ECG signal.
  • sampling_rate (int, float, optional) – Sampling frequency (Hz).
Returns:

rpeaks (array) – R-peak location indices.

References

[Hami02]P.S. Hamilton, “Open Source ECG Analysis Software Documentation”, E.P.Limited, 2002
biosppy.signals.ecg.ssf_segmenter(signal=None, sampling_rate=1000.0, threshold=20, before=0.03, after=0.01)

ECG R-peak segmentation based on the Slope Sum Function (SSF).

Parameters:
  • signal (array) – Input filtered ECG signal.
  • sampling_rate (int, float, optional) – Sampling frequency (Hz).
  • threshold (float, optional) – SSF threshold.
  • before (float, optional) – Search window size before R-peak candidate (seconds).
  • after (float, optional) – Search window size after R-peak candidate (seconds).
Returns:

rpeaks (array) – R-peak location indices.

biosppy.signals.eda

This module provides methods to process Electrodermal Activity (EDA) signals, also known as Galvanic Skin Response (GSR).

copyright:
  1. 2015 by Instituto de Telecomunicacoes
license:

BSD 3-clause, see LICENSE for more details.

biosppy.signals.eda.basic_scr(signal=None, sampling_rate=1000.0)

Basic method to extract Skin Conductivity Responses (SCR) from an EDA signal.

Parameters:
  • signal (array) – Input filterd EDA signal.
  • sampling_rate (int, float, optional) – Sampling frequency (Hz).
Returns:

  • onsets (array) – Indices of the SCR onsets.
  • peaks (array) – Indices of the SRC peaks.
  • amplitudes (array) – SCR pulse amplitudes.

biosppy.signals.eda.eda(signal=None, sampling_rate=1000.0, show=True)

Process a raw EDA signal and extract relevant signal features using default parameters.

Parameters:
  • signal (array) – Raw EDA 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 EDA signal.
  • onsets (array) – Indices of SCR pulse onsets.
  • peaks (array) – Indices of the SCR peaks.
  • amplitudes (array) – SCR pulse amplitudes.

biosppy.signals.eda.kbk_scr(signal=None, sampling_rate=1000.0)

KBK method to extract Skin Conductivity Responses (SCR) from an EDA signal.

Follows the approach by Kim et al. [KiBK04].

Parameters:
  • signal (array) – Input filterd EDA signal.
  • sampling_rate (int, float, optional) – Sampling frequency (Hz).
Returns:

  • onsets (array) – Indices of the SCR onsets.
  • peaks (array) – Indices of the SRC peaks.
  • amplitudes (array) – SCR pulse amplitudes.

References

[KiBK04]K.H. Kim, S.W. Bang, and S.R. Kim, “Emotion recognition system using short-term monitoring of physiological signals”, Med. Biol. Eng. Comput., vol. 42, pp. 419-427, 2004

biosppy.signals.eeg

This module provides methods to process Electroencephalographic (EEG) signals.

copyright:
  1. 2015 by Instituto de Telecomunicacoes
license:

BSD 3-clause, see LICENSE for more details.

biosppy.signals.eeg.car_reference(signal=None)

Change signal reference to the Common Average Reference (CAR).

Parameters:signal (array) – Input EEG signal matrix; each column is one EEG channel.
Returns:signal (array) – Re-referenced EEG signal matrix; each column is one EEG channel.
biosppy.signals.eeg.eeg(signal=None, sampling_rate=1000.0, labels=None, show=True)

Process raw EEG signals and extract relevant signal features using default parameters.

Parameters:
  • signal (array) – Raw EEG signal matrix; each column is one EEG channel.
  • sampling_rate (int, float, optional) – Sampling frequency (Hz).
  • labels (list, optional) – Channel labels.
  • show (bool, optional) – If True, show a summary plot.
Returns:

  • ts (array) – Signal time axis reference (seconds).
  • filtered (array) – Filtered BVP signal.
  • features_ts (array) – Features time axis reference (seconds).
  • theta (array) – Average power in the 4 to 8 Hz frequency band; each column is one EEG channel.
  • alpha_low (array) – Average power in the 8 to 10 Hz frequency band; each column is one EEG channel.
  • alpha_high (array) – Average power in the 10 to 13 Hz frequency band; each column is one EEG channel.
  • beta (array) – Average power in the 13 to 25 Hz frequency band; each column is one EEG channel.
  • gamma (array) – Average power in the 25 to 40 Hz frequency band; each column is one EEG channel.
  • plf_pairs (list) – PLF pair indices.
  • plf (array) – PLF matrix; each column is a channel pair.

biosppy.signals.eeg.get_plf_features(signal=None, sampling_rate=1000.0, size=0.25, overlap=0.5)

Extract Phase-Locking Factor (PLF) features from EEG signals between all channel pairs.

Parameters:
  • signal (array) – Filtered EEG signal matrix; each column is one EEG channel.
  • sampling_rate (int, float, optional) – Sampling frequency (Hz).
  • size (float, optional) – Window size (seconds).
  • overlap (float, optional) – Window overlap (0 to 1).
Returns:

  • ts (array) – Features time axis reference (seconds).
  • plf_pairs (list) – PLF pair indices.
  • plf (array) – PLF matrix; each column is a channel pair.

biosppy.signals.eeg.get_power_features(signal=None, sampling_rate=1000.0, size=0.25, overlap=0.5)

Extract band power features from EEG signals.

Computes the average signal power, with overlapping windows, in typical EEG frequency bands: * Theta: from 4 to 8 Hz, * Lower Alpha: from 8 to 10 Hz, * Higher Alpha: from 10 to 13 Hz, * Beta: from 13 to 25 Hz, * Gamma: from 25 to 40 Hz.

Parameters:
  • array (signal) – Filtered EEG signal matrix; each column is one EEG channel.
  • sampling_rate (int, float, optional) – Sampling frequency (Hz).
  • size (float, optional) – Window size (seconds).
  • overlap (float, optional) – Window overlap (0 to 1).
Returns:

  • ts (array) – Features time axis reference (seconds).
  • theta (array) – Average power in the 4 to 8 Hz frequency band; each column is one EEG channel.
  • alpha_low (array) – Average power in the 8 to 10 Hz frequency band; each column is one EEG channel.
  • alpha_high (array) – Average power in the 10 to 13 Hz frequency band; each column is one EEG channel.
  • beta (array) – Average power in the 13 to 25 Hz frequency band; each column is one EEG channel.
  • gamma (array) – Average power in the 25 to 40 Hz frequency band; each column is one EEG channel.

biosppy.signals.emg

This module provides methods to process Electromyographic (EMG) signals.

copyright:
  1. 2015 by Instituto de Telecomunicacoes
license:

BSD 3-clause, see LICENSE for more details.

biosppy.signals.emg.emg(signal=None, sampling_rate=1000.0, show=True)

Process a raw EMG signal and extract relevant signal features using default parameters.

Parameters:
  • signal (array) – Raw EMG 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 EMG signal.
  • onsets (array) – Indices of EMG pulse onsets.

biosppy.signals.emg.find_onsets(signal=None, sampling_rate=1000.0, size=0.05, threshold=None)

Determine onsets of EMG pulses.

Skips corrupted signal parts.

Parameters:
  • signal (array) – Input filtered BVP signal.
  • sampling_rate (int, float, optional) – Sampling frequency (Hz).
  • size (float, optional) – Detection window size (seconds).
  • threshold (float, optional) – Detection threshold.
Returns:

onsets (array) – Indices of BVP pulse onsets.

biosppy.signals.resp

This module provides methods to process Respiration (Resp) signals.

copyright:
  1. 2015 by Instituto de Telecomunicacoes
license:

BSD 3-clause, see LICENSE for more details.

biosppy.signals.resp.resp(signal=None, sampling_rate=1000.0, show=True)

Process a raw Respiration signal and extract relevant signal features using default parameters.

Parameters:
  • signal (array) – Raw Respiration 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 Respiration signal.
  • zeros (array) – Indices of Respiration zero crossings.
  • resp_rate_ts (array) – Respiration rate time axis reference (seconds).
  • resp_rate (array) – Instantaneous respiration rate (Hz).

biosppy.signals.tools

This module provides various signal analysis methods in the time and frequency domains.

copyright:
  1. 2015 by Instituto de Telecomunicacoes
license:

BSD 3-clause, see LICENSE for more details.

biosppy.signals.tools.analytic_signal(signal=None, N=None)

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)

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.filter_signal(signal=None, ftype='FIR', band='lowpass', order=None, frequency=None, sampling_rate=1000.0, **kwargs)

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’).
  • 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.
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')

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 extrama 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)

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.get_filter(ftype='FIR', band='lowpass', order=None, frequency=None, sampling_rate=1000.0, **kwargs)

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’).
  • 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.
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)

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 190 bpm.
biosppy.signals.tools.normalize(signal=None)

Normalize a signal.

Parameters:signal (array) – Input signal.
Returns:signal (array) – Normalized signal.
biosppy.signals.tools.phase_locking(signal1=None, signal2=None, N=None)

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)

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).
  • 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.signal_stats(signal=None)

Compute various metrics describing the signal.

Parameters:signal (array) – Input signal.
Returns:
  • mean (float) – Mean of the signal.
  • median (float) – Median of the signal.
  • max (float) – Maximum signal amplitude.
  • var (float) – Signal variance (unbiased).
  • std_dev (float) – Standard signal deviation (unbiased).
  • abs_dev (float) – Absolute signal deviation.
  • kurtosis (float) – Signal kurtosis (unbiased).
  • skew (float) – Signal skewness (unbiased).
biosppy.signals.tools.smoother(signal=None, kernel='boxzen', size=10, mirror=True, **kwargs)

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(signal1=None, signal2=None)

Align two signals based on cross-correlation.

Parameters:
  • signal1 (array) – First input signal.
  • signal2 (array) – Second input signal.
Returns:

  • delay (int) – Delay (number of samples) of ‘signal1’ in relation to ‘signal2’; if ‘delay’ < 0 , ‘signal1’ is ahead in relation to ‘signal2’; if ‘delay’ > 0 , ‘signal1’ is delayed in relation to ‘signal2’.
  • corr (float) – Value of maximum correlation.
  • synch1 (array) – Biggest possible portion of ‘signal1’ in synchronization.
  • synch2 (array) – Biggest possible portion of ‘signal2’ in synchronization.

biosppy.signals.tools.windower(signal=None, size=None, step=None, fcn=None, fcn_kwargs=None, kernel='boxcar', kernel_kwargs=None)

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)

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.