spectral.py 72 KB

12345678910111213141516171819202122232425262728293031323334353637383940414243444546474849505152535455565758596061626364656667686970717273747576777879808182838485868788899091929394959697989910010110210310410510610710810911011111211311411511611711811912012112212312412512612712812913013113213313413513613713813914014114214314414514614714814915015115215315415515615715815916016116216316416516616716816917017117217317417517617717817918018118218318418518618718818919019119219319419519619719819920020120220320420520620720820921021121221321421521621721821922022122222322422522622722822923023123223323423523623723823924024124224324424524624724824925025125225325425525625725825926026126226326426526626726826927027127227327427527627727827928028128228328428528628728828929029129229329429529629729829930030130230330430530630730830931031131231331431531631731831932032132232332432532632732832933033133233333433533633733833934034134234334434534634734834935035135235335435535635735835936036136236336436536636736836937037137237337437537637737837938038138238338438538638738838939039139239339439539639739839940040140240340440540640740840941041141241341441541641741841942042142242342442542642742842943043143243343443543643743843944044144244344444544644744844945045145245345445545645745845946046146246346446546646746846947047147247347447547647747847948048148248348448548648748848949049149249349449549649749849950050150250350450550650750850951051151251351451551651751851952052152252352452552652752852953053153253353453553653753853954054154254354454554654754854955055155255355455555655755855956056156256356456556656756856957057157257357457557657757857958058158258358458558658758858959059159259359459559659759859960060160260360460560660760860961061161261361461561661761861962062162262362462562662762862963063163263363463563663763863964064164264364464564664764864965065165265365465565665765865966066166266366466566666766866967067167267367467567667767867968068168268368468568668768868969069169269369469569669769869970070170270370470570670770870971071171271371471571671771871972072172272372472572672772872973073173273373473573673773873974074174274374474574674774874975075175275375475575675775875976076176276376476576676776876977077177277377477577677777877978078178278378478578678778878979079179279379479579679779879980080180280380480580680780880981081181281381481581681781881982082182282382482582682782882983083183283383483583683783883984084184284384484584684784884985085185285385485585685785885986086186286386486586686786886987087187287387487587687787887988088188288388488588688788888989089189289389489589689789889990090190290390490590690790890991091191291391491591691791891992092192292392492592692792892993093193293393493593693793893994094194294394494594694794894995095195295395495595695795895996096196296396496596696796896997097197297397497597697797897998098198298398498598698798898999099199299399499599699799899910001001100210031004100510061007100810091010101110121013101410151016101710181019102010211022102310241025102610271028102910301031103210331034103510361037103810391040104110421043104410451046104710481049105010511052105310541055105610571058105910601061106210631064106510661067106810691070107110721073107410751076107710781079108010811082108310841085108610871088108910901091109210931094109510961097109810991100110111021103110411051106110711081109111011111112111311141115111611171118111911201121112211231124112511261127112811291130113111321133113411351136113711381139114011411142114311441145114611471148114911501151115211531154115511561157115811591160116111621163116411651166116711681169117011711172117311741175117611771178117911801181118211831184118511861187118811891190119111921193119411951196119711981199120012011202120312041205120612071208120912101211121212131214121512161217121812191220122112221223122412251226122712281229123012311232123312341235123612371238123912401241124212431244124512461247124812491250125112521253125412551256125712581259126012611262126312641265126612671268126912701271127212731274127512761277127812791280128112821283128412851286128712881289129012911292129312941295129612971298129913001301130213031304130513061307130813091310131113121313131413151316131713181319132013211322132313241325132613271328132913301331133213331334133513361337133813391340134113421343134413451346134713481349135013511352135313541355135613571358135913601361136213631364136513661367136813691370137113721373137413751376137713781379138013811382138313841385138613871388138913901391139213931394139513961397139813991400140114021403140414051406140714081409141014111412141314141415141614171418141914201421142214231424142514261427142814291430143114321433143414351436143714381439144014411442144314441445144614471448144914501451145214531454145514561457145814591460146114621463146414651466146714681469147014711472147314741475147614771478147914801481148214831484148514861487148814891490149114921493149414951496149714981499150015011502150315041505150615071508150915101511151215131514151515161517151815191520152115221523152415251526152715281529153015311532153315341535153615371538153915401541154215431544154515461547154815491550155115521553155415551556155715581559156015611562156315641565156615671568156915701571157215731574157515761577157815791580158115821583158415851586158715881589159015911592159315941595159615971598159916001601160216031604160516061607160816091610161116121613161416151616161716181619162016211622162316241625162616271628162916301631163216331634163516361637163816391640164116421643164416451646164716481649165016511652165316541655165616571658165916601661166216631664166516661667166816691670167116721673167416751676167716781679168016811682168316841685168616871688168916901691169216931694169516961697169816991700170117021703170417051706170717081709171017111712171317141715171617171718171917201721172217231724172517261727172817291730173117321733173417351736173717381739174017411742174317441745174617471748174917501751175217531754175517561757175817591760176117621763176417651766176717681769177017711772177317741775177617771778177917801781178217831784178517861787178817891790179117921793179417951796179717981799180018011802180318041805180618071808180918101811181218131814181518161817181818191820182118221823182418251826182718281829183018311832183318341835183618371838183918401841184218431844184518461847184818491850185118521853185418551856185718581859186018611862186318641865186618671868186918701871187218731874187518761877187818791880188118821883188418851886188718881889189018911892189318941895189618971898189919001901190219031904190519061907190819091910191119121913191419151916191719181919192019211922192319241925192619271928192919301931193219331934193519361937193819391940194119421943194419451946194719481949195019511952195319541955195619571958195919601961196219631964196519661967196819691970197119721973197419751976197719781979198019811982198319841985198619871988198919901991199219931994199519961997199819992000200120022003200420052006
  1. """Tools for spectral analysis.
  2. """
  3. from __future__ import division, print_function, absolute_import
  4. import numpy as np
  5. from scipy import fftpack
  6. from . import signaltools
  7. from .windows import get_window
  8. from ._spectral import _lombscargle
  9. from ._arraytools import const_ext, even_ext, odd_ext, zero_ext
  10. import warnings
  11. from scipy._lib.six import string_types
  12. __all__ = ['periodogram', 'welch', 'lombscargle', 'csd', 'coherence',
  13. 'spectrogram', 'stft', 'istft', 'check_COLA', 'check_NOLA']
  14. def lombscargle(x,
  15. y,
  16. freqs,
  17. precenter=False,
  18. normalize=False):
  19. """
  20. lombscargle(x, y, freqs)
  21. Computes the Lomb-Scargle periodogram.
  22. The Lomb-Scargle periodogram was developed by Lomb [1]_ and further
  23. extended by Scargle [2]_ to find, and test the significance of weak
  24. periodic signals with uneven temporal sampling.
  25. When *normalize* is False (default) the computed periodogram
  26. is unnormalized, it takes the value ``(A**2) * N/4`` for a harmonic
  27. signal with amplitude A for sufficiently large N.
  28. When *normalize* is True the computed periodogram is normalized by
  29. the residuals of the data around a constant reference model (at zero).
  30. Input arrays should be one-dimensional and will be cast to float64.
  31. Parameters
  32. ----------
  33. x : array_like
  34. Sample times.
  35. y : array_like
  36. Measurement values.
  37. freqs : array_like
  38. Angular frequencies for output periodogram.
  39. precenter : bool, optional
  40. Pre-center amplitudes by subtracting the mean.
  41. normalize : bool, optional
  42. Compute normalized periodogram.
  43. Returns
  44. -------
  45. pgram : array_like
  46. Lomb-Scargle periodogram.
  47. Raises
  48. ------
  49. ValueError
  50. If the input arrays `x` and `y` do not have the same shape.
  51. Notes
  52. -----
  53. This subroutine calculates the periodogram using a slightly
  54. modified algorithm due to Townsend [3]_ which allows the
  55. periodogram to be calculated using only a single pass through
  56. the input arrays for each frequency.
  57. The algorithm running time scales roughly as O(x * freqs) or O(N^2)
  58. for a large number of samples and frequencies.
  59. References
  60. ----------
  61. .. [1] N.R. Lomb "Least-squares frequency analysis of unequally spaced
  62. data", Astrophysics and Space Science, vol 39, pp. 447-462, 1976
  63. .. [2] J.D. Scargle "Studies in astronomical time series analysis. II -
  64. Statistical aspects of spectral analysis of unevenly spaced data",
  65. The Astrophysical Journal, vol 263, pp. 835-853, 1982
  66. .. [3] R.H.D. Townsend, "Fast calculation of the Lomb-Scargle
  67. periodogram using graphics processing units.", The Astrophysical
  68. Journal Supplement Series, vol 191, pp. 247-253, 2010
  69. See Also
  70. --------
  71. istft: Inverse Short Time Fourier Transform
  72. check_COLA: Check whether the Constant OverLap Add (COLA) constraint is met
  73. welch: Power spectral density by Welch's method
  74. spectrogram: Spectrogram by Welch's method
  75. csd: Cross spectral density by Welch's method
  76. Examples
  77. --------
  78. >>> import matplotlib.pyplot as plt
  79. First define some input parameters for the signal:
  80. >>> A = 2.
  81. >>> w = 1.
  82. >>> phi = 0.5 * np.pi
  83. >>> nin = 1000
  84. >>> nout = 100000
  85. >>> frac_points = 0.9 # Fraction of points to select
  86. Randomly select a fraction of an array with timesteps:
  87. >>> r = np.random.rand(nin)
  88. >>> x = np.linspace(0.01, 10*np.pi, nin)
  89. >>> x = x[r >= frac_points]
  90. Plot a sine wave for the selected times:
  91. >>> y = A * np.sin(w*x+phi)
  92. Define the array of frequencies for which to compute the periodogram:
  93. >>> f = np.linspace(0.01, 10, nout)
  94. Calculate Lomb-Scargle periodogram:
  95. >>> import scipy.signal as signal
  96. >>> pgram = signal.lombscargle(x, y, f, normalize=True)
  97. Now make a plot of the input data:
  98. >>> plt.subplot(2, 1, 1)
  99. >>> plt.plot(x, y, 'b+')
  100. Then plot the normalized periodogram:
  101. >>> plt.subplot(2, 1, 2)
  102. >>> plt.plot(f, pgram)
  103. >>> plt.show()
  104. """
  105. x = np.asarray(x, dtype=np.float64)
  106. y = np.asarray(y, dtype=np.float64)
  107. freqs = np.asarray(freqs, dtype=np.float64)
  108. assert x.ndim == 1
  109. assert y.ndim == 1
  110. assert freqs.ndim == 1
  111. if precenter:
  112. pgram = _lombscargle(x, y - y.mean(), freqs)
  113. else:
  114. pgram = _lombscargle(x, y, freqs)
  115. if normalize:
  116. pgram *= 2 / np.dot(y, y)
  117. return pgram
  118. def periodogram(x, fs=1.0, window='boxcar', nfft=None, detrend='constant',
  119. return_onesided=True, scaling='density', axis=-1):
  120. """
  121. Estimate power spectral density using a periodogram.
  122. Parameters
  123. ----------
  124. x : array_like
  125. Time series of measurement values
  126. fs : float, optional
  127. Sampling frequency of the `x` time series. Defaults to 1.0.
  128. window : str or tuple or array_like, optional
  129. Desired window to use. If `window` is a string or tuple, it is
  130. passed to `get_window` to generate the window values, which are
  131. DFT-even by default. See `get_window` for a list of windows and
  132. required parameters. If `window` is array_like it will be used
  133. directly as the window and its length must be nperseg. Defaults
  134. to 'boxcar'.
  135. nfft : int, optional
  136. Length of the FFT used. If `None` the length of `x` will be
  137. used.
  138. detrend : str or function or `False`, optional
  139. Specifies how to detrend each segment. If `detrend` is a
  140. string, it is passed as the `type` argument to the `detrend`
  141. function. If it is a function, it takes a segment and returns a
  142. detrended segment. If `detrend` is `False`, no detrending is
  143. done. Defaults to 'constant'.
  144. return_onesided : bool, optional
  145. If `True`, return a one-sided spectrum for real data. If
  146. `False` return a two-sided spectrum. Note that for complex
  147. data, a two-sided spectrum is always returned.
  148. scaling : { 'density', 'spectrum' }, optional
  149. Selects between computing the power spectral density ('density')
  150. where `Pxx` has units of V**2/Hz and computing the power
  151. spectrum ('spectrum') where `Pxx` has units of V**2, if `x`
  152. is measured in V and `fs` is measured in Hz. Defaults to
  153. 'density'
  154. axis : int, optional
  155. Axis along which the periodogram is computed; the default is
  156. over the last axis (i.e. ``axis=-1``).
  157. Returns
  158. -------
  159. f : ndarray
  160. Array of sample frequencies.
  161. Pxx : ndarray
  162. Power spectral density or power spectrum of `x`.
  163. Notes
  164. -----
  165. .. versionadded:: 0.12.0
  166. See Also
  167. --------
  168. welch: Estimate power spectral density using Welch's method
  169. lombscargle: Lomb-Scargle periodogram for unevenly sampled data
  170. Examples
  171. --------
  172. >>> from scipy import signal
  173. >>> import matplotlib.pyplot as plt
  174. >>> np.random.seed(1234)
  175. Generate a test signal, a 2 Vrms sine wave at 1234 Hz, corrupted by
  176. 0.001 V**2/Hz of white noise sampled at 10 kHz.
  177. >>> fs = 10e3
  178. >>> N = 1e5
  179. >>> amp = 2*np.sqrt(2)
  180. >>> freq = 1234.0
  181. >>> noise_power = 0.001 * fs / 2
  182. >>> time = np.arange(N) / fs
  183. >>> x = amp*np.sin(2*np.pi*freq*time)
  184. >>> x += np.random.normal(scale=np.sqrt(noise_power), size=time.shape)
  185. Compute and plot the power spectral density.
  186. >>> f, Pxx_den = signal.periodogram(x, fs)
  187. >>> plt.semilogy(f, Pxx_den)
  188. >>> plt.ylim([1e-7, 1e2])
  189. >>> plt.xlabel('frequency [Hz]')
  190. >>> plt.ylabel('PSD [V**2/Hz]')
  191. >>> plt.show()
  192. If we average the last half of the spectral density, to exclude the
  193. peak, we can recover the noise power on the signal.
  194. >>> np.mean(Pxx_den[25000:])
  195. 0.00099728892368242854
  196. Now compute and plot the power spectrum.
  197. >>> f, Pxx_spec = signal.periodogram(x, fs, 'flattop', scaling='spectrum')
  198. >>> plt.figure()
  199. >>> plt.semilogy(f, np.sqrt(Pxx_spec))
  200. >>> plt.ylim([1e-4, 1e1])
  201. >>> plt.xlabel('frequency [Hz]')
  202. >>> plt.ylabel('Linear spectrum [V RMS]')
  203. >>> plt.show()
  204. The peak height in the power spectrum is an estimate of the RMS
  205. amplitude.
  206. >>> np.sqrt(Pxx_spec.max())
  207. 2.0077340678640727
  208. """
  209. x = np.asarray(x)
  210. if x.size == 0:
  211. return np.empty(x.shape), np.empty(x.shape)
  212. if window is None:
  213. window = 'boxcar'
  214. if nfft is None:
  215. nperseg = x.shape[axis]
  216. elif nfft == x.shape[axis]:
  217. nperseg = nfft
  218. elif nfft > x.shape[axis]:
  219. nperseg = x.shape[axis]
  220. elif nfft < x.shape[axis]:
  221. s = [np.s_[:]]*len(x.shape)
  222. s[axis] = np.s_[:nfft]
  223. x = x[tuple(s)]
  224. nperseg = nfft
  225. nfft = None
  226. return welch(x, fs=fs, window=window, nperseg=nperseg, noverlap=0,
  227. nfft=nfft, detrend=detrend, return_onesided=return_onesided,
  228. scaling=scaling, axis=axis)
  229. def welch(x, fs=1.0, window='hann', nperseg=None, noverlap=None, nfft=None,
  230. detrend='constant', return_onesided=True, scaling='density',
  231. axis=-1, average='mean'):
  232. r"""
  233. Estimate power spectral density using Welch's method.
  234. Welch's method [1]_ computes an estimate of the power spectral
  235. density by dividing the data into overlapping segments, computing a
  236. modified periodogram for each segment and averaging the
  237. periodograms.
  238. Parameters
  239. ----------
  240. x : array_like
  241. Time series of measurement values
  242. fs : float, optional
  243. Sampling frequency of the `x` time series. Defaults to 1.0.
  244. window : str or tuple or array_like, optional
  245. Desired window to use. If `window` is a string or tuple, it is
  246. passed to `get_window` to generate the window values, which are
  247. DFT-even by default. See `get_window` for a list of windows and
  248. required parameters. If `window` is array_like it will be used
  249. directly as the window and its length must be nperseg. Defaults
  250. to a Hann window.
  251. nperseg : int, optional
  252. Length of each segment. Defaults to None, but if window is str or
  253. tuple, is set to 256, and if window is array_like, is set to the
  254. length of the window.
  255. noverlap : int, optional
  256. Number of points to overlap between segments. If `None`,
  257. ``noverlap = nperseg // 2``. Defaults to `None`.
  258. nfft : int, optional
  259. Length of the FFT used, if a zero padded FFT is desired. If
  260. `None`, the FFT length is `nperseg`. Defaults to `None`.
  261. detrend : str or function or `False`, optional
  262. Specifies how to detrend each segment. If `detrend` is a
  263. string, it is passed as the `type` argument to the `detrend`
  264. function. If it is a function, it takes a segment and returns a
  265. detrended segment. If `detrend` is `False`, no detrending is
  266. done. Defaults to 'constant'.
  267. return_onesided : bool, optional
  268. If `True`, return a one-sided spectrum for real data. If
  269. `False` return a two-sided spectrum. Note that for complex
  270. data, a two-sided spectrum is always returned.
  271. scaling : { 'density', 'spectrum' }, optional
  272. Selects between computing the power spectral density ('density')
  273. where `Pxx` has units of V**2/Hz and computing the power
  274. spectrum ('spectrum') where `Pxx` has units of V**2, if `x`
  275. is measured in V and `fs` is measured in Hz. Defaults to
  276. 'density'
  277. axis : int, optional
  278. Axis along which the periodogram is computed; the default is
  279. over the last axis (i.e. ``axis=-1``).
  280. average : { 'mean', 'median' }, optional
  281. Method to use when averaging periodograms. Defaults to 'mean'.
  282. .. versionadded:: 1.2.0
  283. Returns
  284. -------
  285. f : ndarray
  286. Array of sample frequencies.
  287. Pxx : ndarray
  288. Power spectral density or power spectrum of x.
  289. See Also
  290. --------
  291. periodogram: Simple, optionally modified periodogram
  292. lombscargle: Lomb-Scargle periodogram for unevenly sampled data
  293. Notes
  294. -----
  295. An appropriate amount of overlap will depend on the choice of window
  296. and on your requirements. For the default Hann window an overlap of
  297. 50% is a reasonable trade off between accurately estimating the
  298. signal power, while not over counting any of the data. Narrower
  299. windows may require a larger overlap.
  300. If `noverlap` is 0, this method is equivalent to Bartlett's method
  301. [2]_.
  302. .. versionadded:: 0.12.0
  303. References
  304. ----------
  305. .. [1] P. Welch, "The use of the fast Fourier transform for the
  306. estimation of power spectra: A method based on time averaging
  307. over short, modified periodograms", IEEE Trans. Audio
  308. Electroacoust. vol. 15, pp. 70-73, 1967.
  309. .. [2] M.S. Bartlett, "Periodogram Analysis and Continuous Spectra",
  310. Biometrika, vol. 37, pp. 1-16, 1950.
  311. Examples
  312. --------
  313. >>> from scipy import signal
  314. >>> import matplotlib.pyplot as plt
  315. >>> np.random.seed(1234)
  316. Generate a test signal, a 2 Vrms sine wave at 1234 Hz, corrupted by
  317. 0.001 V**2/Hz of white noise sampled at 10 kHz.
  318. >>> fs = 10e3
  319. >>> N = 1e5
  320. >>> amp = 2*np.sqrt(2)
  321. >>> freq = 1234.0
  322. >>> noise_power = 0.001 * fs / 2
  323. >>> time = np.arange(N) / fs
  324. >>> x = amp*np.sin(2*np.pi*freq*time)
  325. >>> x += np.random.normal(scale=np.sqrt(noise_power), size=time.shape)
  326. Compute and plot the power spectral density.
  327. >>> f, Pxx_den = signal.welch(x, fs, nperseg=1024)
  328. >>> plt.semilogy(f, Pxx_den)
  329. >>> plt.ylim([0.5e-3, 1])
  330. >>> plt.xlabel('frequency [Hz]')
  331. >>> plt.ylabel('PSD [V**2/Hz]')
  332. >>> plt.show()
  333. If we average the last half of the spectral density, to exclude the
  334. peak, we can recover the noise power on the signal.
  335. >>> np.mean(Pxx_den[256:])
  336. 0.0009924865443739191
  337. Now compute and plot the power spectrum.
  338. >>> f, Pxx_spec = signal.welch(x, fs, 'flattop', 1024, scaling='spectrum')
  339. >>> plt.figure()
  340. >>> plt.semilogy(f, np.sqrt(Pxx_spec))
  341. >>> plt.xlabel('frequency [Hz]')
  342. >>> plt.ylabel('Linear spectrum [V RMS]')
  343. >>> plt.show()
  344. The peak height in the power spectrum is an estimate of the RMS
  345. amplitude.
  346. >>> np.sqrt(Pxx_spec.max())
  347. 2.0077340678640727
  348. If we now introduce a discontinuity in the signal, by increasing the
  349. amplitude of a small portion of the signal by 50, we can see the
  350. corruption of the mean average power spectral density, but using a
  351. median average better estimates the normal behaviour.
  352. >>> x[int(N//2):int(N//2)+10] *= 50.
  353. >>> f, Pxx_den = signal.welch(x, fs, nperseg=1024)
  354. >>> f_med, Pxx_den_med = signal.welch(x, fs, nperseg=1024, average='median')
  355. >>> plt.semilogy(f, Pxx_den, label='mean')
  356. >>> plt.semilogy(f_med, Pxx_den_med, label='median')
  357. >>> plt.ylim([0.5e-3, 1])
  358. >>> plt.xlabel('frequency [Hz]')
  359. >>> plt.ylabel('PSD [V**2/Hz]')
  360. >>> plt.legend()
  361. >>> plt.show()
  362. """
  363. freqs, Pxx = csd(x, x, fs=fs, window=window, nperseg=nperseg,
  364. noverlap=noverlap, nfft=nfft, detrend=detrend,
  365. return_onesided=return_onesided, scaling=scaling,
  366. axis=axis, average=average)
  367. return freqs, Pxx.real
  368. def csd(x, y, fs=1.0, window='hann', nperseg=None, noverlap=None, nfft=None,
  369. detrend='constant', return_onesided=True, scaling='density',
  370. axis=-1, average='mean'):
  371. r"""
  372. Estimate the cross power spectral density, Pxy, using Welch's
  373. method.
  374. Parameters
  375. ----------
  376. x : array_like
  377. Time series of measurement values
  378. y : array_like
  379. Time series of measurement values
  380. fs : float, optional
  381. Sampling frequency of the `x` and `y` time series. Defaults
  382. to 1.0.
  383. window : str or tuple or array_like, optional
  384. Desired window to use. If `window` is a string or tuple, it is
  385. passed to `get_window` to generate the window values, which are
  386. DFT-even by default. See `get_window` for a list of windows and
  387. required parameters. If `window` is array_like it will be used
  388. directly as the window and its length must be nperseg. Defaults
  389. to a Hann window.
  390. nperseg : int, optional
  391. Length of each segment. Defaults to None, but if window is str or
  392. tuple, is set to 256, and if window is array_like, is set to the
  393. length of the window.
  394. noverlap: int, optional
  395. Number of points to overlap between segments. If `None`,
  396. ``noverlap = nperseg // 2``. Defaults to `None`.
  397. nfft : int, optional
  398. Length of the FFT used, if a zero padded FFT is desired. If
  399. `None`, the FFT length is `nperseg`. Defaults to `None`.
  400. detrend : str or function or `False`, optional
  401. Specifies how to detrend each segment. If `detrend` is a
  402. string, it is passed as the `type` argument to the `detrend`
  403. function. If it is a function, it takes a segment and returns a
  404. detrended segment. If `detrend` is `False`, no detrending is
  405. done. Defaults to 'constant'.
  406. return_onesided : bool, optional
  407. If `True`, return a one-sided spectrum for real data. If
  408. `False` return a two-sided spectrum. Note that for complex
  409. data, a two-sided spectrum is always returned.
  410. scaling : { 'density', 'spectrum' }, optional
  411. Selects between computing the cross spectral density ('density')
  412. where `Pxy` has units of V**2/Hz and computing the cross spectrum
  413. ('spectrum') where `Pxy` has units of V**2, if `x` and `y` are
  414. measured in V and `fs` is measured in Hz. Defaults to 'density'
  415. axis : int, optional
  416. Axis along which the CSD is computed for both inputs; the
  417. default is over the last axis (i.e. ``axis=-1``).
  418. average : { 'mean', 'median' }, optional
  419. Method to use when averaging periodograms. Defaults to 'mean'.
  420. .. versionadded:: 1.2.0
  421. Returns
  422. -------
  423. f : ndarray
  424. Array of sample frequencies.
  425. Pxy : ndarray
  426. Cross spectral density or cross power spectrum of x,y.
  427. See Also
  428. --------
  429. periodogram: Simple, optionally modified periodogram
  430. lombscargle: Lomb-Scargle periodogram for unevenly sampled data
  431. welch: Power spectral density by Welch's method. [Equivalent to
  432. csd(x,x)]
  433. coherence: Magnitude squared coherence by Welch's method.
  434. Notes
  435. --------
  436. By convention, Pxy is computed with the conjugate FFT of X
  437. multiplied by the FFT of Y.
  438. If the input series differ in length, the shorter series will be
  439. zero-padded to match.
  440. An appropriate amount of overlap will depend on the choice of window
  441. and on your requirements. For the default Hann window an overlap of
  442. 50% is a reasonable trade off between accurately estimating the
  443. signal power, while not over counting any of the data. Narrower
  444. windows may require a larger overlap.
  445. .. versionadded:: 0.16.0
  446. References
  447. ----------
  448. .. [1] P. Welch, "The use of the fast Fourier transform for the
  449. estimation of power spectra: A method based on time averaging
  450. over short, modified periodograms", IEEE Trans. Audio
  451. Electroacoust. vol. 15, pp. 70-73, 1967.
  452. .. [2] Rabiner, Lawrence R., and B. Gold. "Theory and Application of
  453. Digital Signal Processing" Prentice-Hall, pp. 414-419, 1975
  454. Examples
  455. --------
  456. >>> from scipy import signal
  457. >>> import matplotlib.pyplot as plt
  458. Generate two test signals with some common features.
  459. >>> fs = 10e3
  460. >>> N = 1e5
  461. >>> amp = 20
  462. >>> freq = 1234.0
  463. >>> noise_power = 0.001 * fs / 2
  464. >>> time = np.arange(N) / fs
  465. >>> b, a = signal.butter(2, 0.25, 'low')
  466. >>> x = np.random.normal(scale=np.sqrt(noise_power), size=time.shape)
  467. >>> y = signal.lfilter(b, a, x)
  468. >>> x += amp*np.sin(2*np.pi*freq*time)
  469. >>> y += np.random.normal(scale=0.1*np.sqrt(noise_power), size=time.shape)
  470. Compute and plot the magnitude of the cross spectral density.
  471. >>> f, Pxy = signal.csd(x, y, fs, nperseg=1024)
  472. >>> plt.semilogy(f, np.abs(Pxy))
  473. >>> plt.xlabel('frequency [Hz]')
  474. >>> plt.ylabel('CSD [V**2/Hz]')
  475. >>> plt.show()
  476. """
  477. freqs, _, Pxy = _spectral_helper(x, y, fs, window, nperseg, noverlap, nfft,
  478. detrend, return_onesided, scaling, axis,
  479. mode='psd')
  480. # Average over windows.
  481. if len(Pxy.shape) >= 2 and Pxy.size > 0:
  482. if Pxy.shape[-1] > 1:
  483. if average == 'median':
  484. Pxy = np.median(Pxy, axis=-1) / _median_bias(Pxy.shape[-1])
  485. elif average == 'mean':
  486. Pxy = Pxy.mean(axis=-1)
  487. else:
  488. raise ValueError('average must be "median" or "mean", got %s'
  489. % (average,))
  490. else:
  491. Pxy = np.reshape(Pxy, Pxy.shape[:-1])
  492. return freqs, Pxy
  493. def spectrogram(x, fs=1.0, window=('tukey', .25), nperseg=None, noverlap=None,
  494. nfft=None, detrend='constant', return_onesided=True,
  495. scaling='density', axis=-1, mode='psd'):
  496. """
  497. Compute a spectrogram with consecutive Fourier transforms.
  498. Spectrograms can be used as a way of visualizing the change of a
  499. nonstationary signal's frequency content over time.
  500. Parameters
  501. ----------
  502. x : array_like
  503. Time series of measurement values
  504. fs : float, optional
  505. Sampling frequency of the `x` time series. Defaults to 1.0.
  506. window : str or tuple or array_like, optional
  507. Desired window to use. If `window` is a string or tuple, it is
  508. passed to `get_window` to generate the window values, which are
  509. DFT-even by default. See `get_window` for a list of windows and
  510. required parameters. If `window` is array_like it will be used
  511. directly as the window and its length must be nperseg.
  512. Defaults to a Tukey window with shape parameter of 0.25.
  513. nperseg : int, optional
  514. Length of each segment. Defaults to None, but if window is str or
  515. tuple, is set to 256, and if window is array_like, is set to the
  516. length of the window.
  517. noverlap : int, optional
  518. Number of points to overlap between segments. If `None`,
  519. ``noverlap = nperseg // 8``. Defaults to `None`.
  520. nfft : int, optional
  521. Length of the FFT used, if a zero padded FFT is desired. If
  522. `None`, the FFT length is `nperseg`. Defaults to `None`.
  523. detrend : str or function or `False`, optional
  524. Specifies how to detrend each segment. If `detrend` is a
  525. string, it is passed as the `type` argument to the `detrend`
  526. function. If it is a function, it takes a segment and returns a
  527. detrended segment. If `detrend` is `False`, no detrending is
  528. done. Defaults to 'constant'.
  529. return_onesided : bool, optional
  530. If `True`, return a one-sided spectrum for real data. If
  531. `False` return a two-sided spectrum. Note that for complex
  532. data, a two-sided spectrum is always returned.
  533. scaling : { 'density', 'spectrum' }, optional
  534. Selects between computing the power spectral density ('density')
  535. where `Sxx` has units of V**2/Hz and computing the power
  536. spectrum ('spectrum') where `Sxx` has units of V**2, if `x`
  537. is measured in V and `fs` is measured in Hz. Defaults to
  538. 'density'.
  539. axis : int, optional
  540. Axis along which the spectrogram is computed; the default is over
  541. the last axis (i.e. ``axis=-1``).
  542. mode : str, optional
  543. Defines what kind of return values are expected. Options are
  544. ['psd', 'complex', 'magnitude', 'angle', 'phase']. 'complex' is
  545. equivalent to the output of `stft` with no padding or boundary
  546. extension. 'magnitude' returns the absolute magnitude of the
  547. STFT. 'angle' and 'phase' return the complex angle of the STFT,
  548. with and without unwrapping, respectively.
  549. Returns
  550. -------
  551. f : ndarray
  552. Array of sample frequencies.
  553. t : ndarray
  554. Array of segment times.
  555. Sxx : ndarray
  556. Spectrogram of x. By default, the last axis of Sxx corresponds
  557. to the segment times.
  558. See Also
  559. --------
  560. periodogram: Simple, optionally modified periodogram
  561. lombscargle: Lomb-Scargle periodogram for unevenly sampled data
  562. welch: Power spectral density by Welch's method.
  563. csd: Cross spectral density by Welch's method.
  564. Notes
  565. -----
  566. An appropriate amount of overlap will depend on the choice of window
  567. and on your requirements. In contrast to welch's method, where the
  568. entire data stream is averaged over, one may wish to use a smaller
  569. overlap (or perhaps none at all) when computing a spectrogram, to
  570. maintain some statistical independence between individual segments.
  571. It is for this reason that the default window is a Tukey window with
  572. 1/8th of a window's length overlap at each end.
  573. .. versionadded:: 0.16.0
  574. References
  575. ----------
  576. .. [1] Oppenheim, Alan V., Ronald W. Schafer, John R. Buck
  577. "Discrete-Time Signal Processing", Prentice Hall, 1999.
  578. Examples
  579. --------
  580. >>> from scipy import signal
  581. >>> import matplotlib.pyplot as plt
  582. Generate a test signal, a 2 Vrms sine wave whose frequency is slowly
  583. modulated around 3kHz, corrupted by white noise of exponentially
  584. decreasing magnitude sampled at 10 kHz.
  585. >>> fs = 10e3
  586. >>> N = 1e5
  587. >>> amp = 2 * np.sqrt(2)
  588. >>> noise_power = 0.01 * fs / 2
  589. >>> time = np.arange(N) / float(fs)
  590. >>> mod = 500*np.cos(2*np.pi*0.25*time)
  591. >>> carrier = amp * np.sin(2*np.pi*3e3*time + mod)
  592. >>> noise = np.random.normal(scale=np.sqrt(noise_power), size=time.shape)
  593. >>> noise *= np.exp(-time/5)
  594. >>> x = carrier + noise
  595. Compute and plot the spectrogram.
  596. >>> f, t, Sxx = signal.spectrogram(x, fs)
  597. >>> plt.pcolormesh(t, f, Sxx)
  598. >>> plt.ylabel('Frequency [Hz]')
  599. >>> plt.xlabel('Time [sec]')
  600. >>> plt.show()
  601. Note, if using output that is not one sided, then use the following:
  602. >>> f, t, Sxx = signal.spectrogram(x, fs, return_onesided=False)
  603. >>> plt.pcolormesh(t, np.fft.fftshift(f), np.fft.fftshift(Sxx, axes=0))
  604. >>> plt.ylabel('Frequency [Hz]')
  605. >>> plt.xlabel('Time [sec]')
  606. >>> plt.show()
  607. """
  608. modelist = ['psd', 'complex', 'magnitude', 'angle', 'phase']
  609. if mode not in modelist:
  610. raise ValueError('unknown value for mode {}, must be one of {}'
  611. .format(mode, modelist))
  612. # need to set default for nperseg before setting default for noverlap below
  613. window, nperseg = _triage_segments(window, nperseg,
  614. input_length=x.shape[axis])
  615. # Less overlap than welch, so samples are more statisically independent
  616. if noverlap is None:
  617. noverlap = nperseg // 8
  618. if mode == 'psd':
  619. freqs, time, Sxx = _spectral_helper(x, x, fs, window, nperseg,
  620. noverlap, nfft, detrend,
  621. return_onesided, scaling, axis,
  622. mode='psd')
  623. else:
  624. freqs, time, Sxx = _spectral_helper(x, x, fs, window, nperseg,
  625. noverlap, nfft, detrend,
  626. return_onesided, scaling, axis,
  627. mode='stft')
  628. if mode == 'magnitude':
  629. Sxx = np.abs(Sxx)
  630. elif mode in ['angle', 'phase']:
  631. Sxx = np.angle(Sxx)
  632. if mode == 'phase':
  633. # Sxx has one additional dimension for time strides
  634. if axis < 0:
  635. axis -= 1
  636. Sxx = np.unwrap(Sxx, axis=axis)
  637. # mode =='complex' is same as `stft`, doesn't need modification
  638. return freqs, time, Sxx
  639. def check_COLA(window, nperseg, noverlap, tol=1e-10):
  640. r"""
  641. Check whether the Constant OverLap Add (COLA) constraint is met
  642. Parameters
  643. ----------
  644. window : str or tuple or array_like
  645. Desired window to use. If `window` is a string or tuple, it is
  646. passed to `get_window` to generate the window values, which are
  647. DFT-even by default. See `get_window` for a list of windows and
  648. required parameters. If `window` is array_like it will be used
  649. directly as the window and its length must be nperseg.
  650. nperseg : int
  651. Length of each segment.
  652. noverlap : int
  653. Number of points to overlap between segments.
  654. tol : float, optional
  655. The allowed variance of a bin's weighted sum from the median bin
  656. sum.
  657. Returns
  658. -------
  659. verdict : bool
  660. `True` if chosen combination satisfies COLA within `tol`,
  661. `False` otherwise
  662. See Also
  663. --------
  664. check_NOLA: Check whether the Nonzero Overlap Add (NOLA) constraint is met
  665. stft: Short Time Fourier Transform
  666. istft: Inverse Short Time Fourier Transform
  667. Notes
  668. -----
  669. In order to enable inversion of an STFT via the inverse STFT in
  670. `istft`, it is sufficient that the signal windowing obeys the constraint of
  671. "Constant OverLap Add" (COLA). This ensures that every point in the input
  672. data is equally weighted, thereby avoiding aliasing and allowing full
  673. reconstruction.
  674. Some examples of windows that satisfy COLA:
  675. - Rectangular window at overlap of 0, 1/2, 2/3, 3/4, ...
  676. - Bartlett window at overlap of 1/2, 3/4, 5/6, ...
  677. - Hann window at 1/2, 2/3, 3/4, ...
  678. - Any Blackman family window at 2/3 overlap
  679. - Any window with ``noverlap = nperseg-1``
  680. A very comprehensive list of other windows may be found in [2]_,
  681. wherein the COLA condition is satisfied when the "Amplitude
  682. Flatness" is unity.
  683. .. versionadded:: 0.19.0
  684. References
  685. ----------
  686. .. [1] Julius O. Smith III, "Spectral Audio Signal Processing", W3K
  687. Publishing, 2011,ISBN 978-0-9745607-3-1.
  688. .. [2] G. Heinzel, A. Ruediger and R. Schilling, "Spectrum and
  689. spectral density estimation by the Discrete Fourier transform
  690. (DFT), including a comprehensive list of window functions and
  691. some new at-top windows", 2002,
  692. http://hdl.handle.net/11858/00-001M-0000-0013-557A-5
  693. Examples
  694. --------
  695. >>> from scipy import signal
  696. Confirm COLA condition for rectangular window of 75% (3/4) overlap:
  697. >>> signal.check_COLA(signal.boxcar(100), 100, 75)
  698. True
  699. COLA is not true for 25% (1/4) overlap, though:
  700. >>> signal.check_COLA(signal.boxcar(100), 100, 25)
  701. False
  702. "Symmetrical" Hann window (for filter design) is not COLA:
  703. >>> signal.check_COLA(signal.hann(120, sym=True), 120, 60)
  704. False
  705. "Periodic" or "DFT-even" Hann window (for FFT analysis) is COLA for
  706. overlap of 1/2, 2/3, 3/4, etc.:
  707. >>> signal.check_COLA(signal.hann(120, sym=False), 120, 60)
  708. True
  709. >>> signal.check_COLA(signal.hann(120, sym=False), 120, 80)
  710. True
  711. >>> signal.check_COLA(signal.hann(120, sym=False), 120, 90)
  712. True
  713. """
  714. nperseg = int(nperseg)
  715. if nperseg < 1:
  716. raise ValueError('nperseg must be a positive integer')
  717. if noverlap >= nperseg:
  718. raise ValueError('noverlap must be less than nperseg.')
  719. noverlap = int(noverlap)
  720. if isinstance(window, string_types) or type(window) is tuple:
  721. win = get_window(window, nperseg)
  722. else:
  723. win = np.asarray(window)
  724. if len(win.shape) != 1:
  725. raise ValueError('window must be 1-D')
  726. if win.shape[0] != nperseg:
  727. raise ValueError('window must have length of nperseg')
  728. step = nperseg - noverlap
  729. binsums = sum(win[ii*step:(ii+1)*step] for ii in range(nperseg//step))
  730. if nperseg % step != 0:
  731. binsums[:nperseg % step] += win[-(nperseg % step):]
  732. deviation = binsums - np.median(binsums)
  733. return np.max(np.abs(deviation)) < tol
  734. def check_NOLA(window, nperseg, noverlap, tol=1e-10):
  735. r"""
  736. Check whether the Nonzero Overlap Add (NOLA) constraint is met
  737. Parameters
  738. ----------
  739. window : str or tuple or array_like
  740. Desired window to use. If `window` is a string or tuple, it is
  741. passed to `get_window` to generate the window values, which are
  742. DFT-even by default. See `get_window` for a list of windows and
  743. required parameters. If `window` is array_like it will be used
  744. directly as the window and its length must be nperseg.
  745. nperseg : int
  746. Length of each segment.
  747. noverlap : int
  748. Number of points to overlap between segments.
  749. tol : float, optional
  750. The allowed variance of a bin's weighted sum from the median bin
  751. sum.
  752. Returns
  753. -------
  754. verdict : bool
  755. `True` if chosen combination satisfies the NOLA constraint within
  756. `tol`, `False` otherwise
  757. See Also
  758. --------
  759. check_COLA: Check whether the Constant OverLap Add (COLA) constraint is met
  760. stft: Short Time Fourier Transform
  761. istft: Inverse Short Time Fourier Transform
  762. Notes
  763. -----
  764. In order to enable inversion of an STFT via the inverse STFT in
  765. `istft`, the signal windowing must obey the constraint of "nonzero
  766. overlap add" (NOLA):
  767. .. math:: \sum_{t}w^{2}[n-tH] \ne 0
  768. for all :math:`n`, where :math:`w` is the window function, :math:`t` is the
  769. frame index, and :math:`H` is the hop size (:math:`H` = `nperseg` -
  770. `noverlap`).
  771. This ensures that the normalization factors in the denominator of the
  772. overlap-add inversion equation are not zero. Only very pathological windows
  773. will fail the NOLA constraint.
  774. .. versionadded:: 1.2.0
  775. References
  776. ----------
  777. .. [1] Julius O. Smith III, "Spectral Audio Signal Processing", W3K
  778. Publishing, 2011,ISBN 978-0-9745607-3-1.
  779. .. [2] G. Heinzel, A. Ruediger and R. Schilling, "Spectrum and
  780. spectral density estimation by the Discrete Fourier transform
  781. (DFT), including a comprehensive list of window functions and
  782. some new at-top windows", 2002,
  783. http://hdl.handle.net/11858/00-001M-0000-0013-557A-5
  784. Examples
  785. --------
  786. >>> from scipy import signal
  787. Confirm NOLA condition for rectangular window of 75% (3/4) overlap:
  788. >>> signal.check_NOLA(signal.boxcar(100), 100, 75)
  789. True
  790. NOLA is also true for 25% (1/4) overlap:
  791. >>> signal.check_NOLA(signal.boxcar(100), 100, 25)
  792. True
  793. "Symmetrical" Hann window (for filter design) is also NOLA:
  794. >>> signal.check_NOLA(signal.hann(120, sym=True), 120, 60)
  795. True
  796. As long as there is overlap, it takes quite a pathological window to fail
  797. NOLA:
  798. >>> w = np.ones(64, dtype="float")
  799. >>> w[::2] = 0
  800. >>> signal.check_NOLA(w, 64, 32)
  801. False
  802. If there is not enough overlap, a window with zeros at the ends will not
  803. work:
  804. >>> signal.check_NOLA(signal.hann(64), 64, 0)
  805. False
  806. >>> signal.check_NOLA(signal.hann(64), 64, 1)
  807. False
  808. >>> signal.check_NOLA(signal.hann(64), 64, 2)
  809. True
  810. """
  811. nperseg = int(nperseg)
  812. if nperseg < 1:
  813. raise ValueError('nperseg must be a positive integer')
  814. if noverlap >= nperseg:
  815. raise ValueError('noverlap must be less than nperseg')
  816. if noverlap < 0:
  817. raise ValueError('noverlap must be a nonnegative integer')
  818. noverlap = int(noverlap)
  819. if isinstance(window, string_types) or type(window) is tuple:
  820. win = get_window(window, nperseg)
  821. else:
  822. win = np.asarray(window)
  823. if len(win.shape) != 1:
  824. raise ValueError('window must be 1-D')
  825. if win.shape[0] != nperseg:
  826. raise ValueError('window must have length of nperseg')
  827. step = nperseg - noverlap
  828. binsums = sum(win[ii*step:(ii+1)*step]**2 for ii in range(nperseg//step))
  829. if nperseg % step != 0:
  830. binsums[:nperseg % step] += win[-(nperseg % step):]**2
  831. return np.min(binsums) > tol
  832. def stft(x, fs=1.0, window='hann', nperseg=256, noverlap=None, nfft=None,
  833. detrend=False, return_onesided=True, boundary='zeros', padded=True,
  834. axis=-1):
  835. r"""
  836. Compute the Short Time Fourier Transform (STFT).
  837. STFTs can be used as a way of quantifying the change of a
  838. nonstationary signal's frequency and phase content over time.
  839. Parameters
  840. ----------
  841. x : array_like
  842. Time series of measurement values
  843. fs : float, optional
  844. Sampling frequency of the `x` time series. Defaults to 1.0.
  845. window : str or tuple or array_like, optional
  846. Desired window to use. If `window` is a string or tuple, it is
  847. passed to `get_window` to generate the window values, which are
  848. DFT-even by default. See `get_window` for a list of windows and
  849. required parameters. If `window` is array_like it will be used
  850. directly as the window and its length must be nperseg. Defaults
  851. to a Hann window.
  852. nperseg : int, optional
  853. Length of each segment. Defaults to 256.
  854. noverlap : int, optional
  855. Number of points to overlap between segments. If `None`,
  856. ``noverlap = nperseg // 2``. Defaults to `None`. When
  857. specified, the COLA constraint must be met (see Notes below).
  858. nfft : int, optional
  859. Length of the FFT used, if a zero padded FFT is desired. If
  860. `None`, the FFT length is `nperseg`. Defaults to `None`.
  861. detrend : str or function or `False`, optional
  862. Specifies how to detrend each segment. If `detrend` is a
  863. string, it is passed as the `type` argument to the `detrend`
  864. function. If it is a function, it takes a segment and returns a
  865. detrended segment. If `detrend` is `False`, no detrending is
  866. done. Defaults to `False`.
  867. return_onesided : bool, optional
  868. If `True`, return a one-sided spectrum for real data. If
  869. `False` return a two-sided spectrum. Note that for complex
  870. data, a two-sided spectrum is always returned. Defaults to
  871. `True`.
  872. boundary : str or None, optional
  873. Specifies whether the input signal is extended at both ends, and
  874. how to generate the new values, in order to center the first
  875. windowed segment on the first input point. This has the benefit
  876. of enabling reconstruction of the first input point when the
  877. employed window function starts at zero. Valid options are
  878. ``['even', 'odd', 'constant', 'zeros', None]``. Defaults to
  879. 'zeros', for zero padding extension. I.e. ``[1, 2, 3, 4]`` is
  880. extended to ``[0, 1, 2, 3, 4, 0]`` for ``nperseg=3``.
  881. padded : bool, optional
  882. Specifies whether the input signal is zero-padded at the end to
  883. make the signal fit exactly into an integer number of window
  884. segments, so that all of the signal is included in the output.
  885. Defaults to `True`. Padding occurs after boundary extension, if
  886. `boundary` is not `None`, and `padded` is `True`, as is the
  887. default.
  888. axis : int, optional
  889. Axis along which the STFT is computed; the default is over the
  890. last axis (i.e. ``axis=-1``).
  891. Returns
  892. -------
  893. f : ndarray
  894. Array of sample frequencies.
  895. t : ndarray
  896. Array of segment times.
  897. Zxx : ndarray
  898. STFT of `x`. By default, the last axis of `Zxx` corresponds
  899. to the segment times.
  900. See Also
  901. --------
  902. istft: Inverse Short Time Fourier Transform
  903. check_COLA: Check whether the Constant OverLap Add (COLA) constraint
  904. is met
  905. check_NOLA: Check whether the Nonzero Overlap Add (NOLA) constraint is met
  906. welch: Power spectral density by Welch's method.
  907. spectrogram: Spectrogram by Welch's method.
  908. csd: Cross spectral density by Welch's method.
  909. lombscargle: Lomb-Scargle periodogram for unevenly sampled data
  910. Notes
  911. -----
  912. In order to enable inversion of an STFT via the inverse STFT in
  913. `istft`, the signal windowing must obey the constraint of "Nonzero
  914. OverLap Add" (NOLA), and the input signal must have complete
  915. windowing coverage (i.e. ``(x.shape[axis] - nperseg) %
  916. (nperseg-noverlap) == 0``). The `padded` argument may be used to
  917. accomplish this.
  918. Given a time-domain signal :math:`x[n]`, a window :math:`w[n]`, and a hop
  919. size :math:`H` = `nperseg - noverlap`, the windowed frame at time index
  920. :math:`t` is given by
  921. .. math:: x_{t}[n]=x[n]w[n-tH]
  922. The overlap-add (OLA) reconstruction equation is given by
  923. .. math:: x[n]=\frac{\sum_{t}x_{t}[n]w[n-tH]}{\sum_{t}w^{2}[n-tH]}
  924. The NOLA constraint ensures that every normalization term that appears
  925. in the denomimator of the OLA reconstruction equation is nonzero. Whether a
  926. choice of `window`, `nperseg`, and `noverlap` satisfy this constraint can
  927. be tested with `check_NOLA`.
  928. .. versionadded:: 0.19.0
  929. References
  930. ----------
  931. .. [1] Oppenheim, Alan V., Ronald W. Schafer, John R. Buck
  932. "Discrete-Time Signal Processing", Prentice Hall, 1999.
  933. .. [2] Daniel W. Griffin, Jae S. Lim "Signal Estimation from
  934. Modified Short-Time Fourier Transform", IEEE 1984,
  935. 10.1109/TASSP.1984.1164317
  936. Examples
  937. --------
  938. >>> from scipy import signal
  939. >>> import matplotlib.pyplot as plt
  940. Generate a test signal, a 2 Vrms sine wave whose frequency is slowly
  941. modulated around 3kHz, corrupted by white noise of exponentially
  942. decreasing magnitude sampled at 10 kHz.
  943. >>> fs = 10e3
  944. >>> N = 1e5
  945. >>> amp = 2 * np.sqrt(2)
  946. >>> noise_power = 0.01 * fs / 2
  947. >>> time = np.arange(N) / float(fs)
  948. >>> mod = 500*np.cos(2*np.pi*0.25*time)
  949. >>> carrier = amp * np.sin(2*np.pi*3e3*time + mod)
  950. >>> noise = np.random.normal(scale=np.sqrt(noise_power),
  951. ... size=time.shape)
  952. >>> noise *= np.exp(-time/5)
  953. >>> x = carrier + noise
  954. Compute and plot the STFT's magnitude.
  955. >>> f, t, Zxx = signal.stft(x, fs, nperseg=1000)
  956. >>> plt.pcolormesh(t, f, np.abs(Zxx), vmin=0, vmax=amp)
  957. >>> plt.title('STFT Magnitude')
  958. >>> plt.ylabel('Frequency [Hz]')
  959. >>> plt.xlabel('Time [sec]')
  960. >>> plt.show()
  961. """
  962. freqs, time, Zxx = _spectral_helper(x, x, fs, window, nperseg, noverlap,
  963. nfft, detrend, return_onesided,
  964. scaling='spectrum', axis=axis,
  965. mode='stft', boundary=boundary,
  966. padded=padded)
  967. return freqs, time, Zxx
  968. def istft(Zxx, fs=1.0, window='hann', nperseg=None, noverlap=None, nfft=None,
  969. input_onesided=True, boundary=True, time_axis=-1, freq_axis=-2):
  970. r"""
  971. Perform the inverse Short Time Fourier transform (iSTFT).
  972. Parameters
  973. ----------
  974. Zxx : array_like
  975. STFT of the signal to be reconstructed. If a purely real array
  976. is passed, it will be cast to a complex data type.
  977. fs : float, optional
  978. Sampling frequency of the time series. Defaults to 1.0.
  979. window : str or tuple or array_like, optional
  980. Desired window to use. If `window` is a string or tuple, it is
  981. passed to `get_window` to generate the window values, which are
  982. DFT-even by default. See `get_window` for a list of windows and
  983. required parameters. If `window` is array_like it will be used
  984. directly as the window and its length must be nperseg. Defaults
  985. to a Hann window. Must match the window used to generate the
  986. STFT for faithful inversion.
  987. nperseg : int, optional
  988. Number of data points corresponding to each STFT segment. This
  989. parameter must be specified if the number of data points per
  990. segment is odd, or if the STFT was padded via ``nfft >
  991. nperseg``. If `None`, the value depends on the shape of
  992. `Zxx` and `input_onesided`. If `input_onesided` is True,
  993. ``nperseg=2*(Zxx.shape[freq_axis] - 1)``. Otherwise,
  994. ``nperseg=Zxx.shape[freq_axis]``. Defaults to `None`.
  995. noverlap : int, optional
  996. Number of points to overlap between segments. If `None`, half
  997. of the segment length. Defaults to `None`. When specified, the
  998. COLA constraint must be met (see Notes below), and should match
  999. the parameter used to generate the STFT. Defaults to `None`.
  1000. nfft : int, optional
  1001. Number of FFT points corresponding to each STFT segment. This
  1002. parameter must be specified if the STFT was padded via ``nfft >
  1003. nperseg``. If `None`, the default values are the same as for
  1004. `nperseg`, detailed above, with one exception: if
  1005. `input_onesided` is True and
  1006. ``nperseg==2*Zxx.shape[freq_axis] - 1``, `nfft` also takes on
  1007. that value. This case allows the proper inversion of an
  1008. odd-length unpadded STFT using ``nfft=None``. Defaults to
  1009. `None`.
  1010. input_onesided : bool, optional
  1011. If `True`, interpret the input array as one-sided FFTs, such
  1012. as is returned by `stft` with ``return_onesided=True`` and
  1013. `numpy.fft.rfft`. If `False`, interpret the input as a a
  1014. two-sided FFT. Defaults to `True`.
  1015. boundary : bool, optional
  1016. Specifies whether the input signal was extended at its
  1017. boundaries by supplying a non-`None` ``boundary`` argument to
  1018. `stft`. Defaults to `True`.
  1019. time_axis : int, optional
  1020. Where the time segments of the STFT is located; the default is
  1021. the last axis (i.e. ``axis=-1``).
  1022. freq_axis : int, optional
  1023. Where the frequency axis of the STFT is located; the default is
  1024. the penultimate axis (i.e. ``axis=-2``).
  1025. Returns
  1026. -------
  1027. t : ndarray
  1028. Array of output data times.
  1029. x : ndarray
  1030. iSTFT of `Zxx`.
  1031. See Also
  1032. --------
  1033. stft: Short Time Fourier Transform
  1034. check_COLA: Check whether the Constant OverLap Add (COLA) constraint
  1035. is met
  1036. check_NOLA: Check whether the Nonzero Overlap Add (NOLA) constraint is met
  1037. Notes
  1038. -----
  1039. In order to enable inversion of an STFT via the inverse STFT with
  1040. `istft`, the signal windowing must obey the constraint of "nonzero
  1041. overlap add" (NOLA):
  1042. .. math:: \sum_{t}w^{2}[n-tH] \ne 0
  1043. This ensures that the normalization factors that appear in the denominator
  1044. of the overlap-add reconstruction equation
  1045. .. math:: x[n]=\frac{\sum_{t}x_{t}[n]w[n-tH]}{\sum_{t}w^{2}[n-tH]}
  1046. are not zero. The NOLA constraint can be checked with the `check_NOLA`
  1047. function.
  1048. An STFT which has been modified (via masking or otherwise) is not
  1049. guaranteed to correspond to a exactly realizible signal. This
  1050. function implements the iSTFT via the least-squares estimation
  1051. algorithm detailed in [2]_, which produces a signal that minimizes
  1052. the mean squared error between the STFT of the returned signal and
  1053. the modified STFT.
  1054. .. versionadded:: 0.19.0
  1055. References
  1056. ----------
  1057. .. [1] Oppenheim, Alan V., Ronald W. Schafer, John R. Buck
  1058. "Discrete-Time Signal Processing", Prentice Hall, 1999.
  1059. .. [2] Daniel W. Griffin, Jae S. Lim "Signal Estimation from
  1060. Modified Short-Time Fourier Transform", IEEE 1984,
  1061. 10.1109/TASSP.1984.1164317
  1062. Examples
  1063. --------
  1064. >>> from scipy import signal
  1065. >>> import matplotlib.pyplot as plt
  1066. Generate a test signal, a 2 Vrms sine wave at 50Hz corrupted by
  1067. 0.001 V**2/Hz of white noise sampled at 1024 Hz.
  1068. >>> fs = 1024
  1069. >>> N = 10*fs
  1070. >>> nperseg = 512
  1071. >>> amp = 2 * np.sqrt(2)
  1072. >>> noise_power = 0.001 * fs / 2
  1073. >>> time = np.arange(N) / float(fs)
  1074. >>> carrier = amp * np.sin(2*np.pi*50*time)
  1075. >>> noise = np.random.normal(scale=np.sqrt(noise_power),
  1076. ... size=time.shape)
  1077. >>> x = carrier + noise
  1078. Compute the STFT, and plot its magnitude
  1079. >>> f, t, Zxx = signal.stft(x, fs=fs, nperseg=nperseg)
  1080. >>> plt.figure()
  1081. >>> plt.pcolormesh(t, f, np.abs(Zxx), vmin=0, vmax=amp)
  1082. >>> plt.ylim([f[1], f[-1]])
  1083. >>> plt.title('STFT Magnitude')
  1084. >>> plt.ylabel('Frequency [Hz]')
  1085. >>> plt.xlabel('Time [sec]')
  1086. >>> plt.yscale('log')
  1087. >>> plt.show()
  1088. Zero the components that are 10% or less of the carrier magnitude,
  1089. then convert back to a time series via inverse STFT
  1090. >>> Zxx = np.where(np.abs(Zxx) >= amp/10, Zxx, 0)
  1091. >>> _, xrec = signal.istft(Zxx, fs)
  1092. Compare the cleaned signal with the original and true carrier signals.
  1093. >>> plt.figure()
  1094. >>> plt.plot(time, x, time, xrec, time, carrier)
  1095. >>> plt.xlim([2, 2.1])
  1096. >>> plt.xlabel('Time [sec]')
  1097. >>> plt.ylabel('Signal')
  1098. >>> plt.legend(['Carrier + Noise', 'Filtered via STFT', 'True Carrier'])
  1099. >>> plt.show()
  1100. Note that the cleaned signal does not start as abruptly as the original,
  1101. since some of the coefficients of the transient were also removed:
  1102. >>> plt.figure()
  1103. >>> plt.plot(time, x, time, xrec, time, carrier)
  1104. >>> plt.xlim([0, 0.1])
  1105. >>> plt.xlabel('Time [sec]')
  1106. >>> plt.ylabel('Signal')
  1107. >>> plt.legend(['Carrier + Noise', 'Filtered via STFT', 'True Carrier'])
  1108. >>> plt.show()
  1109. """
  1110. # Make sure input is an ndarray of appropriate complex dtype
  1111. Zxx = np.asarray(Zxx) + 0j
  1112. freq_axis = int(freq_axis)
  1113. time_axis = int(time_axis)
  1114. if Zxx.ndim < 2:
  1115. raise ValueError('Input stft must be at least 2d!')
  1116. if freq_axis == time_axis:
  1117. raise ValueError('Must specify differing time and frequency axes!')
  1118. nseg = Zxx.shape[time_axis]
  1119. if input_onesided:
  1120. # Assume even segment length
  1121. n_default = 2*(Zxx.shape[freq_axis] - 1)
  1122. else:
  1123. n_default = Zxx.shape[freq_axis]
  1124. # Check windowing parameters
  1125. if nperseg is None:
  1126. nperseg = n_default
  1127. else:
  1128. nperseg = int(nperseg)
  1129. if nperseg < 1:
  1130. raise ValueError('nperseg must be a positive integer')
  1131. if nfft is None:
  1132. if (input_onesided) and (nperseg == n_default + 1):
  1133. # Odd nperseg, no FFT padding
  1134. nfft = nperseg
  1135. else:
  1136. nfft = n_default
  1137. elif nfft < nperseg:
  1138. raise ValueError('nfft must be greater than or equal to nperseg.')
  1139. else:
  1140. nfft = int(nfft)
  1141. if noverlap is None:
  1142. noverlap = nperseg//2
  1143. else:
  1144. noverlap = int(noverlap)
  1145. if noverlap >= nperseg:
  1146. raise ValueError('noverlap must be less than nperseg.')
  1147. nstep = nperseg - noverlap
  1148. # Rearrange axes if necessary
  1149. if time_axis != Zxx.ndim-1 or freq_axis != Zxx.ndim-2:
  1150. # Turn negative indices to positive for the call to transpose
  1151. if freq_axis < 0:
  1152. freq_axis = Zxx.ndim + freq_axis
  1153. if time_axis < 0:
  1154. time_axis = Zxx.ndim + time_axis
  1155. zouter = list(range(Zxx.ndim))
  1156. for ax in sorted([time_axis, freq_axis], reverse=True):
  1157. zouter.pop(ax)
  1158. Zxx = np.transpose(Zxx, zouter+[freq_axis, time_axis])
  1159. # Get window as array
  1160. if isinstance(window, string_types) or type(window) is tuple:
  1161. win = get_window(window, nperseg)
  1162. else:
  1163. win = np.asarray(window)
  1164. if len(win.shape) != 1:
  1165. raise ValueError('window must be 1-D')
  1166. if win.shape[0] != nperseg:
  1167. raise ValueError('window must have length of {0}'.format(nperseg))
  1168. if input_onesided:
  1169. ifunc = np.fft.irfft
  1170. else:
  1171. ifunc = fftpack.ifft
  1172. xsubs = ifunc(Zxx, axis=-2, n=nfft)[..., :nperseg, :]
  1173. # Initialize output and normalization arrays
  1174. outputlength = nperseg + (nseg-1)*nstep
  1175. x = np.zeros(list(Zxx.shape[:-2])+[outputlength], dtype=xsubs.dtype)
  1176. norm = np.zeros(outputlength, dtype=xsubs.dtype)
  1177. if np.result_type(win, xsubs) != xsubs.dtype:
  1178. win = win.astype(xsubs.dtype)
  1179. xsubs *= win.sum() # This takes care of the 'spectrum' scaling
  1180. # Construct the output from the ifft segments
  1181. # This loop could perhaps be vectorized/strided somehow...
  1182. for ii in range(nseg):
  1183. # Window the ifft
  1184. x[..., ii*nstep:ii*nstep+nperseg] += xsubs[..., ii] * win
  1185. norm[..., ii*nstep:ii*nstep+nperseg] += win**2
  1186. # Remove extension points
  1187. if boundary:
  1188. x = x[..., nperseg//2:-(nperseg//2)]
  1189. norm = norm[..., nperseg//2:-(nperseg//2)]
  1190. # Divide out normalization where non-tiny
  1191. if np.sum(norm > 1e-10) != len(norm):
  1192. warnings.warn("NOLA condition failed, STFT may not be invertible")
  1193. x /= np.where(norm > 1e-10, norm, 1.0)
  1194. if input_onesided:
  1195. x = x.real
  1196. # Put axes back
  1197. if x.ndim > 1:
  1198. if time_axis != Zxx.ndim-1:
  1199. if freq_axis < time_axis:
  1200. time_axis -= 1
  1201. x = np.rollaxis(x, -1, time_axis)
  1202. time = np.arange(x.shape[0])/float(fs)
  1203. return time, x
  1204. def coherence(x, y, fs=1.0, window='hann', nperseg=None, noverlap=None,
  1205. nfft=None, detrend='constant', axis=-1):
  1206. r"""
  1207. Estimate the magnitude squared coherence estimate, Cxy, of
  1208. discrete-time signals X and Y using Welch's method.
  1209. ``Cxy = abs(Pxy)**2/(Pxx*Pyy)``, where `Pxx` and `Pyy` are power
  1210. spectral density estimates of X and Y, and `Pxy` is the cross
  1211. spectral density estimate of X and Y.
  1212. Parameters
  1213. ----------
  1214. x : array_like
  1215. Time series of measurement values
  1216. y : array_like
  1217. Time series of measurement values
  1218. fs : float, optional
  1219. Sampling frequency of the `x` and `y` time series. Defaults
  1220. to 1.0.
  1221. window : str or tuple or array_like, optional
  1222. Desired window to use. If `window` is a string or tuple, it is
  1223. passed to `get_window` to generate the window values, which are
  1224. DFT-even by default. See `get_window` for a list of windows and
  1225. required parameters. If `window` is array_like it will be used
  1226. directly as the window and its length must be nperseg. Defaults
  1227. to a Hann window.
  1228. nperseg : int, optional
  1229. Length of each segment. Defaults to None, but if window is str or
  1230. tuple, is set to 256, and if window is array_like, is set to the
  1231. length of the window.
  1232. noverlap: int, optional
  1233. Number of points to overlap between segments. If `None`,
  1234. ``noverlap = nperseg // 2``. Defaults to `None`.
  1235. nfft : int, optional
  1236. Length of the FFT used, if a zero padded FFT is desired. If
  1237. `None`, the FFT length is `nperseg`. Defaults to `None`.
  1238. detrend : str or function or `False`, optional
  1239. Specifies how to detrend each segment. If `detrend` is a
  1240. string, it is passed as the `type` argument to the `detrend`
  1241. function. If it is a function, it takes a segment and returns a
  1242. detrended segment. If `detrend` is `False`, no detrending is
  1243. done. Defaults to 'constant'.
  1244. axis : int, optional
  1245. Axis along which the coherence is computed for both inputs; the
  1246. default is over the last axis (i.e. ``axis=-1``).
  1247. Returns
  1248. -------
  1249. f : ndarray
  1250. Array of sample frequencies.
  1251. Cxy : ndarray
  1252. Magnitude squared coherence of x and y.
  1253. See Also
  1254. --------
  1255. periodogram: Simple, optionally modified periodogram
  1256. lombscargle: Lomb-Scargle periodogram for unevenly sampled data
  1257. welch: Power spectral density by Welch's method.
  1258. csd: Cross spectral density by Welch's method.
  1259. Notes
  1260. --------
  1261. An appropriate amount of overlap will depend on the choice of window
  1262. and on your requirements. For the default Hann window an overlap of
  1263. 50% is a reasonable trade off between accurately estimating the
  1264. signal power, while not over counting any of the data. Narrower
  1265. windows may require a larger overlap.
  1266. .. versionadded:: 0.16.0
  1267. References
  1268. ----------
  1269. .. [1] P. Welch, "The use of the fast Fourier transform for the
  1270. estimation of power spectra: A method based on time averaging
  1271. over short, modified periodograms", IEEE Trans. Audio
  1272. Electroacoust. vol. 15, pp. 70-73, 1967.
  1273. .. [2] Stoica, Petre, and Randolph Moses, "Spectral Analysis of
  1274. Signals" Prentice Hall, 2005
  1275. Examples
  1276. --------
  1277. >>> from scipy import signal
  1278. >>> import matplotlib.pyplot as plt
  1279. Generate two test signals with some common features.
  1280. >>> fs = 10e3
  1281. >>> N = 1e5
  1282. >>> amp = 20
  1283. >>> freq = 1234.0
  1284. >>> noise_power = 0.001 * fs / 2
  1285. >>> time = np.arange(N) / fs
  1286. >>> b, a = signal.butter(2, 0.25, 'low')
  1287. >>> x = np.random.normal(scale=np.sqrt(noise_power), size=time.shape)
  1288. >>> y = signal.lfilter(b, a, x)
  1289. >>> x += amp*np.sin(2*np.pi*freq*time)
  1290. >>> y += np.random.normal(scale=0.1*np.sqrt(noise_power), size=time.shape)
  1291. Compute and plot the coherence.
  1292. >>> f, Cxy = signal.coherence(x, y, fs, nperseg=1024)
  1293. >>> plt.semilogy(f, Cxy)
  1294. >>> plt.xlabel('frequency [Hz]')
  1295. >>> plt.ylabel('Coherence')
  1296. >>> plt.show()
  1297. """
  1298. freqs, Pxx = welch(x, fs=fs, window=window, nperseg=nperseg,
  1299. noverlap=noverlap, nfft=nfft, detrend=detrend,
  1300. axis=axis)
  1301. _, Pyy = welch(y, fs=fs, window=window, nperseg=nperseg, noverlap=noverlap,
  1302. nfft=nfft, detrend=detrend, axis=axis)
  1303. _, Pxy = csd(x, y, fs=fs, window=window, nperseg=nperseg,
  1304. noverlap=noverlap, nfft=nfft, detrend=detrend, axis=axis)
  1305. Cxy = np.abs(Pxy)**2 / Pxx / Pyy
  1306. return freqs, Cxy
  1307. def _spectral_helper(x, y, fs=1.0, window='hann', nperseg=None, noverlap=None,
  1308. nfft=None, detrend='constant', return_onesided=True,
  1309. scaling='spectrum', axis=-1, mode='psd', boundary=None,
  1310. padded=False):
  1311. """
  1312. Calculate various forms of windowed FFTs for PSD, CSD, etc.
  1313. This is a helper function that implements the commonality between
  1314. the stft, psd, csd, and spectrogram functions. It is not designed to
  1315. be called externally. The windows are not averaged over; the result
  1316. from each window is returned.
  1317. Parameters
  1318. ---------
  1319. x : array_like
  1320. Array or sequence containing the data to be analyzed.
  1321. y : array_like
  1322. Array or sequence containing the data to be analyzed. If this is
  1323. the same object in memory as `x` (i.e. ``_spectral_helper(x,
  1324. x, ...)``), the extra computations are spared.
  1325. fs : float, optional
  1326. Sampling frequency of the time series. Defaults to 1.0.
  1327. window : str or tuple or array_like, optional
  1328. Desired window to use. If `window` is a string or tuple, it is
  1329. passed to `get_window` to generate the window values, which are
  1330. DFT-even by default. See `get_window` for a list of windows and
  1331. required parameters. If `window` is array_like it will be used
  1332. directly as the window and its length must be nperseg. Defaults
  1333. to a Hann window.
  1334. nperseg : int, optional
  1335. Length of each segment. Defaults to None, but if window is str or
  1336. tuple, is set to 256, and if window is array_like, is set to the
  1337. length of the window.
  1338. noverlap : int, optional
  1339. Number of points to overlap between segments. If `None`,
  1340. ``noverlap = nperseg // 2``. Defaults to `None`.
  1341. nfft : int, optional
  1342. Length of the FFT used, if a zero padded FFT is desired. If
  1343. `None`, the FFT length is `nperseg`. Defaults to `None`.
  1344. detrend : str or function or `False`, optional
  1345. Specifies how to detrend each segment. If `detrend` is a
  1346. string, it is passed as the `type` argument to the `detrend`
  1347. function. If it is a function, it takes a segment and returns a
  1348. detrended segment. If `detrend` is `False`, no detrending is
  1349. done. Defaults to 'constant'.
  1350. return_onesided : bool, optional
  1351. If `True`, return a one-sided spectrum for real data. If
  1352. `False` return a two-sided spectrum. Note that for complex
  1353. data, a two-sided spectrum is always returned.
  1354. scaling : { 'density', 'spectrum' }, optional
  1355. Selects between computing the cross spectral density ('density')
  1356. where `Pxy` has units of V**2/Hz and computing the cross
  1357. spectrum ('spectrum') where `Pxy` has units of V**2, if `x`
  1358. and `y` are measured in V and `fs` is measured in Hz.
  1359. Defaults to 'density'
  1360. axis : int, optional
  1361. Axis along which the FFTs are computed; the default is over the
  1362. last axis (i.e. ``axis=-1``).
  1363. mode: str {'psd', 'stft'}, optional
  1364. Defines what kind of return values are expected. Defaults to
  1365. 'psd'.
  1366. boundary : str or None, optional
  1367. Specifies whether the input signal is extended at both ends, and
  1368. how to generate the new values, in order to center the first
  1369. windowed segment on the first input point. This has the benefit
  1370. of enabling reconstruction of the first input point when the
  1371. employed window function starts at zero. Valid options are
  1372. ``['even', 'odd', 'constant', 'zeros', None]``. Defaults to
  1373. `None`.
  1374. padded : bool, optional
  1375. Specifies whether the input signal is zero-padded at the end to
  1376. make the signal fit exactly into an integer number of window
  1377. segments, so that all of the signal is included in the output.
  1378. Defaults to `False`. Padding occurs after boundary extension, if
  1379. `boundary` is not `None`, and `padded` is `True`.
  1380. Returns
  1381. -------
  1382. freqs : ndarray
  1383. Array of sample frequencies.
  1384. t : ndarray
  1385. Array of times corresponding to each data segment
  1386. result : ndarray
  1387. Array of output data, contents dependent on *mode* kwarg.
  1388. Notes
  1389. -----
  1390. Adapted from matplotlib.mlab
  1391. .. versionadded:: 0.16.0
  1392. """
  1393. if mode not in ['psd', 'stft']:
  1394. raise ValueError("Unknown value for mode %s, must be one of: "
  1395. "{'psd', 'stft'}" % mode)
  1396. boundary_funcs = {'even': even_ext,
  1397. 'odd': odd_ext,
  1398. 'constant': const_ext,
  1399. 'zeros': zero_ext,
  1400. None: None}
  1401. if boundary not in boundary_funcs:
  1402. raise ValueError("Unknown boundary option '{0}', must be one of: {1}"
  1403. .format(boundary, list(boundary_funcs.keys())))
  1404. # If x and y are the same object we can save ourselves some computation.
  1405. same_data = y is x
  1406. if not same_data and mode != 'psd':
  1407. raise ValueError("x and y must be equal if mode is 'stft'")
  1408. axis = int(axis)
  1409. # Ensure we have np.arrays, get outdtype
  1410. x = np.asarray(x)
  1411. if not same_data:
  1412. y = np.asarray(y)
  1413. outdtype = np.result_type(x, y, np.complex64)
  1414. else:
  1415. outdtype = np.result_type(x, np.complex64)
  1416. if not same_data:
  1417. # Check if we can broadcast the outer axes together
  1418. xouter = list(x.shape)
  1419. youter = list(y.shape)
  1420. xouter.pop(axis)
  1421. youter.pop(axis)
  1422. try:
  1423. outershape = np.broadcast(np.empty(xouter), np.empty(youter)).shape
  1424. except ValueError:
  1425. raise ValueError('x and y cannot be broadcast together.')
  1426. if same_data:
  1427. if x.size == 0:
  1428. return np.empty(x.shape), np.empty(x.shape), np.empty(x.shape)
  1429. else:
  1430. if x.size == 0 or y.size == 0:
  1431. outshape = outershape + (min([x.shape[axis], y.shape[axis]]),)
  1432. emptyout = np.rollaxis(np.empty(outshape), -1, axis)
  1433. return emptyout, emptyout, emptyout
  1434. if x.ndim > 1:
  1435. if axis != -1:
  1436. x = np.rollaxis(x, axis, len(x.shape))
  1437. if not same_data and y.ndim > 1:
  1438. y = np.rollaxis(y, axis, len(y.shape))
  1439. # Check if x and y are the same length, zero-pad if necessary
  1440. if not same_data:
  1441. if x.shape[-1] != y.shape[-1]:
  1442. if x.shape[-1] < y.shape[-1]:
  1443. pad_shape = list(x.shape)
  1444. pad_shape[-1] = y.shape[-1] - x.shape[-1]
  1445. x = np.concatenate((x, np.zeros(pad_shape)), -1)
  1446. else:
  1447. pad_shape = list(y.shape)
  1448. pad_shape[-1] = x.shape[-1] - y.shape[-1]
  1449. y = np.concatenate((y, np.zeros(pad_shape)), -1)
  1450. if nperseg is not None: # if specified by user
  1451. nperseg = int(nperseg)
  1452. if nperseg < 1:
  1453. raise ValueError('nperseg must be a positive integer')
  1454. # parse window; if array like, then set nperseg = win.shape
  1455. win, nperseg = _triage_segments(window, nperseg, input_length=x.shape[-1])
  1456. if nfft is None:
  1457. nfft = nperseg
  1458. elif nfft < nperseg:
  1459. raise ValueError('nfft must be greater than or equal to nperseg.')
  1460. else:
  1461. nfft = int(nfft)
  1462. if noverlap is None:
  1463. noverlap = nperseg//2
  1464. else:
  1465. noverlap = int(noverlap)
  1466. if noverlap >= nperseg:
  1467. raise ValueError('noverlap must be less than nperseg.')
  1468. nstep = nperseg - noverlap
  1469. # Padding occurs after boundary extension, so that the extended signal ends
  1470. # in zeros, instead of introducing an impulse at the end.
  1471. # I.e. if x = [..., 3, 2]
  1472. # extend then pad -> [..., 3, 2, 2, 3, 0, 0, 0]
  1473. # pad then extend -> [..., 3, 2, 0, 0, 0, 2, 3]
  1474. if boundary is not None:
  1475. ext_func = boundary_funcs[boundary]
  1476. x = ext_func(x, nperseg//2, axis=-1)
  1477. if not same_data:
  1478. y = ext_func(y, nperseg//2, axis=-1)
  1479. if padded:
  1480. # Pad to integer number of windowed segments
  1481. # I.e make x.shape[-1] = nperseg + (nseg-1)*nstep, with integer nseg
  1482. nadd = (-(x.shape[-1]-nperseg) % nstep) % nperseg
  1483. zeros_shape = list(x.shape[:-1]) + [nadd]
  1484. x = np.concatenate((x, np.zeros(zeros_shape)), axis=-1)
  1485. if not same_data:
  1486. zeros_shape = list(y.shape[:-1]) + [nadd]
  1487. y = np.concatenate((y, np.zeros(zeros_shape)), axis=-1)
  1488. # Handle detrending and window functions
  1489. if not detrend:
  1490. def detrend_func(d):
  1491. return d
  1492. elif not hasattr(detrend, '__call__'):
  1493. def detrend_func(d):
  1494. return signaltools.detrend(d, type=detrend, axis=-1)
  1495. elif axis != -1:
  1496. # Wrap this function so that it receives a shape that it could
  1497. # reasonably expect to receive.
  1498. def detrend_func(d):
  1499. d = np.rollaxis(d, -1, axis)
  1500. d = detrend(d)
  1501. return np.rollaxis(d, axis, len(d.shape))
  1502. else:
  1503. detrend_func = detrend
  1504. if np.result_type(win, np.complex64) != outdtype:
  1505. win = win.astype(outdtype)
  1506. if scaling == 'density':
  1507. scale = 1.0 / (fs * (win*win).sum())
  1508. elif scaling == 'spectrum':
  1509. scale = 1.0 / win.sum()**2
  1510. else:
  1511. raise ValueError('Unknown scaling: %r' % scaling)
  1512. if mode == 'stft':
  1513. scale = np.sqrt(scale)
  1514. if return_onesided:
  1515. if np.iscomplexobj(x):
  1516. sides = 'twosided'
  1517. warnings.warn('Input data is complex, switching to '
  1518. 'return_onesided=False')
  1519. else:
  1520. sides = 'onesided'
  1521. if not same_data:
  1522. if np.iscomplexobj(y):
  1523. sides = 'twosided'
  1524. warnings.warn('Input data is complex, switching to '
  1525. 'return_onesided=False')
  1526. else:
  1527. sides = 'twosided'
  1528. if sides == 'twosided':
  1529. freqs = fftpack.fftfreq(nfft, 1/fs)
  1530. elif sides == 'onesided':
  1531. freqs = np.fft.rfftfreq(nfft, 1/fs)
  1532. # Perform the windowed FFTs
  1533. result = _fft_helper(x, win, detrend_func, nperseg, noverlap, nfft, sides)
  1534. if not same_data:
  1535. # All the same operations on the y data
  1536. result_y = _fft_helper(y, win, detrend_func, nperseg, noverlap, nfft,
  1537. sides)
  1538. result = np.conjugate(result) * result_y
  1539. elif mode == 'psd':
  1540. result = np.conjugate(result) * result
  1541. result *= scale
  1542. if sides == 'onesided' and mode == 'psd':
  1543. if nfft % 2:
  1544. result[..., 1:] *= 2
  1545. else:
  1546. # Last point is unpaired Nyquist freq point, don't double
  1547. result[..., 1:-1] *= 2
  1548. time = np.arange(nperseg/2, x.shape[-1] - nperseg/2 + 1,
  1549. nperseg - noverlap)/float(fs)
  1550. if boundary is not None:
  1551. time -= (nperseg/2) / fs
  1552. result = result.astype(outdtype)
  1553. # All imaginary parts are zero anyways
  1554. if same_data and mode != 'stft':
  1555. result = result.real
  1556. # Output is going to have new last axis for time/window index, so a
  1557. # negative axis index shifts down one
  1558. if axis < 0:
  1559. axis -= 1
  1560. # Roll frequency axis back to axis where the data came from
  1561. result = np.rollaxis(result, -1, axis)
  1562. return freqs, time, result
  1563. def _fft_helper(x, win, detrend_func, nperseg, noverlap, nfft, sides):
  1564. """
  1565. Calculate windowed FFT, for internal use by
  1566. scipy.signal._spectral_helper
  1567. This is a helper function that does the main FFT calculation for
  1568. `_spectral helper`. All input validation is performed there, and the
  1569. data axis is assumed to be the last axis of x. It is not designed to
  1570. be called externally. The windows are not averaged over; the result
  1571. from each window is returned.
  1572. Returns
  1573. -------
  1574. result : ndarray
  1575. Array of FFT data
  1576. Notes
  1577. -----
  1578. Adapted from matplotlib.mlab
  1579. .. versionadded:: 0.16.0
  1580. """
  1581. # Created strided array of data segments
  1582. if nperseg == 1 and noverlap == 0:
  1583. result = x[..., np.newaxis]
  1584. else:
  1585. # https://stackoverflow.com/a/5568169
  1586. step = nperseg - noverlap
  1587. shape = x.shape[:-1]+((x.shape[-1]-noverlap)//step, nperseg)
  1588. strides = x.strides[:-1]+(step*x.strides[-1], x.strides[-1])
  1589. result = np.lib.stride_tricks.as_strided(x, shape=shape,
  1590. strides=strides)
  1591. # Detrend each data segment individually
  1592. result = detrend_func(result)
  1593. # Apply window by multiplication
  1594. result = win * result
  1595. # Perform the fft. Acts on last axis by default. Zero-pads automatically
  1596. if sides == 'twosided':
  1597. func = fftpack.fft
  1598. else:
  1599. result = result.real
  1600. func = np.fft.rfft
  1601. result = func(result, n=nfft)
  1602. return result
  1603. def _triage_segments(window, nperseg, input_length):
  1604. """
  1605. Parses window and nperseg arguments for spectrogram and _spectral_helper.
  1606. This is a helper function, not meant to be called externally.
  1607. Parameters
  1608. ----------
  1609. window : string, tuple, or ndarray
  1610. If window is specified by a string or tuple and nperseg is not
  1611. specified, nperseg is set to the default of 256 and returns a window of
  1612. that length.
  1613. If instead the window is array_like and nperseg is not specified, then
  1614. nperseg is set to the length of the window. A ValueError is raised if
  1615. the user supplies both an array_like window and a value for nperseg but
  1616. nperseg does not equal the length of the window.
  1617. nperseg : int
  1618. Length of each segment
  1619. input_length: int
  1620. Length of input signal, i.e. x.shape[-1]. Used to test for errors.
  1621. Returns
  1622. -------
  1623. win : ndarray
  1624. window. If function was called with string or tuple than this will hold
  1625. the actual array used as a window.
  1626. nperseg : int
  1627. Length of each segment. If window is str or tuple, nperseg is set to
  1628. 256. If window is array_like, nperseg is set to the length of the
  1629. 6
  1630. window.
  1631. """
  1632. # parse window; if array like, then set nperseg = win.shape
  1633. if isinstance(window, string_types) or isinstance(window, tuple):
  1634. # if nperseg not specified
  1635. if nperseg is None:
  1636. nperseg = 256 # then change to default
  1637. if nperseg > input_length:
  1638. warnings.warn('nperseg = {0:d} is greater than input length '
  1639. ' = {1:d}, using nperseg = {1:d}'
  1640. .format(nperseg, input_length))
  1641. nperseg = input_length
  1642. win = get_window(window, nperseg)
  1643. else:
  1644. win = np.asarray(window)
  1645. if len(win.shape) != 1:
  1646. raise ValueError('window must be 1-D')
  1647. if input_length < win.shape[-1]:
  1648. raise ValueError('window is longer than input signal')
  1649. if nperseg is None:
  1650. nperseg = win.shape[0]
  1651. elif nperseg is not None:
  1652. if nperseg != win.shape[0]:
  1653. raise ValueError("value specified for nperseg is different"
  1654. " from length of window")
  1655. return win, nperseg
  1656. def _median_bias(n):
  1657. """
  1658. Returns the bias of the median of a set of periodograms relative to
  1659. the mean.
  1660. See arXiv:gr-qc/0509116 Appendix B for details.
  1661. Parameters
  1662. ----------
  1663. n : int
  1664. Numbers of periodograms being averaged.
  1665. Returns
  1666. -------
  1667. bias : float
  1668. Calculated bias.
  1669. """
  1670. ii_2 = 2 * np.arange(1., (n-1) // 2 + 1)
  1671. return 1 + np.sum(1. / (ii_2 + 1) - 1. / ii_2)