detrend#
- scipy.signal.detrend(data, axis=-1, type='linear', bp=0, overwrite_data=False)[source]#
Remove linear or constant trend along axis from data.
- Parameters:
- dataarray_like
The input data.
- axisint, optional
The axis along which to detrend the data. By default this is the last axis (-1).
- type{‘linear’, ‘constant’}, optional
The type of detrending. If
type == 'linear'
(default), the result of a linear least-squares fit to data is subtracted from data. Iftype == 'constant'
, only the mean of data is subtracted.- bparray_like of ints, optional
A sequence of break points. If given, an individual linear fit is performed for each part of data between two break points. Break points are specified as indices into data. This parameter only has an effect when
type == 'linear'
.- overwrite_databool, optional
If True, perform in place detrending and avoid a copy. Default is False
- Returns:
- retndarray
The detrended input data.
See also
numpy.polynomial.polynomial.Polynomial.fit
Create least squares fit polynomial.
Notes
Detrending can be interpreted as subtracting a least squares fit polynomial: Setting the parameter type to ‘constant’ corresponds to fitting a zeroth degree polynomial, ‘linear’ to a first degree polynomial. Consult the example below.
Examples
The following example detrends the function \(x(t) = \sin(\pi t) + 1/4\):
>>> import matplotlib.pyplot as plt >>> import numpy as np >>> from scipy.signal import detrend ... >>> t = np.linspace(-0.5, 0.5, 21) >>> x = np.sin(np.pi*t) + 1/4 ... >>> x_d_const = detrend(x, type='constant') >>> x_d_linear = detrend(x, type='linear') ... >>> fig1, ax1 = plt.subplots() >>> ax1.set_title(r"Detrending $x(t)=\sin(\pi t) + 1/4$") >>> ax1.set(xlabel="t", ylabel="$x(t)$", xlim=(t[0], t[-1])) >>> ax1.axhline(y=0, color='black', linewidth=.5) >>> ax1.axvline(x=0, color='black', linewidth=.5) >>> ax1.plot(t, x, 'C0.-', label="No detrending") >>> ax1.plot(t, x_d_const, 'C1x-', label="type='constant'") >>> ax1.plot(t, x_d_linear, 'C2+-', label="type='linear'") >>> ax1.legend() >>> plt.show()
Alternatively, NumPy’s
Polynomial
can be used for detrending as well:>>> pp0 = np.polynomial.Polynomial.fit(t, x, deg=0) # fit degree 0 polynomial >>> np.allclose(x_d_const, x - pp0(t)) # compare with constant detrend True >>> pp1 = np.polynomial.Polynomial.fit(t, x, deg=1) # fit degree 1 polynomial >>> np.allclose(x_d_linear, x - pp1(t)) # compare with linear detrend True
Note that
Polynomial
also allows fitting higher degree polynomials. Consult its documentation on how to extract the polynomial coefficients.