Signal processing

Basic processing procedures for analog signals (e.g., performing a z-score of a signal, or filtering a signal).

elephant.signal_processing.butter(signal, highpass_freq=None, lowpass_freq=None, order=4, filter_function='filtfilt', fs=1.0, axis=-1)[source]

Butterworth filtering function for neo.AnalogSignal.

Filter type is determined according to how values of highpass_freq and lowpass_freq are given (see “Parameters” section for details).

Parameters:
signalneo.AnalogSignal or pq.Quantity or np.ndarray

Time series data to be filtered. If pq.Quantity or np.ndarray, the sampling frequency should be given through the keyword argument fs.

highpass_freqpq.Quantity of float, optional

High-pass cut-off frequency. If float, the given value is taken as frequency in Hz. Default: None.

lowpass_freqpq.Quantity or float, optional

Low-pass cut-off frequency. If float, the given value is taken as frequency in Hz. Filter type is determined depending on the values of lowpass_freq and highpass_freq:

  • highpass_freq only (lowpass_freq is None): highpass filter
  • lowpass_freq only (highpass_freq is None): lowpass filter
  • highpass_freq < lowpass_freq: bandpass filter
  • highpass_freq > lowpass_freq: bandstop filter

Default: None.

orderint, optional

Order of the Butterworth filter. Default: 4.

filter_function{‘filtfilt’, ‘lfilter’, ‘sosfiltfilt’}, optional

Filtering function to be used. Available filters:

  • ‘filtfilt’: scipy.signal.filtfilt;
  • ‘lfilter’: scipy.signal.lfilter;
  • ‘sosfiltfilt’: scipy.signal.sosfiltfilt.

In most applications ‘filtfilt’ should be used, because it doesn’t bring about phase shift due to filtering. For numerically stable filtering, in particular higher order filters, use ‘sosfiltfilt’ (see [1]). Default: ‘filtfilt’.

fspq.Quantity or float, optional

The sampling frequency of the input time series. When given as float, its value is taken as frequency in Hz. When signal is given as neo.AnalogSignal, its attribute is used to specify the sampling frequency and this parameter is ignored. Default: 1.0.

axisint, optional

Axis along which filter is applied. Default: last axis (-1).

Returns:
filtered_signalneo.AnalogSignal or pq.Quantity or np.ndarray

Filtered input data. The shape and type is identical to those of the input signal.

Raises:
ValueError

If filter_function is not one of ‘lfilter’, ‘filtfilt’, or ‘sosfiltfilt’.

If both highpass_freq and lowpass_freq are None.

References

[1]https://github.com/NeuralEnsemble/elephant/issues/220
elephant.signal_processing.cross_correlation_function(signal, ch_pairs, env=False, nlags=None, scaleopt='unbiased')[source]

Computes unbiased estimator of the cross-correlation function.

The calculations are based on [1]:

R(\tau) = \frac{1}{N-|k|} R'(\tau) \\

where R'(\tau) = \left<x(t)y(t+\tau)\right> in a pairwise manner, i.e.:

signal[ch_pairs[0,0]] vs signal[ch_pairs[0,1]],

signal[ch_pairs[1,0]] vs signal[ch_pairs[1,1]],

and so on.

The input time series are z-scored beforehand. scaleopt controls the choice of R_{xy}(\tau) normalizer. Alternatively, returns the Hilbert envelope of R_{xy}(\tau), which is useful to determine the correlation length of oscillatory signals.

Parameters:
signal(nt, nch) neo.AnalogSignal

Signal with nt number of samples that contains nch LFP channels.

ch_pairslist or (n, 2) np.ndarray

List with n channel pairs for which to compute cross-correlation. Each element of the list must contain 2 channel indices. If np.ndarray, the second axis must have dimension 2.

envbool, optional

If True, returns the Hilbert envelope of cross-correlation function. Default: False.

nlagsint, optional

Defines the number of lags for cross-correlation function. If a float is passed, it will be rounded to the nearest integer. Number of samples of output is 2*nlags+1. If None, the number of samples of the output is equal to the number of samples of the input signal (namely nt). Default: None.

scaleopt{‘none’, ‘biased’, ‘unbiased’, ‘normalized’, ‘coeff’}, optional

Normalization option, equivalent to matlab xcorr(…, scaleopt). Specified as one of the following.

  • ‘none’: raw, unscaled cross-correlation

R_{xy}(\tau)

  • ‘biased’: biased estimate of the cross-correlation:

R_{xy,biased}(\tau) = \frac{1}{N} R_{xy}(\tau)

  • ‘unbiased’: unbiased estimate of the cross-correlation:

