L iCwdZddlZddlmZddlmZddl m Z ddl m Z m Z ddlmZddlmZmZmZmZddlZdd lmZmZgd Z d$dd d d ej2dej2dej2dedeedzdej6dzdedej6fdZ d%dZ d&dZ d&dZ d'dZ d(dZ!d(dZ" d)dZ# d*dZ$ d+dZ% d,d Z&d!Z'd"Z(d#Z)y)-zTools for spectral analysis. N)fft) _signaltools) ShortTimeFFT FFT_MODE_TYPE) get_window) const_exteven_extodd_extzero_ext)castLiteral) periodogramwelch lombscarglecsd coherence spectrogramstftistft check_COLA check_NOLAF)weights floating_meanxyfreqs precenter normalizepowerr amplituderrreturnc |&tj|tj}n%tj|tj}tj|tj}tj|tj}tj|tj}|jdk(rD|j dkDr5|j |j cxk(r|j k(stdtd|jdk(r|j dkDs tdtj|dk\rtj|dkDs tdt|tr|rdnd}|d vr td |d |jz z}|r||jz }|jdd }|jd d}|jd d}|jd d}||z}||z}tj|} tj|} tj |j"|} tj |j"| | z} d | z } tj |j"| | z}|rXtj |j"| }tj |j"| }| ||zz} | ||zz} |||zz}d tj$d|z| | z z}||z }tj|}tj|}tj |j"|}tj |j"|}tj |j"||z} d | z } |r`tj |j"|}tj |j"|}|| |zz}|| |zz}| ||zz} | ||zz} tj&|j(j*}|| | |k<|| | |k<|| z }|| z }d||z||zzz}tj,|}|dk(r |t/|j ddz z}|S|dk(rGtj |j"|}|r|| | zz}|d tj,|z z}|Stj,|}tj,|}tj,|}|d|zztj0d|zz}|S)u Compute the generalized Lomb-Scargle periodogram. The Lomb-Scargle periodogram was developed by Lomb [1]_ and further extended by Scargle [2]_ to find, and test the significance of weak periodic signals with uneven temporal sampling. The algorithm used here is based on a weighted least-squares fit of the form ``y(ω) = a*cos(ω*x) + b*sin(ω*x) + c``, where the fit is calculated for each frequency independently. This algorithm was developed by Zechmeister and Kürster which improves the Lomb-Scargle periodogram by enabling the weighting of individual samples and calculating an unknown y offset (also called a "floating-mean" model) [3]_. For more details, and practical considerations, see the excellent reference on the Lomb-Scargle periodogram [4]_. When *normalize* is False (or "power") (default) the computed periodogram is unnormalized, it takes the value ``(A**2) * N/4`` for a harmonic signal with amplitude A for sufficiently large N. Where N is the length of x or y. When *normalize* is True (or "normalize") the computed periodogram is normalized by the residuals of the data around a constant reference model (at zero). When *normalize* is "amplitude" the computed periodogram is the complex representation of the amplitude and phase. Input arrays should be 1-D of a real floating data type, which are converted into float64 arrays before processing. Parameters ---------- x : array_like Sample times. y : array_like Measurement values. Values are assumed to have a baseline of ``y = 0``. If there is a possibility of a y offset, it is recommended to set `floating_mean` to True. freqs : array_like Angular frequencies (e.g., having unit rad/s=2π/s for `x` having unit s) for output periodogram. Frequencies are normally >= 0, as any peak at ``-freq`` will also exist at ``+freq``. precenter : bool, optional Pre-center measurement values by subtracting the mean, if True. This is a legacy parameter and unnecessary if `floating_mean` is True. normalize : bool | str, optional Compute normalized or complex (amplitude + phase) periodogram. Valid options are: ``False``/``"power"``, ``True``/``"normalize"``, or ``"amplitude"``. weights : array_like, optional Weights for each sample. Weights must be nonnegative. floating_mean : bool, optional Determines a y offset for each frequency independently, if True. Else the y offset is assumed to be `0`. Returns ------- pgram : array_like Lomb-Scargle periodogram. Raises ------ ValueError If any of the input arrays x, y, freqs, or weights are not 1D, or if any are zero length. Or, if the input arrays x, y, and weights do not have the same shape as each other. ValueError If any weight is < 0, or the sum of the weights is <= 0. ValueError If the normalize parameter is not one of the allowed options. See Also -------- periodogram: Power spectral density using a periodogram welch: Power spectral density by Welch's method csd: Cross spectral density by Welch's method Notes ----- The algorithm used will not automatically account for any unknown y offset, unless floating_mean is True. Therefore, for most use cases, if there is a possibility of a y offset, it is recommended to set floating_mean to True. If precenter is True, it performs the operation ``y -= y.mean()``. However, precenter is a legacy parameter, and unnecessary when floating_mean is True. Furthermore, the mean removed by precenter does not account for sample weights, nor will it correct for any bias due to consistently missing observations at peaks and/or troughs. When the normalize parameter is "amplitude", for any frequency in freqs that is below ``(2*pi)/(x.max() - x.min())``, the predicted amplitude will tend towards infinity. The concept of a "Nyquist frequency" limit (see Nyquist-Shannon sampling theorem) is not generally applicable to unevenly sampled data. Therefore, with unevenly sampled data, valid frequencies in freqs can often be much higher than expected. References ---------- .. [1] N.R. Lomb "Least-squares frequency analysis of unequally spaced data", Astrophysics and Space Science, vol 39, pp. 447-462, 1976 :doi:`10.1007/bf00648343` .. [2] J.D. Scargle "Studies in astronomical time series analysis. II - Statistical aspects of spectral analysis of unevenly spaced data", The Astrophysical Journal, vol 263, pp. 835-853, 1982 :doi:`10.1086/160554` .. [3] M. Zechmeister and M. Kürster, "The generalised Lomb-Scargle periodogram. A new formalism for the floating-mean and Keplerian periodograms," Astronomy and Astrophysics, vol. 496, pp. 577-584, 2009 :doi:`10.1051/0004-6361:200811296` .. [4] J.T. VanderPlas, "Understanding the Lomb-Scargle Periodogram," The Astrophysical Journal Supplement Series, vol. 236, no. 1, p. 16, May 2018 :doi:`10.3847/1538-4365/aab766` Examples -------- >>> import numpy as np >>> rng = np.random.default_rng() First define some input parameters for the signal: >>> A = 2. # amplitude >>> c = 2. # offset >>> w0 = 1. # rad/sec >>> nin = 150 >>> nout = 1002 Randomly generate sample times: >>> x = rng.uniform(0, 10*np.pi, nin) Plot a sine wave for the selected times: >>> y = A * np.cos(w0*x) + c Define the array of frequencies for which to compute the periodogram: >>> w = np.linspace(0.25, 10, nout) Calculate Lomb-Scargle periodogram for each of the normalize options: >>> from scipy.signal import lombscargle >>> pgram_power = lombscargle(x, y, w, normalize=False) >>> pgram_norm = lombscargle(x, y, w, normalize=True) >>> pgram_amp = lombscargle(x, y, w, normalize='amplitude') ... >>> pgram_power_f = lombscargle(x, y, w, normalize=False, floating_mean=True) >>> pgram_norm_f = lombscargle(x, y, w, normalize=True, floating_mean=True) >>> pgram_amp_f = lombscargle(x, y, w, normalize='amplitude', floating_mean=True) Now make a plot of the input data: >>> import matplotlib.pyplot as plt >>> fig, (ax_t, ax_p, ax_n, ax_a) = plt.subplots(4, 1, figsize=(5, 6)) >>> ax_t.plot(x, y, 'b+') >>> ax_t.set_xlabel('Time [s]') >>> ax_t.set_ylabel('Amplitude') Then plot the periodogram for each of the normalize options, as well as with and without floating_mean=True: >>> ax_p.plot(w, pgram_power, label='default') >>> ax_p.plot(w, pgram_power_f, label='floating_mean=True') >>> ax_p.set_xlabel('Angular frequency [rad/s]') >>> ax_p.set_ylabel('Power') >>> ax_p.legend(prop={'size': 7}) ... >>> ax_n.plot(w, pgram_norm, label='default') >>> ax_n.plot(w, pgram_norm_f, label='floating_mean=True') >>> ax_n.set_xlabel('Angular frequency [rad/s]') >>> ax_n.set_ylabel('Normalized') >>> ax_n.legend(prop={'size': 7}) ... >>> ax_a.plot(w, np.abs(pgram_amp), label='default') >>> ax_a.plot(w, np.abs(pgram_amp_f), label='floating_mean=True') >>> ax_a.set_xlabel('Angular frequency [rad/s]') >>> ax_a.set_ylabel('Amplitude') >>> ax_a.legend(prop={'size': 7}) ... >>> plt.tight_layout() >>> plt.show() dtyperrzEParameters x, y, weights must be 1-D arrays of equal non-zero length!z7Parameter freqs must be a 1-D array of non-zero length!zTParameter weights must have only non-negative entries which sum to a positive value!rr!r zMNormalize must be: False (or 'power'), True (or 'normalize'), or 'amplitude'.?g?g@g@?)np ones_likefloat64asarrayndimsizeshape ValueErrorallsum isinstanceboolmeanreshapecossindotTarctan2finfor&epsnegsqueezefloatexp)rrrrrrr weights_yfreqstcoswtsinwtYCCSSCSCStau freqst_tau coswt_tau sinwt_tauYCYSr>abpgramYYs _/mnt/ssd/data/python-lab/Trading/venv/lib/python3.12/site-packages/scipy/signal/_spectral_py.pyrrs@,,q 3**WBJJ7 1BJJ'A 1BJJ'A JJuBJJ /E FFaKAFFQJ177agg+N+N23 3,O23 3 JJ!O QRSS FF7a< RVVG_q%8:; ;)T"#,K' ;;   w{{},-G L MM!R E "aA "aAoob!$G! I QYF FF6NE FF6NE wyy!A  55= )B rB  55= )B FF799e $ FF799e $ a!e  a!e  a!e   38R"W- -C#Jz"Iz"I  Y 'B  Y 'B  9y0 1B rB FF799i ( FF799i ( a!e  a!e  a!e  a!e  XXAGG $ + +FBrF{OBrF{O RA RA 1r6AF? #E JJu EG qwwqz"S((. L+ k !VVIKK #  !a%KB rzz"~%% L JJqM JJqMjjoR!Vrvvb3h// Lc tj|}|jdk(r>tj|jtj|jfS|d}||j|}n||j|k(r|}n||j|kDr|j|}nd||j|krRtj ddgt |jz} tj d|| |<|t| }|}d}t|dr|jk7r tdt|||d||||| S)u Estimate power spectral density using a periodogram. Parameters ---------- x : array_like Time series of measurement values fs : float, optional Sampling frequency of the `x` time series. Defaults to 1.0. window : str or tuple or array_like, optional Desired window to use. If `window` is a string or tuple, it is passed to `get_window` to generate the window values, which are DFT-even by default. See `get_window` for a list of windows and required parameters. If `window` is array_like it will be used directly as the window and its length must be equal to the length of the axis over which the periodogram is computed. Defaults to 'boxcar'. nfft : int, optional Length of the FFT used. If `None` the length of `x` will be used. detrend : str or function or `False`, optional Specifies how to detrend each segment. If `detrend` is a string, it is passed as the `type` argument to the `detrend` function. If it is a function, it takes a segment and returns a detrended segment. If `detrend` is `False`, no detrending is done. Defaults to 'constant'. return_onesided : bool, optional If `True`, return a one-sided spectrum for real data. If `False` return a two-sided spectrum. Defaults to `True`, but for complex data, a two-sided spectrum is always returned. scaling : { 'density', 'spectrum' }, optional Selects between computing the power spectral density ('density') where `Pxx` has units of V²/Hz and computing the squared magnitude spectrum ('spectrum') where `Pxx` has units of V², if `x` is measured in V and `fs` is measured in Hz. Defaults to 'density' axis : int, optional Axis along which the periodogram is computed; the default is over the last axis (i.e. ``axis=-1``). Returns ------- f : ndarray Array of sample frequencies. Pxx : ndarray Power spectral density or power spectrum of `x`. See Also -------- welch: Estimate power spectral density using Welch's method lombscargle: Lomb-Scargle periodogram for unevenly sampled data Notes ----- The ratio of the squared magnitude (``scaling='spectrum'``) divided by the spectral power density (``scaling='density'``) is the constant factor of ``sum(abs(window)**2)*fs / abs(sum(window))**2``. If `return_onesided` is ``True``, the values of the negative frequencies are added to values of the corresponding positive ones. Consult the :ref:`tutorial_SpectralAnalysis` section of the :ref:`user_guide` for a discussion of the scalings of the power spectral density and the magnitude (squared) spectrum. .. versionadded:: 0.12.0 Examples -------- >>> import numpy as np >>> from scipy import signal >>> import matplotlib.pyplot as plt >>> rng = np.random.default_rng() Generate a test signal, a 2 Vrms sine wave at 1234 Hz, corrupted by 0.001 V**2/Hz of white noise sampled at 10 kHz. >>> fs = 10e3 >>> N = 1e5 >>> amp = 2*np.sqrt(2) >>> freq = 1234.0 >>> noise_power = 0.001 * fs / 2 >>> time = np.arange(N) / fs >>> x = amp*np.sin(2*np.pi*freq*time) >>> x += rng.normal(scale=np.sqrt(noise_power), size=time.shape) Compute and plot the power spectral density. >>> f, Pxx_den = signal.periodogram(x, fs) >>> plt.semilogy(f, Pxx_den) >>> plt.ylim([1e-7, 1e2]) >>> plt.xlabel('frequency [Hz]') >>> plt.ylabel('PSD [V**2/Hz]') >>> plt.show() If we average the last half of the spectral density, to exclude the peak, we can recover the noise power on the signal. >>> np.mean(Pxx_den[25000:]) 0.000985320699252543 Now compute and plot the power spectrum. >>> f, Pxx_spec = signal.periodogram(x, fs, 'flattop', scaling='spectrum') >>> plt.figure() >>> plt.semilogy(f, np.sqrt(Pxx_spec)) >>> plt.ylim([1e-4, 1e1]) >>> plt.xlabel('frequency [Hz]') >>> plt.ylabel('Linear spectrum [V RMS]') >>> plt.show() The peak height in the power spectrum is an estimate of the RMS amplitude. >>> np.sqrt(Pxx_spec.max()) 2.0077340678640727 rNboxcarr/zOthe size of the window must be the same size of the input on the specified axis) fswindownpersegnoverlapnfftdetrendreturn_onesidedscalingaxis) r*r-r/emptyr0s_lentuplehasattrr1r) rrZr[r^r_r`rarbr\ss rVrr[s6n 1 Avv{xx "((177"333 ~ |''$-    ''$-   UU1XJs177| #%%,$ eAhKvv ;;' !BC C r&'AG_ t --rWc Rt||||||||||| |  \} } | | jfS)a Estimate power spectral density using Welch's method. Welch's method [1]_ computes an estimate of the power spectral density by dividing the data into overlapping segments, computing a modified periodogram for each segment and averaging the periodograms. Parameters ---------- x : array_like Time series of measurement values fs : float, optional Sampling frequency of the `x` time series. Defaults to 1.0. window : str or tuple or array_like, optional Desired window to use. If `window` is a string or tuple, it is passed to `get_window` to generate the window values, which are DFT-even by default. See `get_window` for a list of windows and required parameters. If `window` is array_like it will be used directly as the window and its length must be nperseg. Defaults to a Hann window. nperseg : int, optional Length of each segment. Defaults to None, but if window is str or tuple, is set to 256, and if window is array_like, is set to the length of the window. noverlap : int, optional Number of points to overlap between segments. If `None`, ``noverlap = nperseg // 2``. Defaults to `None`. nfft : int, optional Length of the FFT used, if a zero padded FFT is desired. If `None`, the FFT length is `nperseg`. Defaults to `None`. detrend : str or function or `False`, optional Specifies how to detrend each segment. If `detrend` is a string, it is passed as the `type` argument to the `detrend` function. If it is a function, it takes a segment and returns a detrended segment. If `detrend` is `False`, no detrending is done. Defaults to 'constant'. return_onesided : bool, optional If `True`, return a one-sided spectrum for real data. If `False` return a two-sided spectrum. Defaults to `True`, but for complex data, a two-sided spectrum is always returned. scaling : { 'density', 'spectrum' }, optional Selects between computing the power spectral density ('density') where `Pxx` has units of V**2/Hz and computing the squared magnitude spectrum ('spectrum') where `Pxx` has units of V**2, if `x` is measured in V and `fs` is measured in Hz. Defaults to 'density' axis : int, optional Axis along which the periodogram is computed; the default is over the last axis (i.e. ``axis=-1``). average : { 'mean', 'median' }, optional Method to use when averaging periodograms. Defaults to 'mean'. .. versionadded:: 1.2.0 Returns ------- f : ndarray Array of sample frequencies. Pxx : ndarray Power spectral density or power spectrum of x. See Also -------- csd: Cross power spectral density using Welch's method periodogram: Simple, optionally modified periodogram lombscargle: Lomb-Scargle periodogram for unevenly sampled data Notes ----- An appropriate amount of overlap will depend on the choice of window and on your requirements. For the default Hann window an overlap of 50% is a reasonable trade-off between accurately estimating the signal power, while not over counting any of the data. Narrower windows may require a larger overlap. If `noverlap` is 0, this method is equivalent to Bartlett's method [2]_. The ratio of the squared magnitude (``scaling='spectrum'``) divided by the spectral power density (``scaling='density'``) is the constant factor of ``sum(abs(window)**2)*fs / abs(sum(window))**2``. If `return_onesided` is ``True``, the values of the negative frequencies are added to values of the corresponding positive ones. Consult the :ref:`tutorial_SpectralAnalysis` section of the :ref:`user_guide` for a discussion of the scalings of the power spectral density and the (squared) magnitude spectrum. .. versionadded:: 0.12.0 References ---------- .. [1] P. Welch, "The use of the fast Fourier transform for the estimation of power spectra: A method based on time averaging over short, modified periodograms", IEEE Trans. Audio Electroacoust. vol. 15, pp. 70-73, 1967. .. [2] M.S. Bartlett, "Periodogram Analysis and Continuous Spectra", Biometrika, vol. 37, pp. 1-16, 1950. Examples -------- >>> import numpy as np >>> from scipy import signal >>> import matplotlib.pyplot as plt >>> rng = np.random.default_rng() Generate a test signal, a 2 Vrms sine wave at 1234 Hz, corrupted by 0.001 V**2/Hz of white noise sampled at 10 kHz. >>> fs = 10e3 >>> N = 1e5 >>> amp = 2*np.sqrt(2) >>> freq = 1234.0 >>> noise_power = 0.001 * fs / 2 >>> time = np.arange(N) / fs >>> x = amp*np.sin(2*np.pi*freq*time) >>> x += rng.normal(scale=np.sqrt(noise_power), size=time.shape) Compute and plot the power spectral density. >>> f, Pxx_den = signal.welch(x, fs, nperseg=1024) >>> plt.semilogy(f, Pxx_den) >>> plt.ylim([0.5e-3, 1]) >>> plt.xlabel('frequency [Hz]') >>> plt.ylabel('PSD [V**2/Hz]') >>> plt.show() If we average the last half of the spectral density, to exclude the peak, we can recover the noise power on the signal. >>> np.mean(Pxx_den[256:]) 0.0009924865443739191 Now compute and plot the power spectrum. >>> f, Pxx_spec = signal.welch(x, fs, 'flattop', 1024, scaling='spectrum') >>> plt.figure() >>> plt.semilogy(f, np.sqrt(Pxx_spec)) >>> plt.xlabel('frequency [Hz]') >>> plt.ylabel('Linear spectrum [V RMS]') >>> plt.show() The peak height in the power spectrum is an estimate of the RMS amplitude. >>> np.sqrt(Pxx_spec.max()) 2.0077340678640727 If we now introduce a discontinuity in the signal, by increasing the amplitude of a small portion of the signal by 50, we can see the corruption of the mean average power spectral density, but using a median average better estimates the normal behaviour. >>> x[int(N//2):int(N//2)+10] *= 50. >>> f, Pxx_den = signal.welch(x, fs, nperseg=1024) >>> f_med, Pxx_den_med = signal.welch(x, fs, nperseg=1024, average='median') >>> plt.semilogy(f, Pxx_den, label='mean') >>> plt.semilogy(f_med, Pxx_den_med, label='median') >>> plt.ylim([0.5e-3, 1]) >>> plt.xlabel('frequency [Hz]') >>> plt.ylabel('PSD [V**2/Hz]') >>> plt.legend() >>> plt.show() ) rZr[r\r]r^r_r`rarbaverage)rreal) rrZr[r\r]r^r_r`rarbrjrPxxs rVrrs>NQb&T7%4g1JE3 #((?rWc  ||ut| } } tj|}| stj|}t|jt|j}} | j | |j |  tj | |}|jdk(s|jdk(rY|t|j| |j| gfz}tjtj|d| }||fStj||tj}nq|jdk(r>tj|jtj|jfStj|tj}| r|j| n%t|j| |j| }t|t st|t"r[| t|nd}|dkrtd|d||kr&t%j&d |d d |d |zd |}t)||}n"tj|}| t+|}|t+|k7rtd |dt+|| t|n|}||krtd|d|d| t|n|dz}||k\rtd|d|dtj,|r|rd}|j| |j| krdt|j}|j| |j| z || <tj.|tj0|f| }n|j| |j| krct|j}|j| |j| z || <tj.|tj0|f| }t3t4|rdnd}| dddx}vrtd| d|dt7|||z ||||| d}|j9|||durdn|d||z |j:z|dz| }|rl| dkr|j<dz | zn| }tj||d}|d!d|j>dzdk(rdndfxxdzcc<tj|d|}|jddkDr| d"k(rtA|jd}tj,|rYtjBtjD|dtjBtjF|dd#zz}ntjB|d}||z}nJ| d$k(r|jId}n2td%| d&tjJ||jdd}|jM|}| r |jD}|jN|fS#t$r}td|d}~wwxYw)'u Estimate the cross power spectral density, Pxy, using Welch's method. Parameters ---------- x : array_like Time series of measurement values y : array_like Time series of measurement values fs : float, optional Sampling frequency of the `x` and `y` time series. Defaults to 1.0. window : str or tuple or array_like, optional Desired window to use. If `window` is a string or tuple, it is passed to `get_window` to generate the window values, which are DFT-even by default. See `get_window` for a list of windows and required parameters. If `window` is array_like it will be used directly as the window and its length must be nperseg. Defaults to a Hann window. nperseg : int, optional Length of each segment. Defaults to None, but if window is str or tuple, is set to 256, and if window is array_like, is set to the length of the window. noverlap: int, optional Number of points to overlap between segments. If `None`, ``noverlap = nperseg // 2``. Defaults to `None` and may not be greater than `nperseg`. nfft : int, optional Length of the FFT used, if a zero padded FFT is desired. If `None`, the FFT length is `nperseg`. Defaults to `None`. detrend : str or function or `False`, optional Specifies how to detrend each segment. If `detrend` is a string, it is passed as the `type` argument to the `detrend` function. If it is a function, it takes a segment and returns a detrended segment. If `detrend` is `False`, no detrending is done. Defaults to 'constant'. return_onesided : bool, optional If `True`, return a one-sided spectrum for real data. If `False` return a two-sided spectrum. Defaults to `True`, but for complex data, a two-sided spectrum is always returned. scaling : { 'density', 'spectrum' }, optional Selects between computing the cross spectral density ('density') where `Pxy` has units of V**2/Hz and computing the cross spectrum ('spectrum') where `Pxy` has units of V**2, if `x` and `y` are measured in V and `fs` is measured in Hz. Defaults to 'density' axis : int, optional Axis along which the CSD is computed for both inputs; the default is over the last axis (i.e. ``axis=-1``). average : { 'mean', 'median' }, optional Method to use when averaging periodograms. If the spectrum is complex, the average is computed separately for the real and imaginary parts. Defaults to 'mean'. .. versionadded:: 1.2.0 Returns ------- f : ndarray Array of sample frequencies. Pxy : ndarray Cross spectral density or cross power spectrum of x,y. See Also -------- periodogram: Simple, optionally modified periodogram lombscargle: Lomb-Scargle periodogram for unevenly sampled data welch: Power spectral density by Welch's method. [Equivalent to csd(x,x)] coherence: Magnitude squared coherence by Welch's method. Notes ----- By convention, Pxy is computed with the conjugate FFT of X multiplied by the FFT of Y. If the input series differ in length, the shorter series will be zero-padded to match. An appropriate amount of overlap will depend on the choice of window and on your requirements. For the default Hann window an overlap of 50% is a reasonable trade-off between accurately estimating the signal power, while not over counting any of the data. Narrower windows may require a larger overlap. The ratio of the cross spectrum (``scaling='spectrum'``) divided by the cross spectral density (``scaling='density'``) is the constant factor of ``sum(abs(window)**2)*fs / abs(sum(window))**2``. If `return_onesided` is ``True``, the values of the negative frequencies are added to values of the corresponding positive ones. Consult the :ref:`tutorial_SpectralAnalysis` section of the :ref:`user_guide` for a discussion of the scalings of a spectral density and an (amplitude) spectrum. Welch's method may be interpreted as taking the average over the time slices of a (cross-) spectrogram. Internally, this function utilizes the `ShortTimeFFT` to determine the required (cross-) spectrogram. An example below illustrates that it is straightforward to calculate `Pxy` directly with the `ShortTimeFFT`. However, there are some notable differences in the behavior of the `ShortTimeFFT`: * There is no direct `ShortTimeFFT` equivalent for the `csd` parameter combination ``return_onesided=True, scaling='density'``, since ``fft_mode='onesided2X'`` requires ``'psd'`` scaling. The is due to `csd` returning the doubled squared magnitude in this case, which does not have a sensible interpretation. * `ShortTimeFFT` uses `float64` / `complex128` internally, which is due to the behavior of the utilized `~scipy.fft` module. Thus, those are the dtypes being returned. The `csd` function casts the return values to `float32` / `complex64` if the input is `float32` / `complex64` as well. * The `csd` function calculates ``np.conj(Sx[q,p]) * Sy[q,p]``, whereas `~ShortTimeFFT.spectrogram` calculates ``Sx[q,p] * np.conj(Sy[q,p])`` where ``Sx[q,p]``, ``Sy[q,p]`` are the STFTs of `x` and `y`. Also, the window positioning is different. .. versionadded:: 0.16.0 References ---------- .. [1] P. Welch, "The use of the fast Fourier transform for the estimation of power spectra: A method based on time averaging over short, modified periodograms", IEEE Trans. Audio Electroacoust. vol. 15, pp. 70-73, 1967. .. [2] Rabiner, Lawrence R., and B. Gold. "Theory and Application of Digital Signal Processing" Prentice-Hall, pp. 414-419, 1975 Examples -------- The following example plots the cross power spectral density of two signals with some common features: >>> import numpy as np >>> from scipy import signal >>> import matplotlib.pyplot as plt >>> rng = np.random.default_rng() ... ... # Generate two test signals with some common features: >>> N, fs = 100_000, 10e3 # number of samples and sampling frequency >>> amp, freq = 20, 1234.0 # amplitude and frequency of utilized sine signal >>> noise_power = 0.001 * fs / 2 >>> time = np.arange(N) / fs >>> b, a = signal.butter(2, 0.25, 'low') >>> x = rng.normal(scale=np.sqrt(noise_power), size=time.shape) >>> y = signal.lfilter(b, a, x) >>> x += amp*np.sin(2*np.pi*freq*time) >>> y += rng.normal(scale=0.1*np.sqrt(noise_power), size=time.shape) ... ... # Compute and plot the magnitude of the cross spectral density: >>> nperseg, noverlap, win = 1024, 512, 'hann' >>> f, Pxy = signal.csd(x, y, fs, win, nperseg, noverlap) >>> fig0, ax0 = plt.subplots(tight_layout=True) >>> ax0.set_title(f"CSD ({win.title()}-window, {nperseg=}, {noverlap=})") >>> ax0.set(xlabel="Frequency $f$ in kHz", ylabel="CSD Magnitude in V²/Hz") >>> ax0.semilogy(f/1e3, np.abs(Pxy)) >>> ax0.grid() >>> plt.show() The cross spectral density is calculated by taking the average over the time slices of a spectrogram: >>> SFT = signal.ShortTimeFFT.from_window('hann', fs, nperseg, noverlap, ... scale_to='psd', fft_mode='onesided2X', ... phase_shift=None) >>> Sxy1 = SFT.spectrogram(y, x, detr='constant', k_offset=nperseg//2, ... p0=0, p1=(N-noverlap) // SFT.hop) >>> Pxy1 = Sxy1.mean(axis=-1) >>> np.allclose(Pxy, Pxy1) # same result as with csd() True As discussed in the Notes section, the results of using an approach analogous to the code snippet above and the `csd` function may deviate due to implementation details. Note that the code snippet above can be easily adapted to determine other statistical properties than the mean value. %x and y cannot be broadcast together.Nrr(rzParameter nperseg=z is not a positive integer!znperseg=z+ is greater than signal length max(len(x), z len(y)) = , using nperseg =  stacklevelz does not equal len(win)=znfft=z* must be greater than or equal to nperseg=!z noverlap=z must be less than nperseg=Frbonesidedtwosided magnitudepsd)spectrumdensityParameter scaling=z not in )fft_modemfftscale_to phase_shift)detrp0p1k_offsetrb.medianr)r6zParameter average=z must be 'median' or 'mean'!)(intr*r-listr0popbroadcast_shapesr1r/minmoveaxisrc result_type complex64maxr4strrfwarningswarnrre iscomplexobj concatenatezerosr rrrhopr.r _median_biasrrkimagr6r7astypef)rrrZr[r\r]r^r_r`rarbrj same_datax_outery_outer outer_shapee out_shape empty_out out_dtypenwinz_shaper~scalesSFTPxyf_axisbiass rVrrsf1fc$itI 1 A  JJqM M4= D D M--gw?K 66Q;!&&A+#sAGGDM1774=+I'J&LLI BHHY$7TBIi' 'NN1a6 66Q;88AGG$bhhqww&77 7NN1bll3 " AGGDM1774=(IA&#*VU";")"5#g,3 Q;2'3NOP P [ MMXWJ&QR&qc);A3?@LM OG)jj  ?#hG#c(HG:%?c#h[ABB(3t9gD g~ED7"MWJaPQQ ( 4s8}'Q,H7IH;&B'1EFF qowwt}qwwt}$qww-  5 NNArxx01 =  &qww-  5 NNArxx01 =M:jQHke!LLvM.gZxxqABB sGh.XD &wT CC //!QW-=T7ALSWW#` between the implementations can be found in the :ref:`tutorial_stft` section of the :ref:`user_guide`. Parameters ---------- x : array_like Time series of measurement values fs : float, optional Sampling frequency of the `x` time series. Defaults to 1.0. window : str or tuple or array_like, optional Desired window to use. If `window` is a string or tuple, it is passed to `get_window` to generate the window values, which are DFT-even by default. See `get_window` for a list of windows and required parameters. If `window` is array_like it will be used directly as the window and its length must be nperseg. Defaults to a Tukey window with shape parameter of 0.25. nperseg : int, optional Length of each segment. Defaults to None, but if window is str or tuple, is set to 256, and if window is array_like, is set to the length of the window. noverlap : int, optional Number of points to overlap between segments. If `None`, ``noverlap = nperseg // 8``. Defaults to `None`. nfft : int, optional Length of the FFT used, if a zero padded FFT is desired. If `None`, the FFT length is `nperseg`. Defaults to `None`. detrend : str or function or `False`, optional Specifies how to detrend each segment. If `detrend` is a string, it is passed as the `type` argument to the `detrend` function. If it is a function, it takes a segment and returns a detrended segment. If `detrend` is `False`, no detrending is done. Defaults to 'constant'. return_onesided : bool, optional If `True`, return a one-sided spectrum for real data. If `False` return a two-sided spectrum. Defaults to `True`, but for complex data, a two-sided spectrum is always returned. scaling : { 'density', 'spectrum' }, optional Selects between computing the power spectral density ('density') where `Sxx` has units of V**2/Hz and computing the power spectrum ('spectrum') where `Sxx` has units of V**2, if `x` is measured in V and `fs` is measured in Hz. Defaults to 'density'. axis : int, optional Axis along which the spectrogram is computed; the default is over the last axis (i.e. ``axis=-1``). mode : str, optional Defines what kind of return values are expected. Options are ['psd', 'complex', 'magnitude', 'angle', 'phase']. 'complex' is equivalent to the output of `stft` with no padding or boundary extension. 'magnitude' returns the absolute magnitude of the STFT. 'angle' and 'phase' return the complex angle of the STFT, with and without unwrapping, respectively. Returns ------- f : ndarray Array of sample frequencies. t : ndarray Array of segment times. Sxx : ndarray Spectrogram of x. By default, the last axis of Sxx corresponds to the segment times. See Also -------- periodogram: Simple, optionally modified periodogram lombscargle: Lomb-Scargle periodogram for unevenly sampled data welch: Power spectral density by Welch's method. csd: Cross spectral density by Welch's method. ShortTimeFFT: Newer STFT/ISTFT implementation providing more features, which also includes a :meth:`~ShortTimeFFT.spectrogram` method. Notes ----- An appropriate amount of overlap will depend on the choice of window and on your requirements. In contrast to welch's method, where the entire data stream is averaged over, one may wish to use a smaller overlap (or perhaps none at all) when computing a spectrogram, to maintain some statistical independence between individual segments. It is for this reason that the default window is a Tukey window with 1/8th of a window's length overlap at each end. .. versionadded:: 0.16.0 References ---------- .. [1] Oppenheim, Alan V., Ronald W. Schafer, John R. Buck "Discrete-Time Signal Processing", Prentice Hall, 1999. Examples -------- >>> import numpy as np >>> from scipy import signal >>> from scipy.fft import fftshift >>> import matplotlib.pyplot as plt >>> rng = np.random.default_rng() Generate a test signal, a 2 Vrms sine wave whose frequency is slowly modulated around 3kHz, corrupted by white noise of exponentially decreasing magnitude sampled at 10 kHz. >>> fs = 10e3 >>> N = 1e5 >>> amp = 2 * np.sqrt(2) >>> noise_power = 0.01 * fs / 2 >>> time = np.arange(N) / float(fs) >>> mod = 500*np.cos(2*np.pi*0.25*time) >>> carrier = amp * np.sin(2*np.pi*3e3*time + mod) >>> noise = rng.normal(scale=np.sqrt(noise_power), size=time.shape) >>> noise *= np.exp(-time/5) >>> x = carrier + noise Compute and plot the spectrogram. >>> f, t, Sxx = signal.spectrogram(x, fs) >>> plt.pcolormesh(t, f, Sxx, shading='gouraud') >>> plt.ylabel('Frequency [Hz]') >>> plt.xlabel('Time [sec]') >>> plt.show() Note, if using output that is not one sided, then use the following: >>> f, t, Sxx = signal.spectrogram(x, fs, return_onesided=False) >>> plt.pcolormesh(t, fftshift(f), fftshift(Sxx, axes=0), shading='gouraud') >>> plt.ylabel('Frequency [Hz]') >>> plt.xlabel('Time [sec]') >>> plt.show() )rzcomplexryanglephasezunknown value for mode z, must be one of  input_lengthrz)moderry)rrrrrrv)r1_triage_segmentsr0_spectral_helperr*absrunwrap)rrZr[r\r]r^r_r`rarbrmodelistrtimeSxxs rVrrs"`AH 824&8I(TUU'vw45GGDMCOFGa< u}+Aq"fg,4dG,;Wd168tS ,Aq"fg,4dG,;Wd179tS ; &&+C ' '((3-Cw!8AIDii$/ $ rWct|}|dkr td||k\r tdt|}t|tst |t ur t ||nUtj|tjdk7r tdjd|k7r td||z tfdt|zD}|zdk7r|d|zxxx|z dz ccc|tj|z }tjtj||kS) a6Check whether the Constant OverLap Add (COLA) constraint is met (legacy function). .. legacy:: function The COLA constraint is equivalent of having a constant dual window, i.e., ``all(ShortTimeFFT.dual_win == ShortTimeFFT.dual_win[0])``. Hence, `closest_STFT_dual_window` generalizes this function, as the following example shows: >>> import numpy as np >>> from scipy.signal import check_COLA, closest_STFT_dual_window, windows ... >>> w, w_rect, hop = windows.hann(12, sym=False), np.ones(12), 6 >>> dual_win, alpha = closest_STFT_dual_window(w, hop, w_rect, scaled=True) >>> np.allclose(dual_win/alpha, w_rect, atol=1e-10, rtol=0) True >>> check_COLA(w, len(w), len(w) - hop) # equivalent legacy function call True Parameters ---------- window : str or tuple or array_like Desired window to use. If `window` is a string or tuple, it is passed to `get_window` to generate the window values, which are DFT-even by default. See `get_window` for a list of windows and required parameters. If `window` is array_like it will be used directly as the window and its length must be nperseg. nperseg : int Length of each segment. noverlap : int Number of points to overlap between segments. tol : float, optional The allowed variance of a bin's weighted sum from the median bin sum. Returns ------- verdict : bool `True` if chosen combination satisfies COLA within `tol`, `False` otherwise See Also -------- closest_STFT_dual_window: Allows determining the closest window meeting the COLA constraint for a given window check_NOLA: Check whether the Nonzero Overlap Add (NOLA) constraint is met ShortTimeFFT: Provide short-time Fourier transform and its inverse stft: Short-time Fourier transform (legacy) istft: Inverse Short-time Fourier transform (legacy) Notes ----- In order to invert a short-time Fourier transfrom (STFT) with the so-called "overlap-add method", the signal windowing must obey the constraint of "Constant OverLap Add" (COLA). This ensures that every point in the input data is equally weighted, thereby avoiding aliasing and allowing full reconstruction. Note that the algorithms implemented in `ShortTimeFFT.istft` and in `istft` (legacy) only require that the weaker "nonzero overlap-add" condition (as in `check_NOLA`) is met. Some examples of windows that satisfy COLA: - Rectangular window at overlap of 0, 1/2, 2/3, 3/4, ... - Bartlett window at overlap of 1/2, 3/4, 5/6, ... - Hann window at 1/2, 2/3, 3/4, ... - Any Blackman family window at 2/3 overlap - Any window with ``noverlap = nperseg-1`` A very comprehensive list of other windows may be found in [2]_, wherein the COLA condition is satisfied when the "Amplitude Flatness" is unity. .. versionadded:: 0.19.0 References ---------- .. [1] Julius O. Smith III, "Spectral Audio Signal Processing", W3K Publishing, 2011,ISBN 978-0-9745607-3-1. .. [2] G. Heinzel, A. Ruediger and R. Schilling, "Spectrum and spectral density estimation by the Discrete Fourier transform (DFT), including a comprehensive list of window functions and some new at-top windows", 2002, http://hdl.handle.net/11858/00-001M-0000-0013-557A-5 Examples -------- >>> from scipy import signal Confirm COLA condition for rectangular window of 75% (3/4) overlap: >>> signal.check_COLA(signal.windows.boxcar(100), 100, 75) True COLA is not true for 25% (1/4) overlap, though: >>> signal.check_COLA(signal.windows.boxcar(100), 100, 25) False "Symmetrical" Hann window (for filter design) is not COLA: >>> signal.check_COLA(signal.windows.hann(120, sym=True), 120, 60) False "Periodic" or "DFT-even" Hann window (for FFT analysis) is COLA for overlap of 1/2, 2/3, 3/4, etc.: >>> signal.check_COLA(signal.windows.hann(120, sym=False), 120, 60) True >>> signal.check_COLA(signal.windows.hann(120, sym=False), 120, 80) True >>> signal.check_COLA(signal.windows.hann(120, sym=False), 120, 90) True r"nperseg must be a positive integer#noverlap must be less than nperseg.window must be 1-Dr"window must have length of npersegc3:K|]}|z|dzzyw)rN.0iisteprs rV zcheck_COLA..s$Jr#bgr!tTk*JsN)rr1r4rtyperfrr*r-rer0r3rangerrr)r[r\r]tolbinsums deviationrrs @@rVrrts*l'lG{=>>7>??8}H&#$v,%"7)jj  syy>Q 12 2 99Q<7 "AB B X DJU7D=5IJJG~4 C'D.(9(:$;; "))G,,I 66"&&# $s **rWcZt|}|dkr td||k\r td|dkr tdt|}t|tst |t ur t ||nUtj|tjdk7r tdjd|k7r td||z tfdt|zD}|zdk7r|d |zxxx|z d d zz ccctj||kDS) aH Check whether the Nonzero Overlap Add (NOLA) constraint is met. Parameters ---------- window : str or tuple or array_like Desired window to use. If `window` is a string or tuple, it is passed to `get_window` to generate the window values, which are DFT-even by default. See `get_window` for a list of windows and required parameters. If `window` is array_like it will be used directly as the window and its length must be nperseg. nperseg : int Length of each segment. noverlap : int Number of points to overlap between segments. tol : float, optional The allowed variance of a bin's weighted sum from the median bin sum. Returns ------- verdict : bool `True` if chosen combination satisfies the NOLA constraint within `tol`, `False` otherwise See Also -------- check_COLA: Check whether the Constant OverLap Add (COLA) constraint is met stft: Short Time Fourier Transform istft: Inverse Short Time Fourier Transform Notes ----- In order to enable inversion of an STFT via the inverse STFT in `istft`, the signal windowing must obey the constraint of "nonzero overlap add" (NOLA): .. math:: \sum_{t}w^{2}[n-tH] \ne 0 for all :math:`n`, where :math:`w` is the window function, :math:`t` is the frame index, and :math:`H` is the hop size (:math:`H` = `nperseg` - `noverlap`). This ensures that the normalization factors in the denominator of the overlap-add inversion equation are not zero. Only very pathological windows will fail the NOLA constraint. .. versionadded:: 1.2.0 References ---------- .. [1] Julius O. Smith III, "Spectral Audio Signal Processing", W3K Publishing, 2011,ISBN 978-0-9745607-3-1. .. [2] G. Heinzel, A. Ruediger and R. Schilling, "Spectrum and spectral density estimation by the Discrete Fourier transform (DFT), including a comprehensive list of window functions and some new at-top windows", 2002, http://hdl.handle.net/11858/00-001M-0000-0013-557A-5 Examples -------- >>> import numpy as np >>> from scipy import signal Confirm NOLA condition for rectangular window of 75% (3/4) overlap: >>> signal.check_NOLA(signal.windows.boxcar(100), 100, 75) True NOLA is also true for 25% (1/4) overlap: >>> signal.check_NOLA(signal.windows.boxcar(100), 100, 25) True "Symmetrical" Hann window (for filter design) is also NOLA: >>> signal.check_NOLA(signal.windows.hann(120, sym=True), 120, 60) True As long as there is overlap, it takes quite a pathological window to fail NOLA: >>> w = np.ones(64, dtype="float") >>> w[::2] = 0 >>> signal.check_NOLA(w, 64, 32) False If there is not enough overlap, a window with zeros at the ends will not work: >>> signal.check_NOLA(signal.windows.hann(64), 64, 0) False >>> signal.check_NOLA(signal.windows.hann(64), 64, 1) False >>> signal.check_NOLA(signal.windows.hann(64), 64, 2) True rrz"noverlap must be less than npersegrz&noverlap must be a nonnegative integerrrc3@K|]}|z|dzzdzyw)rruNrrs rVrzcheck_NOLA..}s)M"#bgr!tTk*A-MsNru)rr1r4rrrfrr*r-rer0r3rr)r[r\r]rrrrs @@rVrrs%D'lG{=>>7=>>!|ABB8}H&#$v,%"7)jj  syy>Q 12 2 99Q<7 "AB B X DMgtm8LMMG~4 C'D.(9(:$;Q$>> 66'?S  rWc ~| dk(rd} n| dk7rtd| dt|||||||||| | d|| \} } }| | |fS)aCompute the Short Time Fourier Transform (legacy function). STFTs can be used as a way of quantifying the change of a nonstationary signal's frequency and phase content over time. .. legacy:: function `ShortTimeFFT` is a newer STFT / ISTFT implementation with more features. A :ref:`comparison ` between the implementations can be found in the :ref:`tutorial_stft` section of the :ref:`user_guide`. Parameters ---------- x : array_like Time series of measurement values fs : float, optional Sampling frequency of the `x` time series. Defaults to 1.0. window : str or tuple or array_like, optional Desired window to use. If `window` is a string or tuple, it is passed to `get_window` to generate the window values, which are DFT-even by default. See `get_window` for a list of windows and required parameters. If `window` is array_like it will be used directly as the window and its length must be nperseg. Defaults to a Hann window. nperseg : int, optional Length of each segment. Defaults to 256. noverlap : int, optional Number of points to overlap between segments. If `None`, ``noverlap = nperseg // 2``. Defaults to `None`. When specified, the COLA constraint must be met (see Notes below). nfft : int, optional Length of the FFT used, if a zero padded FFT is desired. If `None`, the FFT length is `nperseg`. Defaults to `None`. detrend : str or function or `False`, optional Specifies how to detrend each segment. If `detrend` is a string, it is passed as the `type` argument to the `detrend` function. If it is a function, it takes a segment and returns a detrended segment. If `detrend` is `False`, no detrending is done. Defaults to `False`. return_onesided : bool, optional If `True`, return a one-sided spectrum for real data. If `False` return a two-sided spectrum. Defaults to `True`, but for complex data, a two-sided spectrum is always returned. boundary : str or None, optional Specifies whether the input signal is extended at both ends, and how to generate the new values, in order to center the first windowed segment on the first input point. This has the benefit of enabling reconstruction of the first input point when the employed window function starts at zero. Valid options are ``['even', 'odd', 'constant', 'zeros', None]``. Defaults to 'zeros', for zero padding extension. I.e. ``[1, 2, 3, 4]`` is extended to ``[0, 1, 2, 3, 4, 0]`` for ``nperseg=3``. padded : bool, optional Specifies whether the input signal is zero-padded at the end to make the signal fit exactly into an integer number of window segments, so that all of the signal is included in the output. Defaults to `True`. Padding occurs after boundary extension, if `boundary` is not `None`, and `padded` is `True`, as is the default. axis : int, optional Axis along which the STFT is computed; the default is over the last axis (i.e. ``axis=-1``). scaling: {'spectrum', 'psd'} The default 'spectrum' scaling allows each frequency line of `Zxx` to be interpreted as a magnitude spectrum. The 'psd' option scales each line to a power spectral density - it allows to calculate the signal's energy by numerically integrating over ``abs(Zxx)**2``. .. versionadded:: 1.9.0 Returns ------- f : ndarray Array of sample frequencies. t : ndarray Array of segment times. Zxx : ndarray STFT of `x`. By default, the last axis of `Zxx` corresponds to the segment times. See Also -------- istft: Inverse Short Time Fourier Transform ShortTimeFFT: Newer STFT/ISTFT implementation providing more features. check_COLA: Check whether the Constant OverLap Add (COLA) constraint is met check_NOLA: Check whether the Nonzero Overlap Add (NOLA) constraint is met welch: Power spectral density by Welch's method. spectrogram: Spectrogram by Welch's method. csd: Cross spectral density by Welch's method. lombscargle: Lomb-Scargle periodogram for unevenly sampled data Notes ----- In order to enable inversion of an STFT via the inverse STFT in `istft`, the signal windowing must obey the constraint of "Nonzero OverLap Add" (NOLA), and the input signal must have complete windowing coverage (i.e. ``(x.shape[axis] - nperseg) % (nperseg-noverlap) == 0``). The `padded` argument may be used to accomplish this. Given a time-domain signal :math:`x[n]`, a window :math:`w[n]`, and a hop size :math:`H` = `nperseg - noverlap`, the windowed frame at time index :math:`t` is given by .. math:: x_{t}[n]=x[n]w[n-tH] The overlap-add (OLA) reconstruction equation is given by .. math:: x[n]=\frac{\sum_{t}x_{t}[n]w[n-tH]}{\sum_{t}w^{2}[n-tH]} The NOLA constraint ensures that every normalization term that appears in the denominator of the OLA reconstruction equation is nonzero. Whether a choice of `window`, `nperseg`, and `noverlap` satisfy this constraint can be tested with `check_NOLA`. .. versionadded:: 0.19.0 References ---------- .. [1] Oppenheim, Alan V., Ronald W. Schafer, John R. Buck "Discrete-Time Signal Processing", Prentice Hall, 1999. .. [2] Daniel W. Griffin, Jae S. Lim "Signal Estimation from Modified Short-Time Fourier Transform", IEEE 1984, 10.1109/TASSP.1984.1164317 Examples -------- >>> import numpy as np >>> from scipy import signal >>> import matplotlib.pyplot as plt >>> rng = np.random.default_rng() Generate a test signal, a 2 Vrms sine wave whose frequency is slowly modulated around 3kHz, corrupted by white noise of exponentially decreasing magnitude sampled at 10 kHz. >>> fs = 10e3 >>> N = 1e5 >>> amp = 2 * np.sqrt(2) >>> noise_power = 0.01 * fs / 2 >>> time = np.arange(N) / float(fs) >>> mod = 500*np.cos(2*np.pi*0.25*time) >>> carrier = amp * np.sin(2*np.pi*3e3*time + mod) >>> noise = rng.normal(scale=np.sqrt(noise_power), ... size=time.shape) >>> noise *= np.exp(-time/5) >>> x = carrier + noise Compute and plot the STFT's magnitude. >>> f, t, Zxx = signal.stft(x, fs, nperseg=1000) >>> plt.pcolormesh(t, f, np.abs(Zxx), vmin=0, vmax=amp, shading='gouraud') >>> plt.title('STFT Magnitude') >>> plt.ylabel('Frequency [Hz]') >>> plt.xlabel('Time [sec]') >>> plt.show() Compare the energy of the signal `x` with the energy of its STFT: >>> E_x = sum(x**2) / fs # Energy of x >>> # Calculate a two-sided STFT with PSD scaling: >>> f, t, Zxx = signal.stft(x, fs, nperseg=1000, return_onesided=False, ... scaling='psd') >>> # Integrate numerically over abs(Zxx)**2: >>> df, dt = f[1] - f[0], t[1] - t[0] >>> E_Zxx = sum(np.sum(Zxx.real**2 + Zxx.imag**2, axis=0) * df) * dt >>> # The energy is the same, but the numerical errors are quite large: >>> np.isclose(E_x, E_Zxx, rtol=1e-2) True rzr|r{r} not in ['spectrum', 'psd']!r)rarbrboundarypadded)r1r)rrZr[r\r]r^r_r`rrrbrarrZxxs rVrrspb% J .gZ/KLMM'1b&'8(,g07d-3h/5 7E4 $ rWc n tj|dz}t| } t|}|jdkr t d| |k(r t d|j |} |rd|j | dz z} n|j | } || }nt|}|dkr t d||r || dzk(r|}n| }n||kr t dt|}||dz}n t|}||k\r t d ||z } ||jdz k7s| |jdz k7r| d kr|j| z} |d kr|j|z}t t|j}t|| gd D]}|j|tj||| |gz}t|tst|tur t||}nXtj|}t!|j dk7r t d |j d |k7rt d||rt"j$nt"j&}||d|dd|ddf}|| dz | zz}tj(t |j dd|gz|j*}tj(||j*}tj,|||j*k7r|j/|j*}| dk(r||j1z}n<| dk(r(|tj2|t1|dzzz}nt d| dt| D]A}|d|| z|| z|zfxx|d|f|zz cc<|d|| z|| z|zfxx|dzz cc<C|r |d|dz|dz f}|d|dz|dz f}tj0|dkDt!|k7rt5j6d|sdndzd|tj8|dkD|dz}|r |j:}|jdkDr3||jdz k7r!| |kr|dz}tj<|d|}tj>|j d tA|z }||fS)aQPerform the inverse Short Time Fourier transform (legacy function). .. legacy:: function `ShortTimeFFT` is a newer STFT / ISTFT implementation with more features. A :ref:`comparison ` between the implementations can be found in the :ref:`tutorial_stft` section of the :ref:`user_guide`. Parameters ---------- Zxx : array_like STFT of the signal to be reconstructed. If a purely real array is passed, it will be cast to a complex data type. fs : float, optional Sampling frequency of the time series. Defaults to 1.0. window : str or tuple or array_like, optional Desired window to use. If `window` is a string or tuple, it is passed to `get_window` to generate the window values, which are DFT-even by default. See `get_window` for a list of windows and required parameters. If `window` is array_like it will be used directly as the window and its length must be nperseg. Defaults to a Hann window. Must match the window used to generate the STFT for faithful inversion. nperseg : int, optional Number of data points corresponding to each STFT segment. This parameter must be specified if the number of data points per segment is odd, or if the STFT was padded via ``nfft > nperseg``. If `None`, the value depends on the shape of `Zxx` and `input_onesided`. If `input_onesided` is `True`, ``nperseg=2*(Zxx.shape[freq_axis] - 1)``. Otherwise, ``nperseg=Zxx.shape[freq_axis]``. Defaults to `None`. noverlap : int, optional Number of points to overlap between segments. If `None`, half of the segment length. Defaults to `None`. When specified, the COLA constraint must be met (see Notes below), and should match the parameter used to generate the STFT. Defaults to `None`. nfft : int, optional Number of FFT points corresponding to each STFT segment. This parameter must be specified if the STFT was padded via ``nfft > nperseg``. If `None`, the default values are the same as for `nperseg`, detailed above, with one exception: if `input_onesided` is True and ``nperseg==2*Zxx.shape[freq_axis] - 1``, `nfft` also takes on that value. This case allows the proper inversion of an odd-length unpadded STFT using ``nfft=None``. Defaults to `None`. input_onesided : bool, optional If `True`, interpret the input array as one-sided FFTs, such as is returned by `stft` with ``return_onesided=True`` and `numpy.fft.rfft`. If `False`, interpret the input as a a two-sided FFT. Defaults to `True`. boundary : bool, optional Specifies whether the input signal was extended at its boundaries by supplying a non-`None` ``boundary`` argument to `stft`. Defaults to `True`. time_axis : int, optional Where the time segments of the STFT is located; the default is the last axis (i.e. ``axis=-1``). freq_axis : int, optional Where the frequency axis of the STFT is located; the default is the penultimate axis (i.e. ``axis=-2``). scaling: {'spectrum', 'psd'} The default 'spectrum' scaling allows each frequency line of `Zxx` to be interpreted as a magnitude spectrum. The 'psd' option scales each line to a power spectral density - it allows to calculate the signal's energy by numerically integrating over ``abs(Zxx)**2``. Returns ------- t : ndarray Array of output data times. x : ndarray iSTFT of `Zxx`. See Also -------- stft: Short Time Fourier Transform ShortTimeFFT: Newer STFT/ISTFT implementation providing more features. check_COLA: Check whether the Constant OverLap Add (COLA) constraint is met check_NOLA: Check whether the Nonzero Overlap Add (NOLA) constraint is met Notes ----- In order to enable inversion of an STFT via the inverse STFT with `istft`, the signal windowing must obey the constraint of "nonzero overlap add" (NOLA): .. math:: \sum_{t}w^{2}[n-tH] \ne 0 This ensures that the normalization factors that appear in the denominator of the overlap-add reconstruction equation .. math:: x[n]=\frac{\sum_{t}x_{t}[n]w[n-tH]}{\sum_{t}w^{2}[n-tH]} are not zero. The NOLA constraint can be checked with the `check_NOLA` function. An STFT which has been modified (via masking or otherwise) is not guaranteed to correspond to a exactly realizible signal. This function implements the iSTFT via the least-squares estimation algorithm detailed in [2]_, which produces a signal that minimizes the mean squared error between the STFT of the returned signal and the modified STFT. .. versionadded:: 0.19.0 References ---------- .. [1] Oppenheim, Alan V., Ronald W. Schafer, John R. Buck "Discrete-Time Signal Processing", Prentice Hall, 1999. .. [2] Daniel W. Griffin, Jae S. Lim "Signal Estimation from Modified Short-Time Fourier Transform", IEEE 1984, 10.1109/TASSP.1984.1164317 Examples -------- >>> import numpy as np >>> from scipy import signal >>> import matplotlib.pyplot as plt >>> rng = np.random.default_rng() Generate a test signal, a 2 Vrms sine wave at 50Hz corrupted by 0.001 V**2/Hz of white noise sampled at 1024 Hz. >>> fs = 1024 >>> N = 10*fs >>> nperseg = 512 >>> amp = 2 * np.sqrt(2) >>> noise_power = 0.001 * fs / 2 >>> time = np.arange(N) / float(fs) >>> carrier = amp * np.sin(2*np.pi*50*time) >>> noise = rng.normal(scale=np.sqrt(noise_power), ... size=time.shape) >>> x = carrier + noise Compute the STFT, and plot its magnitude >>> f, t, Zxx = signal.stft(x, fs=fs, nperseg=nperseg) >>> plt.figure() >>> plt.pcolormesh(t, f, np.abs(Zxx), vmin=0, vmax=amp, shading='gouraud') >>> plt.ylim([f[1], f[-1]]) >>> plt.title('STFT Magnitude') >>> plt.ylabel('Frequency [Hz]') >>> plt.xlabel('Time [sec]') >>> plt.yscale('log') >>> plt.show() Zero the components that are 10% or less of the carrier magnitude, then convert back to a time series via inverse STFT >>> Zxx = np.where(np.abs(Zxx) >= amp/10, Zxx, 0) >>> _, xrec = signal.istft(Zxx, fs) Compare the cleaned signal with the original and true carrier signals. >>> plt.figure() >>> plt.plot(time, x, time, xrec, time, carrier) >>> plt.xlim([2, 2.1]) >>> plt.xlabel('Time [sec]') >>> plt.ylabel('Signal') >>> plt.legend(['Carrier + Noise', 'Filtered via STFT', 'True Carrier']) >>> plt.show() Note that the cleaned signal does not start as abruptly as the original, since some of the coefficients of the transient were also removed: >>> plt.figure() >>> plt.plot(time, x, time, xrec, time, carrier) >>> plt.xlim([0, 0.1]) >>> plt.xlabel('Time [sec]') >>> plt.ylabel('Signal') >>> plt.legend(['Carrier + Noise', 'Filtered via STFT', 'True Carrier']) >>> plt.show() yruzInput stft must be at least 2d!z/Must specify differing time and frequency axes!rNr.nfft must be greater than or equal to nperseg.rrT)reverserzwindow must have length of )rbr.r%r{rzr}r绽|=z2NOLA condition failed, STFT may not be invertible.z! Possibly due to missing boundaryrrr'r()!r*r-rr.r1r0rrsortedr transposer4rrrfrresp_fftirfftifftrr&rrr3sqrtrrwhererkraranger@)rrZr[r\r]r^input_onesidedr time_axis freq_axisranseg n_defaultnstepzouteraxrifuncxsubs outputlengthrnormrrs rVrrDsl **S/B CIIII xx!|:;;IJKK 99Y Dsyy+a/0 IIi( g, Q;AB B | IM!9DD IJJ4yA:x=7>?? h ECHHQJ)sxxz"9 q=9,I q=9,IeCHHo&)Y/> B JJrN ll3 9'= =>&#$v,%"7)jj  syy>Q 12 2 99Q<7 ":7)DE E*FLL E #B$ 'XgXq(8 9Ed1fe^+L ciin%|n4EKKHA 88L 4D ~~c5!U[[0jj%*  E  c#q&k)**.gZ/KLMMDk7 #r%x5(( ()U37^c-AA) S"U(2e8G++ +,Q6,7  c7A: m++ ,C!gqjM112 vvdUls4y( @:B2 L $,c **A FF vvz  "9$Q  Ar9-A 99QWWQZ r *D 7NrWc t||||||||\} } t||||||||\} } t||||||||| \} } tj| dz| z | z }| |fS)a Estimate the magnitude squared coherence estimate, Cxy, of discrete-time signals X and Y using Welch's method. ``Cxy = abs(Pxy)**2/(Pxx*Pyy)``, where `Pxx` and `Pyy` are power spectral density estimates of X and Y, and `Pxy` is the cross spectral density estimate of X and Y. Parameters ---------- x : array_like Time series of measurement values y : array_like Time series of measurement values fs : float, optional Sampling frequency of the `x` and `y` time series. Defaults to 1.0. window : str or tuple or array_like, optional Desired window to use. If `window` is a string or tuple, it is passed to `get_window` to generate the window values, which are DFT-even by default. See `get_window` for a list of windows and required parameters. If `window` is array_like it will be used directly as the window and its length must be nperseg. Defaults to a Hann window. nperseg : int, optional Length of each segment. Defaults to None, but if window is str or tuple, is set to 256, and if window is array_like, is set to the length of the window. noverlap: int, optional Number of points to overlap between segments. If `None`, ``noverlap = nperseg // 2``. Defaults to `None`. nfft : int, optional Length of the FFT used, if a zero padded FFT is desired. If `None`, the FFT length is `nperseg`. Defaults to `None`. detrend : str or function or `False`, optional Specifies how to detrend each segment. If `detrend` is a string, it is passed as the `type` argument to the `detrend` function. If it is a function, it takes a segment and returns a detrended segment. If `detrend` is `False`, no detrending is done. Defaults to 'constant'. axis : int, optional Axis along which the coherence is computed for both inputs; the default is over the last axis (i.e. ``axis=-1``). Returns ------- f : ndarray Array of sample frequencies. Cxy : ndarray Magnitude squared coherence of x and y. See Also -------- periodogram: Simple, optionally modified periodogram lombscargle: Lomb-Scargle periodogram for unevenly sampled data welch: Power spectral density by Welch's method. csd: Cross spectral density by Welch's method. Notes ----- An appropriate amount of overlap will depend on the choice of window and on your requirements. For the default Hann window an overlap of 50% is a reasonable trade-off between accurately estimating the signal power, while not over counting any of the data. Narrower windows may require a larger overlap. .. versionadded:: 0.16.0 References ---------- .. [1] P. Welch, "The use of the fast Fourier transform for the estimation of power spectra: A method based on time averaging over short, modified periodograms", IEEE Trans. Audio Electroacoust. vol. 15, pp. 70-73, 1967. .. [2] Stoica, Petre, and Randolph Moses, "Spectral Analysis of Signals" Prentice Hall, 2005 Examples -------- >>> import numpy as np >>> from scipy import signal >>> import matplotlib.pyplot as plt >>> rng = np.random.default_rng() Generate two test signals with some common features. >>> fs = 10e3 >>> N = 1e5 >>> amp = 20 >>> freq = 1234.0 >>> noise_power = 0.001 * fs / 2 >>> time = np.arange(N) / fs >>> b, a = signal.butter(2, 0.25, 'low') >>> x = rng.normal(scale=np.sqrt(noise_power), size=time.shape) >>> y = signal.lfilter(b, a, x) >>> x += amp*np.sin(2*np.pi*freq*time) >>> y += rng.normal(scale=0.1*np.sqrt(noise_power), size=time.shape) Compute and plot the coherence. >>> f, Cxy = signal.coherence(x, y, fs, nperseg=1024) >>> plt.semilogy(f, Cxy) >>> plt.xlabel('frequency [Hz]') >>> plt.ylabel('Coherence') >>> plt.show() )rZr[r\r]r^r_rbru)rrr*r)rrrZr[r\r]r^r_rbrrl_PyyrCxys rVrrrsZqR (tW "JE31FGhgD:FAs A"VW"wTKFAs &&+q.3  $C #:rWc  | dvrtd| dttttdd}| |vr(td| dt |j ||u}|s| dk7r td t  tj|}|s;tj|}tj||tj}n$tj|tj}|st |j}t |j}|j |j  tjtj|tj|j}|rk|j d k(rtj|jtj|jtj|jfS|j d k(s|j d k(rZt#|j |j gfz}tj$tj|d }|||fS|j&d kDrD d k7r?tj$| d }|s&|j&d kDrtj$| d }|s|jd |jd k7r|jd |jd krct |j}|jd |jd z |d <tj(|tj*|fd }nbt |j}|jd |jd z |d <tj(|tj*|fd }|t|}|d kr tdt-|||jd \}}||}n||kr tdt|}||dz}n t|}||k\r td||z }| #|| }|||dzd }|s|||dzd }| r|jd |z |z|z}t |jdd |gz}tj(|tj*|fd }|sHt |jdd |gz}tj(|tj*|fd }sd}n t/dsfd}n d k7r fd}n}tj|tj|k7r|j1|}| dk(rd|||zj3zz }n*| dk(rd|j3dzz }ntd| | dk(rtj4|}|rbtj6|rd}t9j:dd n5d!}|s1tj6|rd}t9j:dd nd}|dk(rt=j>|d |z } n|d!k(rt=j@|d |z } tC|||||||}!|s*tC|||||||}"tjD|!|"z}!n| dk(rtjD|!|!z}!|!|z}!|d!k(r-| dk(r(|dzr|!d"d dfxxdzcc<n|!d"d d fxxdzcc<tjF|dz |jd |dz z d z||z tI|z }#| |#|dz |z z}#|!j1|}!|r| dk7r |!jJ}! d kr d z tj$|!d }! |#|!fS#t$r}td |d}~wwxYw)#uCalculate various forms of windowed FFTs for PSD, CSD, etc. .. legacy:: function This function is soley used by the legacy functions `spectrogram` and `stft` (which are also in this same source file `scipy/signal/_spectral_py.py`). This is a helper function that implements the commonality between the stft, psd, csd, and spectrogram functions. It is not designed to be called externally. The windows are not averaged over; the result from each window is returned. Parameters ---------- x : array_like Array or sequence containing the data to be analyzed. y : array_like Array or sequence containing the data to be analyzed. If this is the same object in memory as `x` (i.e. ``_spectral_helper(x, x, ...)``), the extra computations are spared. fs : float, optional Sampling frequency of the time series. Defaults to 1.0. window : str or tuple or array_like, optional Desired window to use. If `window` is a string or tuple, it is passed to `get_window` to generate the window values, which are DFT-even by default. See `get_window` for a list of windows and required parameters. If `window` is array_like it will be used directly as the window and its length must be nperseg. Defaults to a Hann window. nperseg : int, optional Length of each segment. Defaults to None, but if window is str or tuple, is set to 256, and if window is array_like, is set to the length of the window. noverlap : int, optional Number of points to overlap between segments. If `None`, ``noverlap = nperseg // 2``. Defaults to `None`. nfft : int, optional Length of the FFT used, if a zero padded FFT is desired. If `None`, the FFT length is `nperseg`. Defaults to `None`. detrend : str or function or `False`, optional Specifies how to detrend each segment. If `detrend` is a string, it is passed as the `type` argument to the `detrend` function. If it is a function, it takes a segment and returns a detrended segment. If `detrend` is `False`, no detrending is done. Defaults to 'constant'. return_onesided : bool, optional If `True`, return a one-sided spectrum for real data. If `False` return a two-sided spectrum. Defaults to `True`, but for complex data, a two-sided spectrum is always returned. scaling : { 'density', 'spectrum' }, optional Selects between computing the cross spectral density ('density') where `Pxy` has units of V²/Hz and computing the cross spectrum ('spectrum') where `Pxy` has units of V², if `x` and `y` are measured in V and `fs` is measured in Hz. Defaults to 'density' axis : int, optional Axis along which the FFTs are computed; the default is over the last axis (i.e. ``axis=-1``). mode: str {'psd', 'stft'}, optional Defines what kind of return values are expected. Defaults to 'psd'. boundary : str or None, optional Specifies whether the input signal is extended at both ends, and how to generate the new values, in order to center the first windowed segment on the first input point. This has the benefit of enabling reconstruction of the first input point when the employed window function starts at zero. Valid options are ``['even', 'odd', 'constant', 'zeros', None]``. Defaults to `None`. padded : bool, optional Specifies whether the input signal is zero-padded at the end to make the signal fit exactly into an integer number of window segments, so that all of the signal is included in the output. Defaults to `False`. Padding occurs after boundary extension, if `boundary` is not `None`, and `padded` is `True`. Returns ------- freqs : ndarray Array of sample frequencies. t : ndarray Array of times corresponding to each data segment result : ndarray Array of output data, contents dependent on *mode* kwarg. Notes ----- Adapted from matplotlib.mlab .. versionadded:: 0.16.0 )rzrzUnknown value for mode z!, must be one of: {'psd', 'stft'}N)evenoddconstantrNzUnknown boundary option 'z', must be one of: rzz'x and y must be equal if mode is 'stft'rnrr(rrrrrurrvc|S)Nr)ds rV detrend_funcz&_spectral_helper..detrend_funcsHrW__call__c4tj|dS)Nr()rrb)rr_)rr_s rVrz&_spectral_helper..detrend_funcs''bA ArWcptj|d}|}tj|dS)Nr()r*r)rrbr_s rVrz&_spectral_helper..detrend_funcs1 Ar4(A A;;q$+ +rWr|r'r{zUnknown scaling: rrxz9Input data is complex, switching to return_onesided=Falserqrrrw.)&r1r r r r rkeysrr*r-rrr0r broadcastrcr/rrr.rrrrgrr3rrrrrfftfreqrfftfreq _fft_helper conjugaterr@rk)$rrrZr[r\r]r^r_r`rarbrrrboundary_funcsroutdtypexouteryouter outershaperoutshapeemptyout pad_shaperrext_funcnadd zeros_shaperscalesidesrresultresult_yrs$ ` ` rVrrs>~ ?"24&9++, ,'$"+' "N ~%4XJ?,,01D1D1F,G+HJK KQI BCC t9D 1 A  JJqM>>!Q 5>>!R\\2 aggagg 4 4 Mbhhv&68HIOOJ 66Q;88AGG$bhhqww&7!''9JJ J 66Q;!&&A+!S!''$-)G%H$JJH{{288H#5r4@HXx/ /vvz 2: AtR(A!KK4,  772;!''"+ %wwr{QWWR[( M ! aggbk 9 " NNArxx ':#;R@ M ! aggbk 9 " NNArxx ':#;R@g, Q;AB B$FG!''"+NLC | IJJ4yA:x=7>?? h E!(+ Q  ,GQJR0A ''"+g%&.'91773B<(D61 NNArxx 45B ?qwws|,v5K288K#89CA   Wj ) B  ,   ~~c2<<(H4jj")rSWMMO+, J cggil",WK899 v~ ??1 E MMU%& (E??1%&EMM#:-.0 tQrT* * ad+Cw$ NF q#|Wh$&f%0 f%. eOF tu} !8 37Oq O 3"9  "  99WQY gai 7! ;x' )).r 3D b   ]]8 $FTV^ ax  [[T *F $ I MDE1 L MsA^$$ ^>- ^99^>cV|dk(r|dk(r|dtjf}n?||z }tjjj ||dd}|ddd|ddf}||}||z}|dk(rt j } n|j}t j} | || }|S) a Calculate windowed FFT, for internal use by `scipy.signal._spectral_helper`. .. legacy:: function This function is solely used by the legacy `_spectral_helper` function, which is located also in this file. This is a helper function that does the main FFT calculation for `_spectral helper`. All input validation is performed there, and the data axis is assumed to be the last axis of x. It is not designed to be called externally. The windows are not averaged over; the result from each window is returned. Returns ------- result : ndarray Array of FFT data Notes ----- Adapted from matplotlib.mlab .. versionadded:: 0.16.0 rr.r(T) window_shaperb writeableNrx)r) r*newaxislib stride_trickssliding_window_viewrrrkrfft) rrrr\r]r^r r rfuncs rVrr s8!|A 3 ?#!%%99 G": QWWa(& !F6\F zz{{ &D !F MrWc t|tst|tr>|d}||kDr%tjd|dd|dd|dd|}t ||}||fSt j|}t|jdk7r td ||jd kr td ||jd }||fS|||jd k7r td ||fS)aD Parses window and nperseg arguments for spectrogram and _spectral_helper. This is a helper function, not meant to be called externally. .. legacy:: function This function is soley used by the legacy functions `spectrogram` and `_spectral_helper` (which are also in this file). Parameters ---------- window : string, tuple, or ndarray If window is specified by a string or tuple and nperseg is not specified, nperseg is set to the default of 256 and returns a window of that length. If instead the window is array_like and nperseg is not specified, then nperseg is set to the length of the window. A ValueError is raised if the user supplies both an array_like window and a value for nperseg but nperseg does not equal the length of the window. nperseg : int Length of each segment input_length: int Length of input signal, i.e. x.shape[-1]. Used to test for errors. Returns ------- win : ndarray window. If function was called with string or tuple than this will hold the actual array used as a window. nperseg : int Length of each segment. If window is str or tuple, nperseg is set to 256. If window is array_like, nperseg is set to the length of the window. roz nperseg = rz! is greater than input length = rprqrrrrr(z"window is longer than input signalrz>value specified for nperseg is different from length of window) r4rrfrrrr*r-rer0r1)r[r\rrs rVrrN sN&#*VU"; ?G \ ! MMJwqk2 ,Q//A,qAQS%& (#G) <jj  syy>Q 12 2 #))B- 'AB B ?iilG <  #))A,& ":;; <rWcdtjd|dz dzdzz}dtjd|dzz d|z z zS)aG Returns the bias of the median of a set of periodograms relative to the mean. See Appendix B from [1]_ for details. Parameters ---------- n : int Numbers of periodograms being averaged. Returns ------- bias : float Calculated bias. References ---------- .. [1] B. Allen, W.G. Anderson, P.R. Brady, D.A. Brown, J.D.E. Creighton. "FINDCHIRP: an algorithm for detection of gravitational waves from inspiraling compact binaries", Physical Review D 85, 2012, :arxiv:`gr-qc/0509116` rur'r)r*rr3)rii_2s rVrr sJ0 ryyacaZ!^, ,D rvvbD1HoT 12 22rW)FF)r'rYNrTr|r() r'hannNNNrTr|r(r6) r')tukeyg?NNNrTr|r(rz)r) r'rroNNFTrTr(r{) r'rNNNTTr(rr{)r'rNNNrr() r'rNNNrTr|r(rzNF)*__doc__numpyr* numpy.typingtypingnptscipyrrrr_short_time_fftrrwindowsr _arraytoolsr r r r rr r__all__ ArrayLiker5NDArrayrrrrrrrrrrrrrrrrWrVr)s8??  GCH F#'F }}F }}F ==F F g?@@ F[[4 FF [[FR @J>@S-lGK