signal_processing.py 36 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893894895896897898899900901902903904905906907908909910911912913914915916917918919920921922923924925926927928929930931932933934935936937938939940941942943944945946947948949950951952953954955956957
  1. # -*- coding: utf-8 -*-
  2. """
  3. Basic processing procedures for analog signals (e.g., performing a z-score of a
  4. signal, or filtering a signal).
  5. :copyright: Copyright 2014-2020 by the Elephant team, see `doc/authors.rst`.
  6. :license: Modified BSD, see LICENSE.txt for details.
  7. """
  8. from __future__ import division, print_function, unicode_literals
  9. import neo
  10. import numpy as np
  11. import quantities as pq
  12. import scipy.signal
  13. from elephant.utils import deprecated_alias, check_same_units
  14. __all__ = [
  15. "zscore",
  16. "cross_correlation_function",
  17. "butter",
  18. "wavelet_transform",
  19. "hilbert",
  20. "rauc",
  21. "derivative"
  22. ]
  23. def zscore(signal, inplace=True):
  24. r"""
  25. Apply a z-score operation to one or several `neo.AnalogSignal` objects.
  26. The z-score operation subtracts the mean :math:`\mu` of the signal, and
  27. divides by its standard deviation :math:`\sigma`:
  28. .. math::
  29. Z(x(t)) = \frac{x(t)-\mu}{\sigma}
  30. If a `neo.AnalogSignal` object containing multiple signals is provided,
  31. the z-transform is always calculated for each signal individually.
  32. If a list of `neo.AnalogSignal` objects is supplied, the mean and standard
  33. deviation are calculated across all objects of the list. Thus, all list
  34. elements are z-transformed by the same values of :math:`\\mu` and
  35. :math:`\sigma`. For a `neo.AnalogSignal` that contains multiple signals,
  36. each signal of the array is treated separately across list elements.
  37. Therefore, the number of signals must be identical for each
  38. `neo.AnalogSignal` object of the list.
  39. Parameters
  40. ----------
  41. signal : neo.AnalogSignal or list of neo.AnalogSignal
  42. Signals for which to calculate the z-score.
  43. inplace : bool, optional
  44. If True, the contents of the input `signal` is replaced by the
  45. z-transformed signal, if possible, i.e when the signal type is float.
  46. If False, a copy of the original `signal` is returned.
  47. Default: True
  48. Returns
  49. -------
  50. signal_ztransofrmed : neo.AnalogSignal or list of neo.AnalogSignal
  51. The output format matches the input format: for each input
  52. `neo.AnalogSignal`, a corresponding `neo.AnalogSignal` is returned,
  53. containing the z-transformed signal with dimensionless unit.
  54. Notes
  55. -----
  56. You may supply a list of `neo.AnalogSignal` objects, where each object in
  57. the list contains the data of one trial of the experiment, and each signal
  58. of the `neo.AnalogSignal` corresponds to the recordings from one specific
  59. electrode in a particular trial. In this scenario, you will z-transform
  60. the signal of each electrode separately, but transform all trials of a
  61. given electrode in the same way.
  62. Examples
  63. --------
  64. Z-transform a single `neo.AnalogSignal`, containing only a single signal.
  65. >>> import neo
  66. >>> import numpy as np
  67. >>> import quantities as pq
  68. ...
  69. >>> a = neo.AnalogSignal(
  70. ... np.array([1, 2, 3, 4, 5, 6]).reshape(-1,1) * pq.mV,
  71. ... t_start=0*pq.s, sampling_rate=1000*pq.Hz)
  72. >>> zscore(a).as_quantity()
  73. [[-1.46385011]
  74. [-0.87831007]
  75. [-0.29277002]
  76. [ 0.29277002]
  77. [ 0.87831007]
  78. [ 1.46385011]] dimensionless
  79. Z-transform a single `neo.AnalogSignal` containing multiple signals.
  80. >>> b = neo.AnalogSignal(
  81. ... np.transpose([[1, 2, 3, 4, 5, 6],
  82. ... [11, 12, 13, 14, 15, 16]]) * pq.mV,
  83. ... t_start=0*pq.s, sampling_rate=1000*pq.Hz)
  84. >>> zscore(b).as_quantity()
  85. [[-1.46385011 -1.46385011]
  86. [-0.87831007 -0.87831007]
  87. [-0.29277002 -0.29277002]
  88. [ 0.29277002 0.29277002]
  89. [ 0.87831007 0.87831007]
  90. [ 1.46385011 1.46385011]] dimensionless
  91. Z-transform a list of `neo.AnalogSignal`, each one containing more than
  92. one signal:
  93. >>> c = neo.AnalogSignal(
  94. ... np.transpose([[21, 22, 23, 24, 25, 26],
  95. ... [31, 32, 33, 34, 35, 36]]) * pq.mV,
  96. ... t_start=0*pq.s, sampling_rate=1000*pq.Hz)
  97. >>> zscore([b, c])
  98. [<AnalogSignal(array([[-1.11669108, -1.08361877],
  99. [-1.0672076 , -1.04878252],
  100. [-1.01772411, -1.01394628],
  101. [-0.96824063, -0.97911003],
  102. [-0.91875714, -0.94427378],
  103. [-0.86927366, -0.90943753]]) * dimensionless, [0.0 s, 0.006 s],
  104. sampling rate: 1000.0 Hz)>,
  105. <AnalogSignal(array([[ 0.78170952, 0.84779261],
  106. [ 0.86621866, 0.90728682],
  107. [ 0.9507278 , 0.96678104],
  108. [ 1.03523694, 1.02627526],
  109. [ 1.11974608, 1.08576948],
  110. [ 1.20425521, 1.1452637 ]]) * dimensionless, [0.0 s, 0.006 s],
  111. sampling rate: 1000.0 Hz)>]
  112. """
  113. # Transform input to a list
  114. if isinstance(signal, neo.AnalogSignal):
  115. signal = [signal]
  116. check_same_units(signal, object_type=neo.AnalogSignal)
  117. # Calculate mean and standard deviation
  118. signal_stacked = np.vstack(signal).magnitude
  119. mean = signal_stacked.mean(axis=0)
  120. std = signal_stacked.std(axis=0)
  121. signal_ztransofrmed = []
  122. for sig in signal:
  123. sig_normalized = sig.magnitude.astype(mean.dtype, copy=not inplace)
  124. sig_normalized -= mean
  125. # items where std is zero are already zero
  126. np.divide(sig_normalized, std, out=sig_normalized, where=std != 0)
  127. sig_dimless = neo.AnalogSignal(signal=sig_normalized,
  128. units=pq.dimensionless,
  129. dtype=sig_normalized.dtype,
  130. copy=False,
  131. t_start=sig.t_start,
  132. sampling_rate=sig.sampling_rate,
  133. name=sig.name,
  134. file_origin=sig.file_origin,
  135. description=sig.description,
  136. array_annotations=sig.array_annotations,
  137. **sig.annotations)
  138. signal_ztransofrmed.append(sig_dimless)
  139. # Return single object, or list of objects
  140. if len(signal_ztransofrmed) == 1:
  141. signal_ztransofrmed = signal_ztransofrmed[0]
  142. return signal_ztransofrmed
  143. @deprecated_alias(ch_pairs='channel_pairs', nlags='n_lags',
  144. env='hilbert_envelope')
  145. def cross_correlation_function(signal, channel_pairs, hilbert_envelope=False,
  146. n_lags=None, scaleopt='unbiased'):
  147. r"""
  148. Computes unbiased estimator of the cross-correlation function.
  149. The calculations are based on [1]_:
  150. .. math::
  151. R(\tau) = \frac{1}{N-|k|} R'(\tau) \\
  152. where :math:`R'(\tau) = \left<x(t)y(t+\tau)\right>` in a pairwise
  153. manner, i.e.:
  154. `signal[channel_pairs[0,0]]` vs `signal[channel_pairs[0,1]]`,
  155. `signal[channel_pairs[1,0]]` vs `signal[channel_pairs[1,1]]`,
  156. and so on.
  157. The input time series are z-scored beforehand. `scaleopt` controls the
  158. choice of :math:`R_{xy}(\tau)` normalizer. Alternatively, returns the
  159. Hilbert envelope of :math:`R_{xy}(\tau)`, which is useful to determine the
  160. correlation length of oscillatory signals.
  161. Parameters
  162. ----------
  163. signal : (nt, nch) neo.AnalogSignal
  164. Signal with `nt` number of samples that contains `nch` LFP channels.
  165. channel_pairs : list or (n, 2) np.ndarray
  166. List with `n` channel pairs for which to compute cross-correlation.
  167. Each element of the list must contain 2 channel indices.
  168. If `np.ndarray`, the second axis must have dimension 2.
  169. hilbert_envelope : bool, optional
  170. If True, returns the Hilbert envelope of cross-correlation function
  171. result.
  172. Default: False.
  173. n_lags : int, optional
  174. Defines the number of lags for cross-correlation function. If a `float`
  175. is passed, it will be rounded to the nearest integer. Number of
  176. samples of output is `2*n_lags+1`.
  177. If None, the number of samples of the output is equal to the number of
  178. samples of the input signal (namely `nt`).
  179. Default: None.
  180. scaleopt : {'none', 'biased', 'unbiased', 'normalized', 'coeff'}, optional
  181. Normalization option, equivalent to matlab `xcorr(..., scaleopt)`.
  182. Specified as one of the following.
  183. * 'none': raw, unscaled cross-correlation
  184. .. math::
  185. R_{xy}(\tau)
  186. * 'biased': biased estimate of the cross-correlation:
  187. .. math::
  188. R_{xy,biased}(\tau) = \frac{1}{N} R_{xy}(\tau)
  189. * 'unbiased': unbiased estimate of the cross-correlation:
  190. .. math::
  191. R_{xy,unbiased}(\tau) = \frac{1}{N-\tau} R_{xy}(\tau)
  192. * 'normalized' or 'coeff': normalizes the sequence so that the
  193. autocorrelations at zero lag equal 1:
  194. .. math::
  195. R_{xy,coeff}(\tau) = \frac{1}{\sqrt{R_{xx}(0) R_{yy}(0)}}
  196. R_{xy}(\tau)
  197. Default: 'unbiased'.
  198. Returns
  199. -------
  200. cross_corr : neo.AnalogSignal
  201. Shape: `[2*n_lags+1, n]`
  202. Pairwise cross-correlation functions for channel pairs given by
  203. `channel_pairs`. If `hilbert_envelope` is True, the output is the
  204. Hilbert envelope of the pairwise cross-correlation function. This is
  205. helpful to compute the correlation length for oscillating
  206. cross-correlation functions.
  207. Raises
  208. ------
  209. ValueError
  210. If input `signal` is not a `neo.AnalogSignal`.
  211. If `channel_pairs` is not a list of channel pair indices with shape
  212. `(n,2)`.
  213. If `hilbert_envelope` is not a boolean.
  214. If `n_lags` is not a positive integer.
  215. If `scaleopt` is not one of the predefined above keywords.
  216. References
  217. ----------
  218. .. [1] Stoica, P., & Moses, R. (2005). Spectral Analysis of Signals.
  219. Prentice Hall. Retrieved from http://user.it.uu.se/~ps/SAS-new.pdf,
  220. Eq. 2.2.3.
  221. Examples
  222. --------
  223. >>> import neo
  224. >>> import quantities as pq
  225. >>> import matplotlib.pyplot as plt
  226. ...
  227. >>> dt = 0.02
  228. >>> N = 2018
  229. >>> f = 0.5
  230. >>> t = np.arange(N)*dt
  231. >>> x = np.zeros((N,2))
  232. >>> x[:,0] = 0.2 * np.sin(2.*np.pi*f*t)
  233. >>> x[:,1] = 5.3 * np.cos(2.*np.pi*f*t)
  234. ...
  235. >>> # Generate neo.AnalogSignals from x and find cross-correlation
  236. >>> signal = neo.AnalogSignal(x, units='mV', t_start=0.*pq.ms,
  237. >>> sampling_rate=1/dt*pq.Hz, dtype=float)
  238. >>> rho = cross_correlation_function(signal, [0,1], n_lags=150)
  239. >>> env = cross_correlation_function(signal, [0,1], n_lags=150,
  240. ... hilbert_envelope=True)
  241. ...
  242. >>> plt.plot(rho.times, rho)
  243. >>> plt.plot(env.times, env) # should be equal to one
  244. >>> plt.show()
  245. """
  246. # Make channel_pairs a 2D array
  247. pairs = np.asarray(channel_pairs)
  248. if pairs.ndim == 1:
  249. pairs = np.expand_dims(pairs, axis=0)
  250. # Check input
  251. if not isinstance(signal, neo.AnalogSignal):
  252. raise ValueError('Input signal must be of type neo.AnalogSignal')
  253. if pairs.shape[1] != 2:
  254. raise ValueError("'channel_pairs' is not a list of channel pair "
  255. "indices. Cannot define pairs for cross-correlation.")
  256. if not isinstance(hilbert_envelope, bool):
  257. raise ValueError("'hilbert_envelope' must be a boolean value")
  258. if n_lags is not None:
  259. if not isinstance(n_lags, int) or n_lags <= 0:
  260. raise ValueError('n_lags must be a non-negative integer')
  261. # z-score analog signal and store channel time series in different arrays
  262. # Cross-correlation will be calculated between xsig and ysig
  263. z_transformed = signal.magnitude - signal.magnitude.mean(axis=0)
  264. z_transformed = np.divide(z_transformed, signal.magnitude.std(axis=0),
  265. out=z_transformed,
  266. where=z_transformed != 0)
  267. # transpose (nch, xy, nt) -> (xy, nt, nch)
  268. xsig, ysig = np.transpose(z_transformed.T[pairs], (1, 2, 0))
  269. # Define vector of lags tau
  270. nt, nch = xsig.shape
  271. tau = np.arange(nt) - nt // 2
  272. # Calculate cross-correlation by taking Fourier transform of signal,
  273. # multiply in Fourier space, and transform back. Correct for bias due
  274. # to zero-padding
  275. xcorr = scipy.signal.fftconvolve(xsig, ysig[::-1], mode='same', axes=0)
  276. if scaleopt == 'biased':
  277. xcorr /= nt
  278. elif scaleopt == 'unbiased':
  279. normalizer = np.expand_dims(nt - np.abs(tau), axis=1)
  280. xcorr /= normalizer
  281. elif scaleopt in ('normalized', 'coeff'):
  282. normalizer = np.sqrt((xsig ** 2).sum(axis=0) * (ysig ** 2).sum(axis=0))
  283. xcorr /= normalizer
  284. elif scaleopt != 'none':
  285. raise ValueError("Invalid scaleopt mode: '{}'".format(scaleopt))
  286. # Calculate envelope of cross-correlation function with Hilbert transform.
  287. # This is useful for transient oscillatory signals.
  288. if hilbert_envelope:
  289. xcorr = np.abs(scipy.signal.hilbert(xcorr, axis=0))
  290. # Cut off lags outside the desired range
  291. if n_lags is not None:
  292. tau0 = np.argwhere(tau == 0).item()
  293. xcorr = xcorr[tau0 - n_lags: tau0 + n_lags + 1, :]
  294. # Return neo.AnalogSignal
  295. cross_corr = neo.AnalogSignal(xcorr,
  296. units='',
  297. t_start=tau[0] * signal.sampling_period,
  298. t_stop=tau[-1] * signal.sampling_period,
  299. sampling_rate=signal.sampling_rate,
  300. dtype=float)
  301. return cross_corr
  302. @deprecated_alias(highpass_freq='highpass_frequency',
  303. lowpass_freq='lowpass_frequency',
  304. fs='sampling_frequency')
  305. def butter(signal, highpass_frequency=None, lowpass_frequency=None, order=4,
  306. filter_function='filtfilt', sampling_frequency=1.0, axis=-1):
  307. """
  308. Butterworth filtering function for `neo.AnalogSignal`.
  309. Filter type is determined according to how values of `highpass_frequency`
  310. and `lowpass_frequency` are given (see "Parameters" section for details).
  311. Parameters
  312. ----------
  313. signal : neo.AnalogSignal or pq.Quantity or np.ndarray
  314. Time series data to be filtered.
  315. If `pq.Quantity` or `np.ndarray`, the sampling frequency should be
  316. given through the keyword argument `fs`.
  317. highpass_frequency : pq.Quantity of float, optional
  318. High-pass cut-off frequency. If `float`, the given value is taken as
  319. frequency in Hz.
  320. Default: None.
  321. lowpass_frequency : pq.Quantity or float, optional
  322. Low-pass cut-off frequency. If `float`, the given value is taken as
  323. frequency in Hz.
  324. Filter type is determined depending on the values of
  325. `lowpass_frequency` and `highpass_frequency`:
  326. * `highpass_frequency` only (`lowpass_frequency` is None):
  327. highpass filter
  328. * `lowpass_frequency` only (`highpass_frequency` is None):
  329. lowpass filter
  330. * `highpass_frequency` < `lowpass_frequency`: bandpass filter
  331. * `highpass_frequency` > `lowpass_frequency`: bandstop filter
  332. Default: None.
  333. order : int, optional
  334. Order of the Butterworth filter.
  335. Default: 4.
  336. filter_function : {'filtfilt', 'lfilter', 'sosfiltfilt'}, optional
  337. Filtering function to be used. Available filters:
  338. * 'filtfilt': `scipy.signal.filtfilt`;
  339. * 'lfilter': `scipy.signal.lfilter`;
  340. * 'sosfiltfilt': `scipy.signal.sosfiltfilt`.
  341. In most applications 'filtfilt' should be used, because it doesn't
  342. bring about phase shift due to filtering. For numerically stable
  343. filtering, in particular higher order filters, use 'sosfiltfilt'
  344. (see [1]_).
  345. Default: 'filtfilt'.
  346. sampling_frequency : pq.Quantity or float, optional
  347. The sampling frequency of the input time series. When given as
  348. `float`, its value is taken as frequency in Hz. When `signal` is given
  349. as `neo.AnalogSignal`, its attribute is used to specify the sampling
  350. frequency and this parameter is ignored.
  351. Default: 1.0.
  352. axis : int, optional
  353. Axis along which filter is applied.
  354. Default: last axis (-1).
  355. Returns
  356. -------
  357. filtered_signal : neo.AnalogSignal or pq.Quantity or np.ndarray
  358. Filtered input data. The shape and type is identical to those of the
  359. input `signal`.
  360. Raises
  361. ------
  362. ValueError
  363. If `filter_function` is not one of 'lfilter', 'filtfilt',
  364. or 'sosfiltfilt'.
  365. If both `highpass_frequency` and `lowpass_frequency` are None.
  366. References
  367. ----------
  368. .. [1] https://github.com/NeuralEnsemble/elephant/issues/220
  369. """
  370. available_filters = 'lfilter', 'filtfilt', 'sosfiltfilt'
  371. if filter_function not in available_filters:
  372. raise ValueError("Invalid `filter_function`: {filter_function}. "
  373. "Available filters: {available_filters}".format(
  374. filter_function=filter_function,
  375. available_filters=available_filters))
  376. # design filter
  377. if hasattr(signal, 'sampling_rate'):
  378. sampling_frequency = signal.sampling_rate.rescale(pq.Hz).magnitude
  379. if isinstance(highpass_frequency, pq.quantity.Quantity):
  380. highpass_frequency = highpass_frequency.rescale(pq.Hz).magnitude
  381. if isinstance(lowpass_frequency, pq.quantity.Quantity):
  382. lowpass_frequency = lowpass_frequency.rescale(pq.Hz).magnitude
  383. Fn = sampling_frequency / 2.
  384. # filter type is determined according to the values of cut-off
  385. # frequencies
  386. if lowpass_frequency and highpass_frequency:
  387. if highpass_frequency < lowpass_frequency:
  388. Wn = (highpass_frequency / Fn, lowpass_frequency / Fn)
  389. btype = 'bandpass'
  390. else:
  391. Wn = (lowpass_frequency / Fn, highpass_frequency / Fn)
  392. btype = 'bandstop'
  393. elif lowpass_frequency:
  394. Wn = lowpass_frequency / Fn
  395. btype = 'lowpass'
  396. elif highpass_frequency:
  397. Wn = highpass_frequency / Fn
  398. btype = 'highpass'
  399. else:
  400. raise ValueError(
  401. "Either highpass_frequency or lowpass_frequency must be given"
  402. )
  403. if filter_function == 'sosfiltfilt':
  404. output = 'sos'
  405. else:
  406. output = 'ba'
  407. designed_filter = scipy.signal.butter(order, Wn, btype=btype,
  408. output=output)
  409. # When the input is AnalogSignal, the axis for time index (i.e. the
  410. # first axis) needs to be rolled to the last
  411. data = np.asarray(signal)
  412. if isinstance(signal, neo.AnalogSignal):
  413. data = np.rollaxis(data, 0, len(data.shape))
  414. # apply filter
  415. if filter_function == 'lfilter':
  416. b, a = designed_filter
  417. filtered_data = scipy.signal.lfilter(b=b, a=a, x=data, axis=axis)
  418. elif filter_function == 'filtfilt':
  419. b, a = designed_filter
  420. filtered_data = scipy.signal.filtfilt(b=b, a=a, x=data, axis=axis)
  421. else:
  422. filtered_data = scipy.signal.sosfiltfilt(sos=designed_filter,
  423. x=data, axis=axis)
  424. if isinstance(signal, neo.AnalogSignal):
  425. filtered_data = np.rollaxis(filtered_data, -1, 0)
  426. signal_out = signal.duplicate_with_new_data(filtered_data)
  427. # todo use flag once is fixed
  428. # https://github.com/NeuralEnsemble/python-neo/issues/752
  429. signal_out.array_annotate(**signal.array_annotations)
  430. return signal_out
  431. elif isinstance(signal, pq.quantity.Quantity):
  432. return filtered_data * signal.units
  433. else:
  434. return filtered_data
  435. @deprecated_alias(nco='n_cycles', freq='frequency', fs='sampling_frequency')
  436. def wavelet_transform(signal, frequency, n_cycles=6.0, sampling_frequency=1.0,
  437. zero_padding=True):
  438. r"""
  439. Compute the wavelet transform of a given signal with Morlet mother
  440. wavelet.
  441. The parametrization of the wavelet is based on [1]_.
  442. Parameters
  443. ----------
  444. signal : (Nt, Nch) neo.AnalogSignal or np.ndarray or list
  445. Time series data to be wavelet-transformed. When multi-dimensional
  446. `np.ndarray` or list is given, the time axis must be the last
  447. dimension. If `neo.AnalogSignal`, `Nt` is the number of time points
  448. and `Nch` is the number of channels.
  449. frequency : float or list of float
  450. Center frequency of the Morlet wavelet in Hz. Multiple center
  451. frequencies can be given as a list, in which case the function
  452. computes the wavelet transforms for all the given frequencies at once.
  453. n_cycles : float, optional
  454. Size of the mother wavelet (approximate number of oscillation cycles
  455. within a wavelet). Corresponds to :math:`nco` in the paper [1]_.
  456. A larger `n_cycles` value leads to a higher frequency resolution and a
  457. lower temporal resolution, and vice versa.
  458. Typically used values are in a range of 3–8, but one should be cautious
  459. when using a value smaller than ~ 6, in which case the admissibility of
  460. the wavelet is not ensured (cf. [2]_).
  461. Default: 6.0.
  462. sampling_frequency : float, optional
  463. Sampling rate of the input data in Hz.
  464. When `signal` is given as a `neo.AnalogSignal`, the sampling frequency
  465. is taken from its attribute and this parameter is ignored.
  466. Default: 1.0.
  467. zero_padding : bool, optional
  468. Specifies whether the data length is extended to the least power of
  469. 2 greater than the original length, by padding zeros to the tail, for
  470. speeding up the computation.
  471. If True, the extended part is cut out from the final result before
  472. returned, so that the output has the same length as the input.
  473. Default: True.
  474. Returns
  475. -------
  476. signal_wt : np.ndarray
  477. Wavelet transform of the input data. When `frequency` was given as a
  478. list, the way how the wavelet transforms for different frequencies are
  479. returned depends on the input type:
  480. * when the input was a `neo.AnalogSignal`, the returned array has
  481. shape (`Nt`, `Nch`, `Nf`), where `Nf` = `len(freq)`, such that the
  482. last dimension indexes the frequencies;
  483. * when the input was a `np.ndarray` or list of shape
  484. (`a`, `b`, ..., `c`, `Nt`), the returned array has a shape
  485. (`a`, `b`, ..., `c`, `Nf`, `Nt`), such that the second last
  486. dimension indexes the frequencies.
  487. To summarize, `signal_wt.ndim` = `signal.ndim` + 1, with the
  488. additional dimension in the last axis (for `neo.AnalogSignal` input)
  489. or the second last axis (`np.ndarray` or list input) indexing the
  490. frequencies.
  491. Raises
  492. ------
  493. ValueError
  494. If `frequency` (or one of the values in `frequency` when it is a list)
  495. is greater than the half of `sampling_frequency`.
  496. If `n_cycles` is not positive.
  497. Notes
  498. -----
  499. `n_cycles` is related to the wavelet number :math:`w` as
  500. :math:`w \sim 2 \pi \frac{n_{\text{cycles}}}{6}`, as defined in [1]_.
  501. References
  502. ----------
  503. .. [1] M. Le Van Quyen, J. Foucher, J. Lachaux, E. Rodriguez, A. Lutz,
  504. J. Martinerie, & F.J. Varela, "Comparison of Hilbert transform and
  505. wavelet methods for the analysis of neuronal synchrony," J Neurosci
  506. Meth, vol. 111, pp. 83–98, 2001.
  507. .. [2] M. Farge, "Wavelet Transforms and their Applications to
  508. Turbulence," Annu Rev Fluid Mech, vol. 24, pp. 395–458, 1992.
  509. """
  510. def _morlet_wavelet_ft(freq, n_cycles, fs, n):
  511. # Generate the Fourier transform of Morlet wavelet as defined
  512. # in Le van Quyen et al. J Neurosci Meth 111:83-98 (2001).
  513. sigma = n_cycles / (6. * freq)
  514. freqs = np.fft.fftfreq(n, 1.0 / fs)
  515. heaviside = np.array(freqs > 0., dtype=np.float)
  516. ft_real = np.sqrt(2 * np.pi * freq) * sigma * np.exp(
  517. -2 * (np.pi * sigma * (freqs - freq)) ** 2) * heaviside * fs
  518. ft_imag = np.zeros_like(ft_real)
  519. return ft_real + 1.0j * ft_imag
  520. data = np.asarray(signal)
  521. # When the input is AnalogSignal, the axis for time index (i.e. the
  522. # first axis) needs to be rolled to the last
  523. if isinstance(signal, neo.AnalogSignal):
  524. data = np.rollaxis(data, 0, data.ndim)
  525. # When the input is AnalogSignal, use its attribute to specify the
  526. # sampling frequency
  527. if hasattr(signal, 'sampling_rate'):
  528. sampling_frequency = signal.sampling_rate
  529. if isinstance(sampling_frequency, pq.quantity.Quantity):
  530. sampling_frequency = sampling_frequency.rescale('Hz').magnitude
  531. if isinstance(frequency, (list, tuple, np.ndarray)):
  532. freqs = np.asarray(frequency)
  533. else:
  534. freqs = np.array([frequency, ])
  535. if isinstance(freqs[0], pq.quantity.Quantity):
  536. freqs = [f.rescale('Hz').magnitude for f in freqs]
  537. # check whether the given central frequencies are less than the
  538. # Nyquist frequency of the signal
  539. if np.any(freqs >= sampling_frequency / 2):
  540. raise ValueError("'frequency' elements must be less than the half of "
  541. "the 'sampling_frequency' ({}) Hz"
  542. .format(sampling_frequency))
  543. # check if n_cycles is positive
  544. if n_cycles <= 0:
  545. raise ValueError("`n_cycles` must be positive")
  546. n_orig = data.shape[-1]
  547. if zero_padding:
  548. n = 2 ** (int(np.log2(n_orig)) + 1)
  549. else:
  550. n = n_orig
  551. # generate Morlet wavelets (in the frequency domain)
  552. wavelet_fts = np.empty([len(freqs), n], dtype=np.complex)
  553. for i, f in enumerate(freqs):
  554. wavelet_fts[i] = _morlet_wavelet_ft(f, n_cycles, sampling_frequency, n)
  555. # perform wavelet transform by convoluting the signal with the wavelets
  556. if data.ndim == 1:
  557. data = np.expand_dims(data, 0)
  558. data = np.expand_dims(data, data.ndim - 1)
  559. data = np.fft.ifft(np.fft.fft(data, n) * wavelet_fts)
  560. signal_wt = data[..., 0:n_orig]
  561. # reshape the result array according to the input
  562. if isinstance(signal, neo.AnalogSignal):
  563. signal_wt = np.rollaxis(signal_wt, -1)
  564. if not isinstance(frequency, (list, tuple, np.ndarray)):
  565. signal_wt = signal_wt[..., 0]
  566. else:
  567. if signal.ndim == 1:
  568. signal_wt = signal_wt[0]
  569. if not isinstance(frequency, (list, tuple, np.ndarray)):
  570. signal_wt = signal_wt[..., 0, :]
  571. return signal_wt
  572. @deprecated_alias(N='padding')
  573. def hilbert(signal, padding='nextpow'):
  574. """
  575. Apply a Hilbert transform to a `neo.AnalogSignal` object in order to
  576. obtain its (complex) analytic signal.
  577. The time series of the instantaneous angle and amplitude can be obtained
  578. as the angle (`np.angle` function) and absolute value (`np.abs` function)
  579. of the complex analytic signal, respectively.
  580. By default, the function will zero-pad the signal to a length
  581. corresponding to the next higher power of 2. This will provide higher
  582. computational efficiency at the expense of memory. In addition, this
  583. circumvents a situation where, for some specific choices of the length of
  584. the input, `scipy.signal.hilbert` function will not terminate.
  585. Parameters
  586. ----------
  587. signal : neo.AnalogSignal
  588. Signal(s) to transform.
  589. padding : int, {'none', 'nextpow'}, or None, optional
  590. Defines whether the signal is zero-padded.
  591. The `padding` argument corresponds to `N` in
  592. `scipy.signal.hilbert(signal, N=padding)` function.
  593. If 'none' or None, no padding.
  594. If 'nextpow', zero-pad to the next length that is a power of 2.
  595. If it is an `int`, directly specify the length to zero-pad to
  596. (indicates the number of Fourier components).
  597. Default: 'nextpow'.
  598. Returns
  599. -------
  600. neo.AnalogSignal
  601. Contains the complex analytic signal(s) corresponding to the input
  602. `signal`. The unit of the returned `neo.AnalogSignal` is
  603. dimensionless.
  604. Raises
  605. ------
  606. ValueError:
  607. If `padding` is not an integer or neither 'nextpow' nor 'none' (None).
  608. Examples
  609. --------
  610. Create a sine signal at 5 Hz with increasing amplitude and calculate the
  611. instantaneous phases:
  612. >>> import numpy as np
  613. >>> import quantities as pq
  614. >>> import neo
  615. >>> import matplotlib.pyplot as plt
  616. ...
  617. >>> t = np.arange(0, 5000) * pq.ms
  618. >>> f = 5. * pq.Hz
  619. >>> a = neo.AnalogSignal(
  620. ... np.array(
  621. ... (1 + t.magnitude / t[-1].magnitude) * np.sin(
  622. ... 2. * np.pi * f * t.rescale(pq.s))).reshape(
  623. ... (-1,1)) * pq.mV,
  624. ... t_start=0*pq.s,
  625. ... sampling_rate=1000*pq.Hz)
  626. ...
  627. >>> analytic_signal = hilbert(a, padding='nextpow')
  628. >>> angles = np.angle(analytic_signal)
  629. >>> amplitudes = np.abs(analytic_signal)
  630. >>> print(angles)
  631. [[-1.57079633]
  632. [-1.51334228]
  633. [-1.46047675]
  634. ...,
  635. [-1.73112977]
  636. [-1.68211683]
  637. [-1.62879501]]
  638. >>> plt.plot(t, angles)
  639. """
  640. # Length of input signals
  641. n_org = signal.shape[0]
  642. # Right-pad signal to desired length using the signal itself
  643. if isinstance(padding, int):
  644. # User defined padding
  645. n = padding
  646. elif padding == 'nextpow':
  647. # To speed up calculation of the Hilbert transform, make sure we change
  648. # the signal to be of a length that is a power of two. Failure to do so
  649. # results in computations of certain signal lengths to not finish (or
  650. # finish in absurd time). This might be a bug in scipy (0.16), e.g.,
  651. # the following code will not terminate for this value of k:
  652. #
  653. # import numpy
  654. # import scipy.signal
  655. # k=679346
  656. # t = np.arange(0, k) / 1000.
  657. # a = (1 + t / t[-1]) * np.sin(2 * np.pi * 5 * t)
  658. # analytic_signal = scipy.signal.hilbert(a)
  659. #
  660. # For this reason, nextpow is the default setting for now.
  661. n = 2 ** (int(np.log2(n_org - 1)) + 1)
  662. elif padding == 'none' or padding is None:
  663. # No padding
  664. n = n_org
  665. else:
  666. raise ValueError("Invalid padding '{}'.".format(padding))
  667. output = signal.duplicate_with_new_data(
  668. scipy.signal.hilbert(signal.magnitude, N=n, axis=0)[:n_org])
  669. # todo use flag once is fixed
  670. # https://github.com/NeuralEnsemble/python-neo/issues/752
  671. output.array_annotate(**signal.array_annotations)
  672. return output / output.units
  673. def rauc(signal, baseline=None, bin_duration=None, t_start=None, t_stop=None):
  674. """
  675. Calculate the rectified area under the curve (RAUC) for a
  676. `neo.AnalogSignal`.
  677. The signal is optionally divided into bins with duration `bin_duration`,
  678. and the rectified signal (absolute value) is integrated within each bin to
  679. find the area under the curve. The mean or median of the signal or an
  680. arbitrary baseline may optionally be subtracted before rectification.
  681. Parameters
  682. ----------
  683. signal : neo.AnalogSignal
  684. The signal to integrate. If `signal` contains more than one channel,
  685. each is integrated separately.
  686. baseline : pq.Quantity or {'mean', 'median'}, optional
  687. A factor to subtract from the signal before rectification.
  688. If 'mean', the mean value of the entire `signal` is subtracted on a
  689. channel-by-channel basis.
  690. If 'median', the median value of the entire `signal` is subtracted on
  691. a channel-by-channel basis.
  692. Default: None.
  693. bin_duration : pq.Quantity, optional
  694. The length of time that each integration should span.
  695. If None, there will be only one bin spanning the entire signal
  696. duration.
  697. If `bin_duration` does not divide evenly into the signal duration, the
  698. end of the signal is padded with zeros to accomodate the final,
  699. overextending bin.
  700. Default: None.
  701. t_start: pq.Quantity, optional
  702. Time to start the algorithm.
  703. If None, starts at the beginning of `signal`.
  704. Default: None.
  705. t_stop : pq.Quantity, optional
  706. Time to end the algorithm.
  707. If None, ends at the last time of `signal`.
  708. The signal is cropped using `signal.time_slice(t_start, t_stop)` after
  709. baseline removal. Useful if you want the RAUC for a short section of
  710. the signal but want the mean or median calculation (`baseline`='mean'
  711. or `baseline`='median') to use the entire signal for better baseline
  712. estimation.
  713. Default: None.
  714. Returns
  715. -------
  716. pq.Quantity or neo.AnalogSignal
  717. If the number of bins is 1, the returned object is a scalar or
  718. vector `pq.Quantity` containing a single RAUC value for each channel.
  719. Otherwise, the returned object is a `neo.AnalogSignal` containing the
  720. RAUC(s) for each bin stored as a sample, with times corresponding to
  721. the center of each bin. The output signal will have the same number
  722. of channels as the input signal.
  723. Raises
  724. ------
  725. ValueError
  726. If `signal` is not `neo.AnalogSignal`.
  727. If `bin_duration` is not None or `pq.Quantity`.
  728. If `baseline` is not None, 'mean', 'median', or `pq.Quantity`.
  729. See Also
  730. --------
  731. neo.AnalogSignal.time_slice : how `t_start` and `t_stop` are used
  732. """
  733. if not isinstance(signal, neo.AnalogSignal):
  734. raise ValueError('Input signal is not a neo.AnalogSignal!')
  735. if baseline is None:
  736. pass
  737. elif baseline == 'mean':
  738. # subtract mean from each channel
  739. signal = signal - signal.mean(axis=0)
  740. elif baseline == 'median':
  741. # subtract median from each channel
  742. signal = signal - np.median(signal.as_quantity(), axis=0)
  743. elif isinstance(baseline, pq.Quantity):
  744. # subtract arbitrary baseline
  745. signal = signal - baseline
  746. else:
  747. raise ValueError("baseline must be either None, 'mean', 'median', or "
  748. "a Quantity. Got {}".format(baseline))
  749. # slice the signal after subtracting baseline
  750. signal = signal.time_slice(t_start, t_stop)
  751. if bin_duration is not None:
  752. # from bin duration, determine samples per bin and number of bins
  753. if isinstance(bin_duration, pq.Quantity):
  754. samples_per_bin = int(
  755. np.round(
  756. bin_duration.rescale('s') /
  757. signal.sampling_period.rescale('s')))
  758. n_bins = int(np.ceil(signal.shape[0] / samples_per_bin))
  759. else:
  760. raise ValueError("bin_duration must be a Quantity. Got {}".format(
  761. bin_duration))
  762. else:
  763. # all samples in one bin
  764. samples_per_bin = signal.shape[0]
  765. n_bins = 1
  766. # store the actual bin duration
  767. bin_duration = samples_per_bin * signal.sampling_period
  768. # reshape into equal size bins, padding the end with zeros if necessary
  769. n_channels = signal.shape[1]
  770. sig_binned = signal.as_quantity().copy()
  771. sig_binned.resize(n_bins * samples_per_bin, n_channels, refcheck=False)
  772. sig_binned = sig_binned.reshape(n_bins, samples_per_bin, n_channels)
  773. # rectify and integrate over each bin
  774. rauc = np.trapz(np.abs(sig_binned), dx=signal.sampling_period, axis=1)
  775. if n_bins == 1:
  776. # return a single value for each channel
  777. return rauc.squeeze()
  778. else:
  779. # return an AnalogSignal with times corresponding to center of each bin
  780. t_start = signal.t_start.rescale(bin_duration.units) + bin_duration / 2
  781. rauc_sig = neo.AnalogSignal(rauc, t_start=t_start,
  782. sampling_period=bin_duration)
  783. return rauc_sig
  784. def derivative(signal):
  785. """
  786. Calculate the derivative of a `neo.AnalogSignal`.
  787. Parameters
  788. ----------
  789. signal : neo.AnalogSignal
  790. The signal to differentiate. If `signal` contains more than one
  791. channel, each is differentiated separately.
  792. Returns
  793. -------
  794. derivative_sig: neo.AnalogSignal
  795. The returned object is a `neo.AnalogSignal` containing the differences
  796. between each successive sample value of the input signal divided by
  797. the sampling period. Times are centered between the successive samples
  798. of the input. The output signal will have the same number of channels
  799. as the input signal.
  800. Raises
  801. ------
  802. TypeError
  803. If `signal` is not a `neo.AnalogSignal`.
  804. """
  805. if not isinstance(signal, neo.AnalogSignal):
  806. raise TypeError('Input signal is not a neo.AnalogSignal!')
  807. derivative_sig = neo.AnalogSignal(
  808. np.diff(signal.as_quantity(), axis=0) / signal.sampling_period,
  809. t_start=signal.t_start + signal.sampling_period / 2,
  810. sampling_period=signal.sampling_period)
  811. return derivative_sig