R_{xy,unbiased}(\tau) = \frac{1}{N-\tau} R_{xy}(\tau)

  • ‘normalized’ or ‘coeff’: normalizes the sequence so that the autocorrelations at zero lag equal 1:

R_{xy,coeff}(\tau) = \frac{1}{\sqrt{R_{xx}(0) R_{yy}(0)}}
                     R_{xy}(\tau)

Default: ‘unbiased’.

Returns:
cross_corrneo.AnalogSignal

Shape: [2*nlags+1, n] Pairwise cross-correlation functions for channel pairs given by ch_pairs. If env is True, the output is the Hilbert envelope of the pairwise cross-correlation function. This is helpful to compute the correlation length for oscillating cross-correlation functions.

Raises:
ValueError

If input signal is not a neo.AnalogSignal.

If ch_pairs is not a list of channel pair indices with shape (n,2).

If env is not a boolean.

If nlags is not a positive integer.

If scaleopt is not one of the predefined above keywords.

References

[1]Stoica, P., & Moses, R. (2005). Spectral Analysis of Signals. Prentice Hall. Retrieved from http://user.it.uu.se/~ps/SAS-new.pdf, Eq. 2.2.3.

Examples

>>> import neo
>>> import quantities as pq
>>> import matplotlib.pyplot as plt
...
>>> dt = 0.02
>>> N = 2018
>>> f = 0.5
>>> t = np.arange(N)*dt
>>> x = np.zeros((N,2))
>>> x[:,0] = 0.2 * np.sin(2.*np.pi*f*t)
>>> x[:,1] = 5.3 * np.cos(2.*np.pi*f*t)
...
>>> # Generate neo.AnalogSignals from x and find cross-correlation
>>> signal = neo.AnalogSignal(x, units='mV', t_start=0.*pq.ms,
>>>     sampling_rate=1/dt*pq.Hz, dtype=float)
>>> rho = cross_correlation_function(signal, [0,1], nlags=150)
>>> env = cross_correlation_function(signal, [0,1], nlags=150,
...     env=True)
...
>>> plt.plot(rho.times, rho)
>>> plt.plot(env.times, env) # should be equal to one
>>> plt.show()
elephant.signal_processing.derivative(signal)[source]

Calculate the derivative of a neo.AnalogSignal.

Parameters:
signalneo.AnalogSignal

The signal to differentiate. If signal contains more than one channel, each is differentiated separately.

Returns:
derivative_sig: neo.AnalogSignal

The returned object is a neo.AnalogSignal containing the differences between each successive sample value of the input signal divided by the sampling period. Times are centered between the successive samples of the input. The output signal will have the same number of channels as the input signal.

Raises:
TypeError

If signal is not a neo.AnalogSignal.

elephant.signal_processing.hilbert(signal, N='nextpow')[source]

Apply a Hilbert transform to a neo.AnalogSignal object in order to obtain its (complex) analytic signal.

The time series of the instantaneous angle and amplitude can be obtained as the angle (np.angle function) and absolute value (np.abs function) of the complex analytic signal, respectively.

By default, the function will zero-pad the signal to a length corresponding to the next higher power of 2. This will provide higher computational efficiency at the expense of memory. In addition, this circumvents a situation where, for some specific choices of the length of the input, scipy.signal.hilbert function will not terminate.

Parameters:
signalneo.AnalogSignal

Signal(s) to transform.

Nint or {‘none’, ‘nextpow’}, optional

Defines whether the signal is zero-padded. If ‘none’, no padding. If ‘nextpow’, zero-pad to the next length that is a power of 2. If it is an int, directly specify the length to zero-pad to (indicates the number of Fourier components). Default: ‘nextpow’.

Returns:
neo.AnalogSignal

Contains the complex analytic signal(s) corresponding to the input signal. The unit of the returned neo.AnalogSignal is dimensionless.

Raises:
ValueError:

If N is not an integer or neither ‘nextpow’ nor ‘none’.

Notes

If N is an integer, this is passed as the parameter N of scipy.signal.hilbert function.

Examples

Create a sine signal at 5 Hz with increasing amplitude and calculate the instantaneous phases:

>>> import numpy as np
>>> import quantities as pq
>>> import neo
>>> import matplotlib.pyplot as plt
...
>>> t = np.arange(0, 5000) * pq.ms
>>> f = 5. * pq.Hz
>>> a = neo.AnalogSignal(
...       np.array(
...           (1 + t.magnitude / t[-1].magnitude) * np.sin(
...               2. * np.pi * f * t.rescale(pq.s))).reshape(
...                   (-1,1)) * pq.mV,
...       t_start=0*pq.s,
...       sampling_rate=1000*pq.Hz)
...
>>> analytic_signal = hilbert(a, N='nextpow')
>>> angles = np.angle(analytic_signal)
>>> amplitudes = np.abs(analytic_signal)
>>> print(angles)
[[-1.57079633]
 [-1.51334228]
 [-1.46047675]
 ...,
 [-1.73112977]
 [-1.68211683]
 [-1.62879501]]
