scipy.signal.lfilter_zi#
- scipy.signal.lfilter_zi(b, a)[source]#
Construct initial conditions for lfilter for step response steady-state.
Compute an initial state zi for the
lfilter
function that corresponds to the steady state of the step response.A typical use of this function is to set the initial state so that the output of the filter starts at the same value as the first element of the signal to be filtered.
- Parameters:
- b, aarray_like (1-D)
The IIR filter coefficients. See
lfilter
for more information.
- Returns:
- zi1-D ndarray
The initial state for the filter.
Notes
A linear filter with order m has a state space representation (A, B, C, D), for which the output y of the filter can be expressed as:
z(n+1) = A*z(n) + B*x(n) y(n) = C*z(n) + D*x(n)
where z(n) is a vector of length m, A has shape (m, m), B has shape (m, 1), C has shape (1, m) and D has shape (1, 1) (assuming x(n) is a scalar). lfilter_zi solves:
zi = A*zi + B
In other words, it finds the initial condition for which the response to an input of all ones is a constant.
Given the filter coefficients a and b, the state space matrices for the transposed direct form II implementation of the linear filter, which is the implementation used by scipy.signal.lfilter, are:
A = scipy.linalg.companion(a).T B = b[1:] - a[1:]*b[0]
assuming a[0] is 1.0; if a[0] is not 1, a and b are first divided by a[0].
Examples
The following code creates a lowpass Butterworth filter. Then it applies that filter to an array whose values are all 1.0; the output is also all 1.0, as expected for a lowpass filter. If the zi argument of
lfilter
had not been given, the output would have shown the transient signal.>>> from numpy import array, ones >>> from scipy.signal import lfilter, lfilter_zi, butter >>> b, a = butter(5, 0.25) >>> zi = lfilter_zi(b, a) >>> y, zo = lfilter(b, a, ones(10), zi=zi) >>> y array([1., 1., 1., 1., 1., 1., 1., 1., 1., 1.])
Another example:
>>> x = array([0.5, 0.5, 0.5, 0.0, 0.0, 0.0, 0.0]) >>> y, zf = lfilter(b, a, x, zi=zi*x[0]) >>> y array([ 0.5 , 0.5 , 0.5 , 0.49836039, 0.48610528, 0.44399389, 0.35505241])
Note that the zi argument to
lfilter
was computed usinglfilter_zi
and scaled by x[0]. Then the output y has no transient until the input drops from 0.5 to 0.0.