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:
barray_like

The numerator coefficient vector as a 1-D sequence.

aarray_like

The denominator coefficient vector as a 1-D sequence. If a[0] is not 1, then both a and b are normalized by a[0]. Hence, a[0] != 0 must hold.

Returns:
zi1-D ndarray

The initial state for the filter.

Raises:
ValueError

If a[0] == 0 (invalid denominator polynomial) or sum(a) == 0 (unstable filter).

See also

lfilter, lfiltic, filtfilt

Notes

The parameters b and a represent a transfer function \(H(z) = Y(z)/X(z)\) which is defined in the Transfer function representation section of the SciPy User Guide. As discussed in [1], the final value of filtering a step response \(X(z) = z / (z-1)\), i.e., steady state, is given by

\[y_\infty := \lim_{k\to\infty} y[k] = \lim_{z\to 1}\ (z-1)\, Y(z) = \frac{\sum_{i=0}^M b_i}{\sum_{j=0}^N a_j} \,.\]

If the denominator is zero, \(H(z)\) has a pole at \(z_\infty=1\), which makes the filter unstable. For the transposed Direct Form II, which is implemented in lfilter, the initialization values \(z_k\) for \(H(z)\) can be determined by the recurrence equation

\[z_k = z_{k+1} + x[0] \big( b_k - y_\infty a_k \big) \,,\]

with \(x[0]\) being the height of the input step function. Note that \(a_0=1\) is assumed here, which is incorporated into this function by performing a normalization step.

Array API Standard Support

lfilter_zi has experimental support for Python Array API Standard compatible backends in addition to NumPy. Please consider testing these features by setting an environment variable SCIPY_ARRAY_API=1 and providing CuPy, PyTorch, JAX, or Dask arrays as array arguments. The following combinations of backend and device (or other capability) are supported.

Library

CPU

GPU

NumPy

n/a

CuPy

n/a

PyTorch

JAX

⚠️ no JIT

Dask

⚠️ computes graph

n/a

See Support for the array API standard for more information.

References

[1]

Boris Likhterov and Norman Kopeika. “Hardware-efficient technique for minimizing startup transients in Direct Form II digital filters”. In: International Journal of Electronics – Volume 90(7), July 2003, pp. 471–479. DOI:10.1080/00207210310001612482

Examples

The following code creates a lowpass Butterworth filter to filter a signal made up of ones. As expected of a lowpass filter, the output is also all ones. If the zi argument of lfilter had not been given, a transient signal would have been produced. The second signal illustrates that using the parameter zi supresses transients at the beginning of the output signal:

>>> import numpy as np
>>> from scipy.signal import lfilter, lfilter_zi, butter
...
>>> b, a = butter(5, 0.25)
>>> zi = lfilter_zi(b, a)
>>> y0, zi0 = lfilter(b, a, np.ones(10), zi=zi)
>>> y0
array([1.,  1.,  1.,  1.,  1.,  1.,  1.,  1.,  1.,  1.])
>>> # Another signal:
>>> x = np.array([0.5, 0.5, 0.5, 0.0, 0.0, 0.0, 0.0])
>>> y1, zi1 = lfilter(b, a, x, zi=zi*x[0])
>>> y1
array([ 0.5       ,  0.5       ,  0.5       ,  0.49836039,  0.48610528,
    0.44399389,  0.35505241])

Note that the zi argument to lfilter is computed using lfilter_zi and scaled by x[0]. As a result, the output y1 has no transient until the input drops from 0.5 to 0.