>>> plt.plot(t, angles)
elephant.signal_processing.rauc(signal, baseline=None, bin_duration=None, t_start=None, t_stop=None)[source]

Calculate the rectified area under the curve (RAUC) for a neo.AnalogSignal.

The signal is optionally divided into bins with duration bin_duration, and the rectified signal (absolute value) is integrated within each bin to find the area under the curve. The mean or median of the signal or an arbitrary baseline may optionally be subtracted before rectification.

Parameters:
signalneo.AnalogSignal

The signal to integrate. If signal contains more than one channel, each is integrated separately.

baselinepq.Quantity or {‘mean’, ‘median’}, optional

A factor to subtract from the signal before rectification. If ‘mean’, the mean value of the entire signal is subtracted on a channel-by-channel basis. If ‘median’, the median value of the entire signal is subtracted on a channel-by-channel basis. Default: None.

bin_durationpq.Quantity, optional

The length of time that each integration should span. If None, there will be only one bin spanning the entire signal duration. If bin_duration does not divide evenly into the signal duration, the end of the signal is padded with zeros to accomodate the final, overextending bin. Default: None.

t_start: pq.Quantity, optional

Time to start the algorithm. If None, starts at the beginning of signal. Default: None.

t_stoppq.Quantity, optional

Time to end the algorithm. If None, ends at the last time of signal. The signal is cropped using signal.time_slice(t_start, t_stop) after baseline removal. Useful if you want the RAUC for a short section of the signal but want the mean or median calculation (`baseline`=’mean’ or `baseline`=’median’) to use the entire signal for better baseline estimation. Default: None.

Returns:
pq.Quantity or neo.AnalogSignal

If the number of bins is 1, the returned object is a scalar or vector pq.Quantity containing a single RAUC value for each channel. Otherwise, the returned object is a neo.AnalogSignal containing the RAUC(s) for each bin stored as a sample, with times corresponding to the center of each bin. The output signal will have the same number of channels as the input signal.

Raises:
ValueError

If signal is not neo.AnalogSignal.

If bin_duration is not None or pq.Quantity.

If baseline is not None, ‘mean’, ‘median’, or pq.Quantity.

See also

neo.AnalogSignal.time_slice
how t_start and t_stop are used
elephant.signal_processing.wavelet_transform(signal, freq, nco=6.0, fs=1.0, zero_padding=True)[source]

Compute the wavelet transform of a given signal with Morlet mother wavelet.

The parametrization of the wavelet is based on [1].

Parameters:
signal(Nt, Nch) neo.AnalogSignal or np.ndarray or list

Time series data to be wavelet-transformed. When multi-dimensional np.ndarray or list is given, the time axis must be the last dimension. If neo.AnalogSignal, Nt is the number of time points and Nch is the number of channels.

freqfloat or list of float

Center frequency of the Morlet wavelet in Hz. Multiple center frequencies can be given as a list, in which case the function computes the wavelet transforms for all the given frequencies at once.

ncofloat, optional

Size of the mother wavelet (approximate number of oscillation cycles within a wavelet). A larger nco value leads to a higher frequency resolution and a lower temporal resolution, and vice versa. Typically used values are in a range of 3–8, but one should be cautious when using a value smaller than ~ 6, in which case the admissibility of the wavelet is not ensured (cf. [2]). Default: 6.0.

fsfloat, optional

Sampling rate of the input data in Hz. When signal is given as a neo.AnalogSignal, the sampling frequency is taken from its attribute and this parameter is ignored. Default: 1.0.

zero_paddingbool, optional

Specifies whether the data length is extended to the least power of 2 greater than the original length, by padding zeros to the tail, for speeding up the computation. If True, the extended part is cut out from the final result before returned, so that the output has the same length as the input. Default: True.

Returns:
signal_wtnp.ndarray

Wavelet transform of the input data. When freq was given as a list, the way how the wavelet transforms for different frequencies are returned depends on the input type:

  • when the input was a neo.AnalogSignal, the returned array has shape (Nt, Nch, Nf), where Nf = len(freq), such that the last dimension indexes the frequencies;
  • when the input was a np.ndarray or list of shape (a, b, …, c, Nt), the returned array has a shape (a, b, …, c, Nf, Nt), such that the second last dimension indexes the frequencies.

To summarize, signal_wt.ndim = signal.ndim + 1, with the additional dimension in the last axis (for neo.AnalogSignal input) or the second last axis (np.ndarray or list input) indexing the frequencies.

Raises:
ValueError

