scipy.signal.

resample#

scipy.signal.resample(x, num, t=None, axis=0, window=None, domain='time')[source]#

Resample x to num samples using the Fourier method along the given axis.

The resampling is performed by shortening or zero-padding the FFT of x. This has the advantages of providing an ideal antialiasing filter and allowing arbitrary up- or down-sampling ratios. The main drawback is the requirement of assuming x to be a periodic signal.

Parameters:
xarray_like

The input signal made up of equidistant samples. If x is a multidimensional array, the parameter axis specifies the time/frequency axis. It is assumed here that n_x = x.shape[axis] specifies the number of samples and T the sampling interval.

numint

The number of samples of the resampled output signal. It may be larger or smaller than n_x.

tarray_like, optional

If t is not None, then the timestamps of the resampled signal are also returned. t must contain at least the first two timestamps of the input signal x (all others are ignored). The timestamps of the output signal are determined by t[0] + T * n_x / num * np.arange(num) with T = t[1] - t[0]. Default is None.

axisint, optional

The time/frequency axis of x along which the resampling take place. The Default is 0.

windowarray_like, callable, string, float, or tuple, optional

If not None, it specifies a filter in the Fourier domain, which is applied before resampling. I.e., the FFT X of x is calculated by X = W * fft(x, axis=axis). W may be interpreted as a spectral windowing function W(f_X) which consumes the frequencies f_X = fftfreq(n_x, T).

If window is a 1d array of length n_x then W=window. If window is a callable then W = window(f_X). Otherwise, window is passed to get_window, i.e., W = fftshift(signal.get_window(window, n_x)). Default is None.

domain‘time’ | ‘freq’, optional

If set to 'time' (default) then an FFT is applied to x, otherwise ('freq') it is asssmued that an FFT was already applied, i.e., x = fft(x_t, axis=axis) with x_t being the input signal in the time domain.

Returns:
x_rndarray

The resampled signal made up of num samples and sampling interval T * n_x / num.

t_rndarray, optional

The num equidistant timestamps of x_r. This is only returned if paramater t is not None.

See also

decimate

Downsample a (periodic/non-periodic) signal after applying an FIR or IIR filter.

resample_poly

Resample a (periodic/non-periodic) signal using polyphase filtering and an FIR filter.

Notes

This function uses the more efficient one-sided FFT, i.e. rfft / irfft, if x is real-valued and in the time domain. Else, the two-sided FFT, i.e., fft / ifft, is used (all FFT functions are taken from the scipy.fft module).

If a window is applied to a real-valued x, the one-sided spectral windowing function is determined by taking the average of the negative and the positive frequency component. This ensures that real-valued signals and complex signals with zero imaginary part are treated identically. I.e., passing x or passing x.astype(np.complex128) produce the same numeric result.

If the number of input or output samples are prime or have few prime factors, this function may be slow due to utilizing FFTs. Consult prev_fast_len and next_fast_len for determining efficient signals lengths. Alternatively, utilizing resample_poly to calculate an intermediate signal (as illustrated in the example below) can result in significant speed increases.

resample is intended to be used for periodic signals with equidistant sampling intervals. For non-periodic signals, resample_poly may be a better choice. Consult the scipy.interpolate module for methods of resampling signals with non-constant sampling intervals.

Examples

The following example depicts a signal being up-sampled from 20 samples to 100 samples. The ringing at the beginning of the up-sampled signal is due to interpreting the signal being periodic. The red square in the plot illustrates that periodictiy by showing the first sample of the next cycle of the signal.

>>> import numpy as np
>>> import matplotlib.pyplot as plt
>>> from scipy.signal import resample
...
>>> n0, n1 = 20, 100  # number of samples
>>> t0 = np.linspace(0, 10, n0, endpoint=False)  # input time stamps
>>> x0 = np.cos(-t0**2/6)  # input signal
...
>>> x1 = resample(x0, n1)  # resampled signal
>>> t1 = np.linspace(0, 10, n1, endpoint=False)  # timestamps of x1
...
>>> fig0, ax0 = plt.subplots(1, 1, tight_layout=True)
>>> ax0.set_title(f"Resampling $x(t)$ from {n0} samples to {n1} samples")
>>> ax0.set(xlabel="Time $t$", ylabel="Amplitude $x(t)$")
>>> ax0.plot(t1, x1, '.-', alpha=.5, label=f"Resampled")
>>> ax0.plot(t0, x0, 'o-', alpha=.5, label="Original")
>>> ax0.plot(10, x0[0], 'rs', alpha=.5, label="Next Cycle")
>>> ax0.legend(loc='best')
>>> ax0.grid(True)
>>> plt.show()
../../_images/scipy-signal-resample-1_00_00.png

The following example illustrates that a rfft / irfft combination does not always produce the correct resampling result, while resample does:

>>> import numpy as np
>>> from scipy.fft import irfft, fft, fftfreq, rfft
>>> from scipy.signal import resample
...
>>> n0, n1 = 8, 4  # number of samples
>>> t0, t1 = (np.linspace(0, 2, n_, endpoint=False) for n_ in (n0, n1))  # in s
>>> x0, x1 = (np.cos(2*np.pi*t_) for t_ in (t0, t1))  # 1 Hz cosine signal
>>> y1_r = resample(x0, num=n1)
>>> np.allclose(y1_r, x1)  # correct result, as expected
True
>>> y1_f = irfft(rfft(x0), n=n1) * n1/n0
>>> np.allclose(y1_f, x1)  # wrong result
False
>>> # Compare FFTs:
>>> fftfreq(n0, t0[1]-t0[0])  #  frequencies of fft(x0)
array([ 0. ,  0.5,  1. ,  1.5, -2. , -1.5, -1. , -0.5])
>>> np.round(fft(x0), 3)  # FFT of x0
array([-0.-0.j, -0.+0.j,  4.-0.j,  0.+0.j,  0.-0.j,  0.-0.j,  4.+0.j, -0.-0.j])
>>> fftfreq(n1, t1[1]-t1[0])  # frequencies of fft(x1) and fft(y1_f)
array([ 0. ,  0.5, -1. , -0.5])
>>> np.round(fft(x1), 3)  # reference FFT
array([0.+0.j, 0.+0.j, 4.+0.j, 0.+0.j])
>>> np.round(fft(y1_f), 3)  # irfft/rfft off by factor 2 in the Nyuist frequency
array([0.+0.j, 0.+0.j, 2.+0.j, 0.+0.j])

The reason for the different results lies in resample correctly treating unpaired frequency bins. I.e., the input x1 has a bin pair ±1 Hz, whereas the output has only one unpaired bin at -1 Hz, which demands rescaling of that bin. Special treatment is required if n_x != num and min(n_x, num) is even. If the bin values at ±m are zero, the, obviously, no special treatment is needed. Consult the source code of resample for details.

The following code snippet shows how to use resample_poly to speed up the down-sampling:

>>> import numpy as np
>>> from scipy.signal import resample, resample_poly
...
>>> n0 = 19937 # number of input samples - prime
>>> n1 = 128  # number of output samples - fast FFT length
>>> x0 = np.random.rand(n0)  # input signal
...
>>> y1 = resample(x0, n1)  # slow due to n0 being prime
>>> # This is faster:
>>> y1_p = resample(resample_poly(x0, 1, n0 // n1, padtype='wrap'), n1)

Note that the y1 and y1_p are not identical due to the differences in the utilized antialiasing filter in each function. In both functions the antialiasing filter can be customized by modifying the window parameter, though resample interprets it as a spectral window function, whereas resample_poly assumes it to be some form of impulse response.