istft — SciPy v1.15.3 Manual (original) (raw)
scipy.signal.
scipy.signal.istft(Zxx, fs=1.0, window='hann', nperseg=None, noverlap=None, nfft=None, input_onesided=True, boundary=True, time_axis=-1, freq_axis=-2, scaling='spectrum')[source]#
Perform the inverse Short Time Fourier transform (legacy function).
Legacy
This function is considered legacy and will no longer receive updates. While we currently have no plans to remove it, we recommend that new code uses more modern alternatives instead. ShortTimeFFT is a newer STFT / ISTFT implementation with more features. A comparison between the implementations can be found in the Short-Time Fourier Transform section of theSciPy User Guide.
Parameters:
Zxxarray_like
STFT of the signal to be reconstructed. If a purely real array is passed, it will be cast to a complex data type.
fsfloat, optional
Sampling frequency of the time series. Defaults to 1.0.
windowstr 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.
npersegint, 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.
noverlapint, 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.
nfftint, 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 andnperseg==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_onesidedbool, optional
If True, interpret the input array as one-sided FFTs, such as is returned by stft with return_onesided=True
andnumpy.fft.rfft. If False, interpret the input as a a two-sided FFT. Defaults to True.
boundarybool, optional
Specifies whether the input signal was extended at its boundaries by supplying a non-None boundary
argument tostft. Defaults to True.
time_axisint, optional
Where the time segments of the STFT is located; the default is the last axis (i.e. axis=-1
).
freq_axisint, 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:
tndarray
Array of output data times.
xndarray
iSTFT of Zxx.
See also
Short Time Fourier Transform
Newer STFT/ISTFT implementation providing more features.
Check whether the Constant OverLap Add (COLA) constraint is met
Check whether the Nonzero Overlap Add (NOLA) constraint is met
Notes
In order to enable inversion of an STFT via the inverse STFT withistft, the signal windowing must obey the constraint of “nonzero overlap add” (NOLA):
\[\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
\[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_NOLAfunction.
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.
Added in version 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 = 10fs nperseg = 512 amp = 2 * np.sqrt(2) noise_power = 0.001 * fs / 2 time = np.arange(N) / float(fs) carrier = amp * np.sin(2np.pi50time) 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()