If freq (or one of the values in freq when it is a list) is greater than the half of fs.

If nco is not positive.

Notes

nco is related to the wavelet number w as w \sim 2 \pi \frac{'nco'}{6}, as defined in [1].

References

[1](1, 2) M. Le Van Quyen, J. Foucher, J. Lachaux, E. Rodriguez, A. Lutz, J. Martinerie, & F.J. Varela, “Comparison of Hilbert transform and wavelet methods for the analysis of neuronal synchrony,” J Neurosci Meth, vol. 111, pp. 83–98, 2001.
[2]M. Farge, “Wavelet Transforms and their Applications to Turbulence,” Annu Rev Fluid Mech, vol. 24, pp. 395–458, 1992.
elephant.signal_processing.zscore(signal, inplace=True)[source]

Apply a z-score operation to one or several neo.AnalogSignal objects.

The z-score operation subtracts the mean \mu of the signal, and divides by its standard deviation \sigma:

Z(x(t)) = \frac{x(t)-\mu}{\sigma}

If a neo.AnalogSignal object containing multiple signals is provided, the z-transform is always calculated for each signal individually.

If a list of neo.AnalogSignal objects is supplied, the mean and standard deviation are calculated across all objects of the list. Thus, all list elements are z-transformed by the same values of \\mu and \sigma. For a neo.AnalogSignal that contains multiple signals, each signal of the array is treated separately across list elements. Therefore, the number of signals must be identical for each neo.AnalogSignal object of the list.

Parameters:
signalneo.AnalogSignal or list of neo.AnalogSignal

Signals for which to calculate the z-score.

inplacebool, optional

If True, the contents of the input signal is replaced by the z-transformed signal. If False, a copy of the original signal is returned. Default: True.

Returns:
signal_ztransofrmedneo.AnalogSignal or list of neo.AnalogSignal

The output format matches the input format: for each input neo.AnalogSignal, a corresponding neo.AnalogSignal is returned, containing the z-transformed signal with dimensionless unit.

Notes

You may supply a list of neo.AnalogSignal objects, where each object in the list contains the data of one trial of the experiment, and each signal of the neo.AnalogSignal corresponds to the recordings from one specific electrode in a particular trial. In this scenario, you will z-transform the signal of each electrode separately, but transform all trials of a given electrode in the same way.

Examples

Z-transform a single neo.AnalogSignal, containing only a single signal.

>>> import neo
>>> import numpy as np
>>> import quantities as pq
...
>>> a = neo.AnalogSignal(
...       np.array([1, 2, 3, 4, 5, 6]).reshape(-1,1) * pq.mV,
...       t_start=0*pq.s, sampling_rate=1000*pq.Hz)
>>> zscore(a).as_quantity()
[[-1.46385011]
 [-0.87831007]
 [-0.29277002]
 [ 0.29277002]
 [ 0.87831007]
 [ 1.46385011]] dimensionless

Z-transform a single neo.AnalogSignal containing multiple signals.

>>> b = neo.AnalogSignal(
...       np.transpose([[1, 2, 3, 4, 5, 6],
...                     [11, 12, 13, 14, 15, 16]]) * pq.mV,
...       t_start=0*pq.s, sampling_rate=1000*pq.Hz)
>>> zscore(b).as_quantity()
[[-1.46385011 -1.46385011]
 [-0.87831007 -0.87831007]
 [-0.29277002 -0.29277002]
 [ 0.29277002  0.29277002]
 [ 0.87831007  0.87831007]
 [ 1.46385011  1.46385011]] dimensionless

Z-transform a list of neo.AnalogSignal, each one containing more than one signal:

>>> c = neo.AnalogSignal(
...       np.transpose([[21, 22, 23, 24, 25, 26],
...                     [31, 32, 33, 34, 35, 36]]) * pq.mV,
...       t_start=0*pq.s, sampling_rate=1000*pq.Hz)
>>> zscore([b, c])
[<AnalogSignal(array([[-1.11669108, -1.08361877],
   [-1.0672076 , -1.04878252],
   [-1.01772411, -1.01394628],
   [-0.96824063, -0.97911003],
   [-0.91875714, -0.94427378],
   [-0.86927366, -0.90943753]]) * dimensionless, [0.0 s, 0.006 s],
   sampling rate: 1000.0 Hz)>,
   <AnalogSignal(array([[ 0.78170952,  0.84779261],
   [ 0.86621866,  0.90728682],
   [ 0.9507278 ,  0.96678104],
   [ 1.03523694,  1.02627526],
   [ 1.11974608,  1.08576948],
   [ 1.20425521,  1.1452637 ]]) * dimensionless, [0.0 s, 0.006 s],
   sampling rate: 1000.0 Hz)>]