signal.filter_design

Filter design.

Module Contents

Classes

BadCoefficients() Warning about badly conditioned filter coefficients

Functions

findfreqs(num,den,N,kind=”ba”) Find array of frequencies for computing the response of an analog filter.
freqs(b,a,worN=None,plot=None) Compute frequency response of analog filter.
freqs_zpk(z,p,k,worN=None) Compute frequency response of analog filter.
freqz(b,a=1,worN=None,whole=False,plot=None) Compute the frequency response of a digital filter.
freqz_zpk(z,p,k,worN=None,whole=False) r
group_delay(system,w=None,whole=False) rCompute the group delay of a digital filter.
_validate_sos(sos) Helper to validate a SOS input
sosfreqz(sos,worN=None,whole=False) Compute the frequency response of a digital filter in SOS format.
_cplxreal(z,tol=None) Split into complex and real parts, combining conjugate pairs.
_cplxpair(z,tol=None) Sort into pairs of complex conjugates.
tf2zpk(b,a) rReturn zero, pole, gain (z, p, k) representation from a numerator,
zpk2tf(z,p,k) Return polynomial transfer function representation from zeros and poles
tf2sos(b,a,pairing=”nearest”) Return second-order sections from transfer function representation
sos2tf(sos) Return a single transfer function from a series of second-order sections
sos2zpk(sos) Return zeros, poles, and gain of a series of second-order sections
_nearest_real_complex_idx(fro,to,which) Get the next closest real or complex element based on distance
zpk2sos(z,p,k,pairing=”nearest”) Return second-order sections from zeros, poles, and gain of a system
_align_nums(nums) Aligns the shapes of multiple numerators.
normalize(b,a) Normalize numerator/denominator of a continuous-time transfer function.
lp2lp(b,a,wo=1.0) Transform a lowpass filter prototype to a different frequency.
lp2hp(b,a,wo=1.0) Transform a lowpass filter prototype to a highpass filter.
lp2bp(b,a,wo=1.0,bw=1.0) Transform a lowpass filter prototype to a bandpass filter.
lp2bs(b,a,wo=1.0,bw=1.0) Transform a lowpass filter prototype to a bandstop filter.
bilinear(b,a,fs=1.0) Return a digital filter from an analog one using a bilinear transform.
iirdesign(wp,ws,gpass,gstop,analog=False,ftype=”ellip”,output=”ba”) Complete IIR digital and analog filter design.
iirfilter(N,Wn,rp=None,rs=None,btype=”band”,analog=False,ftype=”butter”,output=”ba”) IIR digital and analog filter design given order and critical points.
_relative_degree(z,p) Return relative degree of transfer function from zeros and poles
_zpkbilinear(z,p,k,fs) Return a digital filter from an analog one using a bilinear transform.
_zpklp2lp(z,p,k,wo=1.0) r
_zpklp2hp(z,p,k,wo=1.0) r
_zpklp2bp(z,p,k,wo=1.0,bw=1.0) r
_zpklp2bs(z,p,k,wo=1.0,bw=1.0) r
butter(N,Wn,btype=”low”,analog=False,output=”ba”) Butterworth digital and analog filter design.
cheby1(N,rp,Wn,btype=”low”,analog=False,output=”ba”) Chebyshev type I digital and analog filter design.
cheby2(N,rs,Wn,btype=”low”,analog=False,output=”ba”) Chebyshev type II digital and analog filter design.
ellip(N,rp,rs,Wn,btype=”low”,analog=False,output=”ba”) Elliptic (Cauer) digital and analog filter design.
bessel(N,Wn,btype=”low”,analog=False,output=”ba”,norm=”phase”) Bessel/Thomson digital and analog filter design.
maxflat()
yulewalk()
band_stop_obj(wp,ind,passb,stopb,gpass,gstop,type) Band Stop Objective Function for order minimization.
buttord(wp,ws,gpass,gstop,analog=False) Butterworth filter order selection.
cheb1ord(wp,ws,gpass,gstop,analog=False) Chebyshev type I filter order selection.
cheb2ord(wp,ws,gpass,gstop,analog=False) Chebyshev type II filter order selection.
ellipord(wp,ws,gpass,gstop,analog=False) Elliptic (Cauer) filter order selection.
buttap(N) Return (z,p,k) for analog prototype of Nth-order Butterworth filter.
cheb1ap(N,rp) Return (z,p,k) for Nth-order Chebyshev type I analog lowpass filter.
cheb2ap(N,rs) Return (z,p,k) for Nth-order Chebyshev type I analog lowpass filter.
_vratio(u,ineps,mp)
_kratio(m,k_ratio)
ellipap(N,rp,rs) Return (z,p,k) of Nth-order elliptic analog lowpass filter.
_falling_factorial(x,n) r
_bessel_poly(n,reverse=False) Return the coefficients of Bessel polynomial of degree n
_campos_zeros(n) Return approximate zero locations of Bessel polynomials y_n(x) for order
_aberth(f,fp,x0,tol=1e-15,maxiter=50) Given a function f, its first derivative fp, and a set of initial
_bessel_zeros(N) Find zeros of ordinary Bessel polynomial of order N, by root-finding of
_norm_factor(p,k) Numerically find frequency shift to apply to delay-normalized filter such
besselap(N,norm=”phase”) Return (z,p,k) for analog prototype of an Nth-order Bessel filter.
iirnotch(w0,Q) Design second-order IIR notch digital filter.
iirpeak(w0,Q) Design second-order IIR peak (resonant) digital filter.
_design_notch_peak_filter(w0,Q,ftype) Design notch or peak digital filter.
class BadCoefficients

Warning about badly conditioned filter coefficients

findfreqs(num, den, N, kind="ba")

Find array of frequencies for computing the response of an analog filter.

num, den : array_like, 1-D
The polynomial coefficients of the numerator and denominator of the transfer function of the filter or LTI system, where the coefficients are ordered from highest to lowest degree. Or, the roots of the transfer function numerator and denominator (i.e. zeroes and poles).
N : int
The length of the array to be computed.
kind : str {‘ba’, ‘zp’}, optional
Specifies whether the numerator and denominator are specified by their polynomial coefficients (‘ba’), or their roots (‘zp’).
w : (N,) ndarray
A 1-D array of frequencies, logarithmically spaced.

Find a set of nine frequencies that span the “interesting part” of the frequency response for the filter with the transfer function

H(s) = s / (s^2 + 8s + 25)
>>> from scipy import signal
>>> signal.findfreqs([1, 0], [1, 8, 25], N=9)
array([  1.00000000e-02,   3.16227766e-02,   1.00000000e-01,
         3.16227766e-01,   1.00000000e+00,   3.16227766e+00,
         1.00000000e+01,   3.16227766e+01,   1.00000000e+02])
freqs(b, a, worN=None, plot=None)

Compute frequency response of analog filter.

Given the M-order numerator b and N-order denominator a of an analog filter, compute its frequency response:

        b[0]*(jw)**M + b[1]*(jw)**(M-1) + ... + b[M]
H(w) = ----------------------------------------------
        a[0]*(jw)**N + a[1]*(jw)**(N-1) + ... + a[N]
b : array_like
Numerator of a linear filter.
a : array_like
Denominator of a linear filter.
worN : {None, int, array_like}, optional
If None, then compute at 200 frequencies around the interesting parts of the response curve (determined by pole-zero locations). If a single integer, then compute at that many frequencies. Otherwise, compute the response at the angular frequencies (e.g. rad/s) given in worN.
plot : callable, optional
A callable that takes two arguments. If given, the return parameters w and h are passed to plot. Useful for plotting the frequency response inside freqs.
w : ndarray
The angular frequencies at which h was computed.
h : ndarray
The frequency response.

freqz : Compute the frequency response of a digital filter.

Using Matplotlib’s “plot” function as the callable for plot produces unexpected results, this plots the real part of the complex transfer function, not the magnitude. Try lambda w, h: plot(w, abs(h)).

>>> from scipy.signal import freqs, iirfilter
>>> b, a = iirfilter(4, [1, 10], 1, 60, analog=True, ftype='cheby1')
>>> w, h = freqs(b, a, worN=np.logspace(-1, 2, 1000))
>>> import matplotlib.pyplot as plt
>>> plt.semilogx(w, 20 * np.log10(abs(h)))
>>> plt.xlabel('Frequency')
>>> plt.ylabel('Amplitude response [dB]')
>>> plt.grid()
>>> plt.show()
freqs_zpk(z, p, k, worN=None)

Compute frequency response of analog filter.

Given the zeros z, poles p, and gain k of a filter, compute its frequency response:

           (jw-z[0]) * (jw-z[1]) * ... * (jw-z[-1])
H(w) = k * ----------------------------------------
           (jw-p[0]) * (jw-p[1]) * ... * (jw-p[-1])
z : array_like
Zeroes of a linear filter
p : array_like
Poles of a linear filter
k : scalar
Gain of a linear filter
worN : {None, int, array_like}, optional
If None, then compute at 200 frequencies around the interesting parts of the response curve (determined by pole-zero locations). If a single integer, then compute at that many frequencies. Otherwise, compute the response at the angular frequencies (e.g. rad/s) given in worN.
w : ndarray
The angular frequencies at which h was computed.
h : ndarray
The frequency response.

freqs : Compute the frequency response of an analog filter in TF form freqz : Compute the frequency response of a digital filter in TF form freqz_zpk : Compute the frequency response of a digital filter in ZPK form

>>> from scipy.signal import freqs_zpk, iirfilter
>>> z, p, k = iirfilter(4, [1, 10], 1, 60, analog=True, ftype='cheby1',
...                     output='zpk')
>>> w, h = freqs_zpk(z, p, k, worN=np.logspace(-1, 2, 1000))
>>> import matplotlib.pyplot as plt
>>> plt.semilogx(w, 20 * np.log10(abs(h)))
>>> plt.xlabel('Frequency')
>>> plt.ylabel('Amplitude response [dB]')
>>> plt.grid()
>>> plt.show()
freqz(b, a=1, worN=None, whole=False, plot=None)

Compute the frequency response of a digital filter.

Given the M-order numerator b and N-order denominator a of a digital filter, compute its frequency response:

            jw                 -jw              -jwM
   jw    B(e  )    b[0] + b[1]e    + ... + b[M]e
H(e  ) = ------ = -----------------------------------
            jw                 -jw              -jwN
         A(e  )    a[0] + a[1]e    + ... + a[N]e
b : array_like
Numerator of a linear filter. If b has dimension greater than 1, it is assumed that the coefficients are stored in the first dimension, and b.shape[1:], a.shape[1:], and the shape of the frequencies array must be compatible for broadcasting.
a : array_like
Denominator of a linear filter. If b has dimension greater than 1, it is assumed that the coefficients are stored in the first dimension, and b.shape[1:], a.shape[1:], and the shape of the frequencies array must be compatible for broadcasting.
worN : {None, int, array_like}, optional

If None (default), then compute at 512 equally spaced frequencies. If a single integer, then compute at that many frequencies. This is a convenient alternative to:

np.linspace(0, 2*pi if whole else pi, N, endpoint=False)

Using a number that is fast for FFT computations can result in faster computations (see Notes). If an array_like, compute the response at the frequencies given (in radians/sample).

whole : bool, optional
Normally, frequencies are computed from 0 to the Nyquist frequency, pi radians/sample (upper-half of unit-circle). If whole is True, compute frequencies from 0 to 2*pi radians/sample.
plot : callable
A callable that takes two arguments. If given, the return parameters w and h are passed to plot. Useful for plotting the frequency response inside freqz.
w : ndarray
The normalized frequencies at which h was computed, in radians/sample.
h : ndarray
The frequency response, as complex numbers.

freqz_zpk sosfreqz

Using Matplotlib’s matplotlib.pyplot.plot() function as the callable for plot produces unexpected results, as this plots the real part of the complex transfer function, not the magnitude. Try lambda w, h: plot(w, np.abs(h)).

A direct computation via (R)FFT is used to compute the frequency response when the following conditions are met:

  1. An integer value is given for worN.
  2. worN is fast to compute via FFT (i.e., next_fast_len(worN) <scipy.fftpack.next_fast_len> equals worN).
  3. The denominator coefficients are a single value (a.shape[0] == 1).
  4. worN is at least as long as the numerator coefficients (worN >= b.shape[0]).
  5. If b.ndim > 1, then b.shape[-1] == 1.

For long FIR filters, the FFT approach can have lower error and be much faster than the equivalent direct polynomial calculation.

>>> from scipy import signal
>>> b = signal.firwin(80, 0.5, window=('kaiser', 8))
>>> w, h = signal.freqz(b)
>>> import matplotlib.pyplot as plt
>>> fig = plt.figure()
>>> plt.title('Digital filter frequency response')
>>> ax1 = fig.add_subplot(111)
>>> plt.plot(w, 20 * np.log10(abs(h)), 'b')
>>> plt.ylabel('Amplitude [dB]', color='b')
>>> plt.xlabel('Frequency [rad/sample]')
>>> ax2 = ax1.twinx()
>>> angles = np.unwrap(np.angle(h))
>>> plt.plot(w, angles, 'g')
>>> plt.ylabel('Angle (radians)', color='g')
>>> plt.grid()
>>> plt.axis('tight')
>>> plt.show()

Broadcasting Examples

Suppose we have two FIR filters whose coefficients are stored in the rows of an array with shape (2, 25). For this demonstration we’ll use random data:

>>> np.random.seed(42)
>>> b = np.random.rand(2, 25)

To compute the frequency response for these two filters with one call to freqz, we must pass in b.T, because freqz expects the first axis to hold the coefficients. We must then extend the shape with a trivial dimension of length 1 to allow broadcasting with the array of frequencies. That is, we pass in b.T[..., np.newaxis], which has shape (25, 2, 1):

>>> w, h = signal.freqz(b.T[..., np.newaxis], worN=1024)
>>> w.shape
(1024,)
>>> h.shape
(2, 1024)

Now suppose we have two transfer functions, with the same numerator coefficients b = [0.5, 0.5]. The coefficients for the two denominators are stored in the first dimension of the two-dimensional array a:

a = [   1      1  ]
    [ -0.25, -0.5 ]
>>> b = np.array([0.5, 0.5])
>>> a = np.array([[1, 1], [-0.25, -0.5]])

Only a is more than one-dimensional. To make it compatible for broadcasting with the frequencies, we extend it with a trivial dimension in the call to freqz:

>>> w, h = signal.freqz(b, a[..., np.newaxis], worN=1024)
>>> w.shape
(1024,)
>>> h.shape
(2, 1024)
freqz_zpk(z, p, k, worN=None, whole=False)

r Compute the frequency response of a digital filter in ZPK form.

Given the Zeros, Poles and Gain of a digital filter, compute its frequency response:

:math:`H(z)=k \prod_i (z - Z[i]) / \prod_j (z - P[j])`

where is the gain, are the zeros and are the poles.

z : array_like
Zeroes of a linear filter
p : array_like
Poles of a linear filter
k : scalar
Gain of a linear filter
worN : {None, int, array_like}, optional
If single integer (default 512, same as None), then compute at worN frequencies equally spaced around the unit circle. If an array_like, compute the response at the frequencies given (in radians/sample).
whole : bool, optional
Normally, frequencies are computed from 0 to the Nyquist frequency, pi radians/sample (upper-half of unit-circle). If whole is True, compute frequencies from 0 to 2*pi radians/sample.
w : ndarray
The normalized frequencies at which h was computed, in radians/sample.
h : ndarray
The frequency response.

freqs : Compute the frequency response of an analog filter in TF form freqs_zpk : Compute the frequency response of an analog filter in ZPK form freqz : Compute the frequency response of a digital filter in TF form

>>> from scipy import signal
>>> z, p, k = signal.butter(4, 0.2, output='zpk')
>>> w, h = signal.freqz_zpk(z, p, k)
>>> import matplotlib.pyplot as plt
>>> fig = plt.figure()
>>> plt.title('Digital filter frequency response')
>>> ax1 = fig.add_subplot(111)
>>> plt.plot(w, 20 * np.log10(abs(h)), 'b')
>>> plt.ylabel('Amplitude [dB]', color='b')
>>> plt.xlabel('Frequency [rad/sample]')
>>> ax2 = ax1.twinx()
>>> angles = np.unwrap(np.angle(h))
>>> plt.plot(w, angles, 'g')
>>> plt.ylabel('Angle (radians)', color='g')
>>> plt.grid()
>>> plt.axis('tight')
>>> plt.show()
group_delay(system, w=None, whole=False)

rCompute the group delay of a digital filter.

The group delay measures by how many samples amplitude envelopes of various spectral components of a signal are delayed by a filter. It is formally defined as the derivative of continuous (unwrapped) phase:

          d        jw
D(w) = - -- arg H(e)
         dw
system : tuple of array_like (b, a)
Numerator and denominator coefficients of a filter transfer function.
w : {None, int, array-like}, optional
If None (default), then compute at 512 frequencies equally spaced around the unit circle. If a single integer, then compute at that many frequencies. If array, compute the delay at the frequencies given (in radians/sample).
whole : bool, optional
Normally, frequencies are computed from 0 to the Nyquist frequency, pi radians/sample (upper-half of unit-circle). If whole is True, compute frequencies from 0 to 2*pi radians/sample.
w : ndarray
The normalized frequencies at which the group delay was computed, in radians/sample.
gd : ndarray
The group delay.

The similar function in MATLAB is called grpdelay.

If the transfer function has zeros or poles on the unit circle, the group delay at corresponding frequencies is undefined. When such a case arises the warning is raised and the group delay is set to 0 at those frequencies.

For the details of numerical computation of the group delay refer to [1]_.

freqz : Frequency response of a digital filter

[1]Richard G. Lyons, “Understanding Digital Signal Processing, 3rd edition”, p. 830.
>>> from scipy import signal
>>> b, a = signal.iirdesign(0.1, 0.3, 5, 50, ftype='cheby1')
>>> w, gd = signal.group_delay((b, a))
>>> import matplotlib.pyplot as plt
>>> plt.title('Digital filter group delay')
>>> plt.plot(w, gd)
>>> plt.ylabel('Group delay [samples]')
>>> plt.xlabel('Frequency [rad/sample]')
>>> plt.show()
_validate_sos(sos)

Helper to validate a SOS input

sosfreqz(sos, worN=None, whole=False)

Compute the frequency response of a digital filter in SOS format.

Given sos, an array with shape (n, 6) of second order sections of a digital filter, compute the frequency response of the system function:

       B0(z)   B1(z)         B{n-1}(z)
H(z) = ----- * ----- * ... * ---------
       A0(z)   A1(z)         A{n-1}(z)

for z = exp(omega*1j), where B{k}(z) and A{k}(z) are numerator and denominator of the transfer function of the k-th second order section.

sos : array_like
Array of second-order filter coefficients, must have shape (n_sections, 6). Each row corresponds to a second-order section, with the first three columns providing the numerator coefficients and the last three providing the denominator coefficients.
worN : {None, int, array_like}, optional
If None (default), then compute at 512 frequencies equally spaced around the unit circle. If a single integer, then compute at that many frequencies. Using a number that is fast for FFT computations can result in faster computations (see Notes of freqz). If an array_like, compute the response at the frequencies given (in radians/sample; must be 1D).
whole : bool, optional
Normally, frequencies are computed from 0 to the Nyquist frequency, pi radians/sample (upper-half of unit-circle). If whole is True, compute frequencies from 0 to 2*pi radians/sample.
w : ndarray
The normalized frequencies at which h was computed, in radians/sample.
h : ndarray
The frequency response, as complex numbers.

freqz, sosfilt

New in version 0.19.0.

Design a 15th-order bandpass filter in SOS format.

>>> from scipy import signal
>>> sos = signal.ellip(15, 0.5, 60, (0.2, 0.4), btype='bandpass',
...                    output='sos')

Compute the frequency response at 1500 points from DC to Nyquist.

>>> w, h = signal.sosfreqz(sos, worN=1500)

Plot the response.

>>> import matplotlib.pyplot as plt
>>> plt.subplot(2, 1, 1)
>>> db = 20*np.log10(np.abs(h))
>>> plt.plot(w/np.pi, db)
>>> plt.ylim(-75, 5)
>>> plt.grid(True)
>>> plt.yticks([0, -20, -40, -60])
>>> plt.ylabel('Gain [dB]')
>>> plt.title('Frequency Response')
>>> plt.subplot(2, 1, 2)
>>> plt.plot(w/np.pi, np.angle(h))
>>> plt.grid(True)
>>> plt.yticks([-np.pi, -0.5*np.pi, 0, 0.5*np.pi, np.pi],
...            [r'$-\\pi$', r'$-\\pi/2$', '0', r'$\\pi/2$', r'$\\pi$'])
>>> plt.ylabel('Phase [rad]')
>>> plt.xlabel('Normalized frequency (1.0 = Nyquist)')
>>> plt.show()

If the same filter is implemented as a single transfer function, numerical error corrupts the frequency response:

>>> b, a = signal.ellip(15, 0.5, 60, (0.2, 0.4), btype='bandpass',
...                    output='ba')
>>> w, h = signal.freqz(b, a, worN=1500)
>>> plt.subplot(2, 1, 1)
>>> db = 20*np.log10(np.abs(h))
>>> plt.plot(w/np.pi, db)
>>> plt.subplot(2, 1, 2)
>>> plt.plot(w/np.pi, np.angle(h))
>>> plt.show()
_cplxreal(z, tol=None)

Split into complex and real parts, combining conjugate pairs.

The 1D input vector z is split up into its complex (zc) and real (zr) elements. Every complex element must be part of a complex-conjugate pair, which are combined into a single number (with positive imaginary part) in the output. Two complex numbers are considered a conjugate pair if their real and imaginary parts differ in magnitude by less than tol * abs(z).

z : array_like
Vector of complex numbers to be sorted and split
tol : float, optional
Relative tolerance for testing realness and conjugate equality. Default is 100 * spacing(1) of z’s data type (i.e. 2e-14 for float64)
zc : ndarray
Complex elements of z, with each pair represented by a single value having positive imaginary part, sorted first by real part, and then by magnitude of imaginary part. The pairs are averaged when combined to reduce error.
zr : ndarray
Real elements of z (those having imaginary part less than tol times their magnitude), sorted by value.
ValueError
If there are any complex numbers in z for which a conjugate cannot be found.

_cplxpair

>>> a = [4, 3, 1, 2-2j, 2+2j, 2-1j, 2+1j, 2-1j, 2+1j, 1+1j, 1-1j]
>>> zc, zr = _cplxreal(a)
>>> print(zc)
[ 1.+1.j  2.+1.j  2.+1.j  2.+2.j]
>>> print(zr)
[ 1.  3.  4.]
_cplxpair(z, tol=None)

Sort into pairs of complex conjugates.

Complex conjugates in z are sorted by increasing real part. In each pair, the number with negative imaginary part appears first.

If pairs have identical real parts, they are sorted by increasing imaginary magnitude.

Two complex numbers are considered a conjugate pair if their real and imaginary parts differ in magnitude by less than tol * abs(z). The pairs are forced to be exact complex conjugates by averaging the positive and negative values.

Purely real numbers are also sorted, but placed after the complex conjugate pairs. A number is considered real if its imaginary part is smaller than tol times the magnitude of the number.

z : array_like
1-dimensional input array to be sorted.
tol : float, optional
Relative tolerance for testing realness and conjugate equality. Default is 100 * spacing(1) of z’s data type (i.e. 2e-14 for float64)
y : ndarray
Complex conjugate pairs followed by real numbers.
ValueError
If there are any complex numbers in z for which a conjugate cannot be found.

_cplxreal

>>> a = [4, 3, 1, 2-2j, 2+2j, 2-1j, 2+1j, 2-1j, 2+1j, 1+1j, 1-1j]
>>> z = _cplxpair(a)
>>> print(z)
[ 1.-1.j  1.+1.j  2.-1.j  2.+1.j  2.-1.j  2.+1.j  2.-2.j  2.+2.j  1.+0.j
  3.+0.j  4.+0.j]
tf2zpk(b, a)

rReturn zero, pole, gain (z, p, k) representation from a numerator, denominator representation of a linear filter.

b : array_like
Numerator polynomial coefficients.
a : array_like
Denominator polynomial coefficients.
z : ndarray
Zeros of the transfer function.
p : ndarray
Poles of the transfer function.
k : float
System gain.

If some values of b are too close to 0, they are removed. In that case, a BadCoefficients warning is emitted.

The b and a arrays are interpreted as coefficients for positive, descending powers of the transfer function variable. So the inputs and can represent an analog filter of the form:

or a discrete-time filter of the form:

This “positive powers” form is found more commonly in controls engineering. If M and N are equal (which is true for all filters generated by the bilinear transform), then this happens to be equivalent to the “negative powers” discrete-time form preferred in DSP:

Although this is true for common filters, remember that this is not true in the general case. If M and N are not equal, the discrete-time transfer function coefficients must first be converted to the “positive powers” form before finding the poles and zeros.

zpk2tf(z, p, k)

Return polynomial transfer function representation from zeros and poles

z : array_like
Zeros of the transfer function.
p : array_like
Poles of the transfer function.
k : float
System gain.
b : ndarray
Numerator polynomial coefficients.
a : ndarray
Denominator polynomial coefficients.
tf2sos(b, a, pairing="nearest")

Return second-order sections from transfer function representation

b : array_like
Numerator polynomial coefficients.
a : array_like
Denominator polynomial coefficients.
pairing : {‘nearest’, ‘keep_odd’}, optional
The method to use to combine pairs of poles and zeros into sections. See zpk2sos.
sos : ndarray
Array of second-order filter coefficients, with shape (n_sections, 6). See sosfilt for the SOS filter format specification.

zpk2sos, sosfilt

It is generally discouraged to convert from TF to SOS format, since doing so usually will not improve numerical precision errors. Instead, consider designing filters in ZPK format and converting directly to SOS. TF is converted to SOS by first converting to ZPK format, then converting ZPK to SOS.

New in version 0.16.0.

sos2tf(sos)

Return a single transfer function from a series of second-order sections

sos : array_like
Array of second-order filter coefficients, must have shape (n_sections, 6). See sosfilt for the SOS filter format specification.
b : ndarray
Numerator polynomial coefficients.
a : ndarray
Denominator polynomial coefficients.

New in version 0.16.0.

sos2zpk(sos)

Return zeros, poles, and gain of a series of second-order sections

sos : array_like
Array of second-order filter coefficients, must have shape (n_sections, 6). See sosfilt for the SOS filter format specification.
z : ndarray
Zeros of the transfer function.
p : ndarray
Poles of the transfer function.
k : float
System gain.

New in version 0.16.0.

_nearest_real_complex_idx(fro, to, which)

Get the next closest real or complex element based on distance

zpk2sos(z, p, k, pairing="nearest")

Return second-order sections from zeros, poles, and gain of a system

z : array_like
Zeros of the transfer function.
p : array_like
Poles of the transfer function.
k : float
System gain.
pairing : {‘nearest’, ‘keep_odd’}, optional
The method to use to combine pairs of poles and zeros into sections. See Notes below.
sos : ndarray
Array of second-order filter coefficients, with shape (n_sections, 6). See sosfilt for the SOS filter format specification.

sosfilt

The algorithm used to convert ZPK to SOS format is designed to minimize errors due to numerical precision issues. The pairing algorithm attempts to minimize the peak gain of each biquadratic section. This is done by pairing poles with the nearest zeros, starting with the poles closest to the unit circle.

Algorithms

The current algorithms are designed specifically for use with digital filters. (The output coefficents are not correct for analog filters.)

The steps in the pairing='nearest' and pairing='keep_odd' algorithms are mostly shared. The nearest algorithm attempts to minimize the peak gain, while 'keep_odd' minimizes peak gain under the constraint that odd-order systems should retain one section as first order. The algorithm steps and are as follows:

As a pre-processing step, add poles or zeros to the origin as necessary to obtain the same number of poles and zeros for pairing. If pairing == 'nearest' and there are an odd number of poles, add an additional pole and a zero at the origin.

The following steps are then iterated over until no more poles or zeros remain:

  1. Take the (next remaining) pole (complex or real) closest to the unit circle to begin a new filter section.

  2. If the pole is real and there are no other remaining real poles [7], add the closest real zero to the section and leave it as a first order section. Note that after this step we are guaranteed to be left with an even number of real poles, complex poles, real zeros, and complex zeros for subsequent pairing iterations.

  3. Else:

    1. If the pole is complex and the zero is the only remaining real zero*, then pair the pole with the next closest zero (guaranteed to be complex). This is necessary to ensure that there will be a real zero remaining to eventually create a first-order section (thus keeping the odd order).

    2. Else pair the pole with the closest remaining zero (complex or real).

    3. Proceed to complete the second-order section by adding another pole and zero to the current pole and zero in the section:

      1. If the current pole and zero are both complex, add their conjugates.
      2. Else if the pole is complex and the zero is real, add the conjugate pole and the next closest real zero.
      3. Else if the pole is real and the zero is complex, add the conjugate zero and the real pole closest to those zeros.
      4. Else (we must have a real pole and real zero) add the next real pole closest to the unit circle, and then add the real zero closest to that pole.
[7]This conditional can only be met for specific odd-order inputs with the pairing == 'keep_odd' method.

New in version 0.16.0.

Design a 6th order low-pass elliptic digital filter for a system with a sampling rate of 8000 Hz that has a pass-band corner frequency of 1000 Hz. The ripple in the pass-band should not exceed 0.087 dB, and the attenuation in the stop-band should be at least 90 dB.

In the following call to signal.ellip, we could use output='sos', but for this example, we’ll use output='zpk', and then convert to SOS format with zpk2sos:

>>> from scipy import signal
>>> z, p, k = signal.ellip(6, 0.087, 90, 1000/(0.5*8000), output='zpk')

Now convert to SOS format.

>>> sos = signal.zpk2sos(z, p, k)

The coefficients of the numerators of the sections:

>>> sos[:, :3]
array([[ 0.0014154 ,  0.00248707,  0.0014154 ],
       [ 1.        ,  0.72965193,  1.        ],
       [ 1.        ,  0.17594966,  1.        ]])

The symmetry in the coefficients occurs because all the zeros are on the unit circle.

The coefficients of the denominators of the sections:

>>> sos[:, 3:]
array([[ 1.        , -1.32543251,  0.46989499],
       [ 1.        , -1.26117915,  0.6262586 ],
       [ 1.        , -1.25707217,  0.86199667]])

The next example shows the effect of the pairing option. We have a system with three poles and three zeros, so the SOS array will have shape (2, 6). The means there is, in effect, an extra pole and an extra zero at the origin in the SOS representation.

>>> z1 = np.array([-1, -0.5-0.5j, -0.5+0.5j])
>>> p1 = np.array([0.75, 0.8+0.1j, 0.8-0.1j])

With pairing='nearest' (the default), we obtain

>>> signal.zpk2sos(z1, p1, 1)
array([[ 1.  ,  1.  ,  0.5 ,  1.  , -0.75,  0.  ],
       [ 1.  ,  1.  ,  0.  ,  1.  , -1.6 ,  0.65]])

The first section has the zeros {-0.5-0.05j, -0.5+0.5j} and the poles {0, 0.75}, and the second section has the zeros {-1, 0} and poles {0.8+0.1j, 0.8-0.1j}. Note that the extra pole and zero at the origin have been assigned to different sections.

With pairing='keep_odd', we obtain:

>>> signal.zpk2sos(z1, p1, 1, pairing='keep_odd')
array([[ 1.  ,  1.  ,  0.  ,  1.  , -0.75,  0.  ],
       [ 1.  ,  1.  ,  0.5 ,  1.  , -1.6 ,  0.65]])

The extra pole and zero at the origin are in the same section. The first section is, in effect, a first-order section.

_align_nums(nums)

Aligns the shapes of multiple numerators.

Given an array of numerator coefficient arrays [[a_1, a_2,…, a_n],…, [b_1, b_2,…, b_m]], this function pads shorter numerator arrays with zero’s so that all numerators have the same length. Such alignment is necessary for functions like ‘tf2ss’, which needs the alignment when dealing with SIMO transfer functions.

nums: array_like
Numerator or list of numerators. Not necessarily with same length.
nums: array
The numerator. If nums input was a list of numerators then a 2d array with padded zeros for shorter numerators is returned. Otherwise returns np.asarray(nums).
normalize(b, a)

Normalize numerator/denominator of a continuous-time transfer function.

If values of b are too close to 0, they are removed. In that case, a BadCoefficients warning is emitted.

b: array_like
Numerator of the transfer function. Can be a 2d array to normalize multiple transfer functions.
a: array_like
Denominator of the transfer function. At most 1d.
num: array
The numerator of the normalized transfer function. At least a 1d array. A 2d-array if the input num is a 2d array.
den: 1d-array
The denominator of the normalized transfer function.

Coefficients for both the numerator and denominator should be specified in descending exponent order (e.g., s^2 + 3s + 5 would be represented as [1, 3, 5]).

lp2lp(b, a, wo=1.0)

Transform a lowpass filter prototype to a different frequency.

Return an analog low-pass filter with cutoff frequency wo from an analog low-pass filter prototype with unity cutoff frequency, in transfer function (‘ba’) representation.

lp2hp(b, a, wo=1.0)

Transform a lowpass filter prototype to a highpass filter.

Return an analog high-pass filter with cutoff frequency wo from an analog low-pass filter prototype with unity cutoff frequency, in transfer function (‘ba’) representation.

lp2bp(b, a, wo=1.0, bw=1.0)

Transform a lowpass filter prototype to a bandpass filter.

Return an analog band-pass filter with center frequency wo and bandwidth bw from an analog low-pass filter prototype with unity cutoff frequency, in transfer function (‘ba’) representation.

lp2bs(b, a, wo=1.0, bw=1.0)

Transform a lowpass filter prototype to a bandstop filter.

Return an analog band-stop filter with center frequency wo and bandwidth bw from an analog low-pass filter prototype with unity cutoff frequency, in transfer function (‘ba’) representation.

bilinear(b, a, fs=1.0)

Return a digital filter from an analog one using a bilinear transform.

The bilinear transform substitutes (z-1) / (z+1) for s.

iirdesign(wp, ws, gpass, gstop, analog=False, ftype="ellip", output="ba")

Complete IIR digital and analog filter design.

Given passband and stopband frequencies and gains, construct an analog or digital IIR filter of minimum order for a given basic type. Return the output in numerator, denominator (‘ba’), pole-zero (‘zpk’) or second order sections (‘sos’) form.

wp, ws : float

Passband and stopband edge frequencies. For digital filters, these are normalized from 0 to 1, where 1 is the Nyquist frequency, pi radians/sample. (wp and ws are thus in half-cycles / sample.) For example:

  • Lowpass: wp = 0.2, ws = 0.3
  • Highpass: wp = 0.3, ws = 0.2
  • Bandpass: wp = [0.2, 0.5], ws = [0.1, 0.6]
  • Bandstop: wp = [0.1, 0.6], ws = [0.2, 0.5]

For analog filters, wp and ws are angular frequencies (e.g. rad/s).

gpass : float
The maximum loss in the passband (dB).
gstop : float
The minimum attenuation in the stopband (dB).
analog : bool, optional
When True, return an analog filter, otherwise a digital filter is returned.
ftype : str, optional

The type of IIR filter to design:

  • Butterworth : ‘butter’
  • Chebyshev I : ‘cheby1’
  • Chebyshev II : ‘cheby2’
  • Cauer/elliptic: ‘ellip’
  • Bessel/Thomson: ‘bessel’
output : {‘ba’, ‘zpk’, ‘sos’}, optional
Type of output: numerator/denominator (‘ba’), pole-zero (‘zpk’), or second-order sections (‘sos’). Default is ‘ba’.
b, a : ndarray, ndarray
Numerator (b) and denominator (a) polynomials of the IIR filter. Only returned if output='ba'.
z, p, k : ndarray, ndarray, float
Zeros, poles, and system gain of the IIR filter transfer function. Only returned if output='zpk'.
sos : ndarray
Second-order sections representation of the IIR filter. Only returned if output=='sos'.

butter : Filter design using order and critical points cheby1, cheby2, ellip, bessel buttord : Find order and critical points from passband and stopband spec cheb1ord, cheb2ord, ellipord iirfilter : General filter design using order and critical frequencies

The 'sos' output parameter was added in 0.16.0.

iirfilter(N, Wn, rp=None, rs=None, btype="band", analog=False, ftype="butter", output="ba")

IIR digital and analog filter design given order and critical points.

Design an Nth-order digital or analog filter and return the filter coefficients.

N : int
The order of the filter.
Wn : array_like
A scalar or length-2 sequence giving the critical frequencies. For digital filters, Wn is normalized from 0 to 1, where 1 is the Nyquist frequency, pi radians/sample. (Wn is thus in half-cycles / sample.) For analog filters, Wn is an angular frequency (e.g. rad/s).
rp : float, optional
For Chebyshev and elliptic filters, provides the maximum ripple in the passband. (dB)
rs : float, optional
For Chebyshev and elliptic filters, provides the minimum attenuation in the stop band. (dB)
btype : {‘bandpass’, ‘lowpass’, ‘highpass’, ‘bandstop’}, optional
The type of filter. Default is ‘bandpass’.
analog : bool, optional
When True, return an analog filter, otherwise a digital filter is returned.
ftype : str, optional

The type of IIR filter to design:

  • Butterworth : ‘butter’
  • Chebyshev I : ‘cheby1’
  • Chebyshev II : ‘cheby2’
  • Cauer/elliptic: ‘ellip’
  • Bessel/Thomson: ‘bessel’
output : {‘ba’, ‘zpk’, ‘sos’}, optional
Type of output: numerator/denominator (‘ba’), pole-zero (‘zpk’), or second-order sections (‘sos’). Default is ‘ba’.
b, a : ndarray, ndarray
Numerator (b) and denominator (a) polynomials of the IIR filter. Only returned if output='ba'.
z, p, k : ndarray, ndarray, float
Zeros, poles, and system gain of the IIR filter transfer function. Only returned if output='zpk'.
sos : ndarray
Second-order sections representation of the IIR filter. Only returned if output=='sos'.

butter : Filter design using order and critical points cheby1, cheby2, ellip, bessel buttord : Find order and critical points from passband and stopband spec cheb1ord, cheb2ord, ellipord iirdesign : General filter design using passband and stopband spec

The 'sos' output parameter was added in 0.16.0.

Generate a 17th-order Chebyshev II bandpass filter and plot the frequency response:

>>> from scipy import signal
>>> import matplotlib.pyplot as plt
>>> b, a = signal.iirfilter(17, [50, 200], rs=60, btype='band',
...                         analog=True, ftype='cheby2')
>>> w, h = signal.freqs(b, a, 1000)
>>> fig = plt.figure()
>>> ax = fig.add_subplot(111)
>>> ax.semilogx(w, 20 * np.log10(abs(h)))
>>> ax.set_title('Chebyshev Type II bandpass frequency response')
>>> ax.set_xlabel('Frequency [radians / second]')
>>> ax.set_ylabel('Amplitude [dB]')
>>> ax.axis((10, 1000, -100, 10))
>>> ax.grid(which='both', axis='both')
>>> plt.show()
_relative_degree(z, p)

Return relative degree of transfer function from zeros and poles

_zpkbilinear(z, p, k, fs)

Return a digital filter from an analog one using a bilinear transform.

Transform a set of poles and zeros from the analog s-plane to the digital z-plane using Tustin’s method, which substitutes (z-1) / (z+1) for s, maintaining the shape of the frequency response.

z : array_like
Zeros of the analog IIR filter transfer function.
p : array_like
Poles of the analog IIR filter transfer function.
k : float
System gain of the analog IIR filter transfer function.
fs : float
Sample rate, as ordinary frequency (e.g. hertz). No prewarping is done in this function.
z : ndarray
Zeros of the transformed digital filter transfer function.
p : ndarray
Poles of the transformed digital filter transfer function.
k : float
System gain of the transformed digital filter.
_zpklp2lp(z, p, k, wo=1.0)

r Transform a lowpass filter prototype to a different frequency.

Return an analog low-pass filter with cutoff frequency wo from an analog low-pass filter prototype with unity cutoff frequency, using zeros, poles, and gain (‘zpk’) representation.

z : array_like
Zeros of the analog IIR filter transfer function.
p : array_like
Poles of the analog IIR filter transfer function.
k : float
System gain of the analog IIR filter transfer function.
wo : float
Desired cutoff, as angular frequency (e.g. rad/s). Defaults to no change.
z : ndarray
Zeros of the transformed low-pass filter transfer function.
p : ndarray
Poles of the transformed low-pass filter transfer function.
k : float
System gain of the transformed low-pass filter.

This is derived from the s-plane substitution

_zpklp2hp(z, p, k, wo=1.0)

r Transform a lowpass filter prototype to a highpass filter.

Return an analog high-pass filter with cutoff frequency wo from an analog low-pass filter prototype with unity cutoff frequency, using zeros, poles, and gain (‘zpk’) representation.

z : array_like
Zeros of the analog IIR filter transfer function.
p : array_like
Poles of the analog IIR filter transfer function.
k : float
System gain of the analog IIR filter transfer function.
wo : float
Desired cutoff, as angular frequency (e.g. rad/s). Defaults to no change.
z : ndarray
Zeros of the transformed high-pass filter transfer function.
p : ndarray
Poles of the transformed high-pass filter transfer function.
k : float
System gain of the transformed high-pass filter.

This is derived from the s-plane substitution

This maintains symmetry of the lowpass and highpass responses on a logarithmic scale.

_zpklp2bp(z, p, k, wo=1.0, bw=1.0)

r Transform a lowpass filter prototype to a bandpass filter.

Return an analog band-pass filter with center frequency wo and bandwidth bw from an analog low-pass filter prototype with unity cutoff frequency, using zeros, poles, and gain (‘zpk’) representation.

z : array_like
Zeros of the analog IIR filter transfer function.
p : array_like
Poles of the analog IIR filter transfer function.
k : float
System gain of the analog IIR filter transfer function.
wo : float
Desired passband center, as angular frequency (e.g. rad/s). Defaults to no change.
bw : float
Desired passband width, as angular frequency (e.g. rad/s). Defaults to 1.
z : ndarray
Zeros of the transformed band-pass filter transfer function.
p : ndarray
Poles of the transformed band-pass filter transfer function.
k : float
System gain of the transformed band-pass filter.

This is derived from the s-plane substitution

This is the “wideband” transformation, producing a passband with geometric (log frequency) symmetry about wo.

_zpklp2bs(z, p, k, wo=1.0, bw=1.0)

r Transform a lowpass filter prototype to a bandstop filter.

Return an analog band-stop filter with center frequency wo and stopband width bw from an analog low-pass filter prototype with unity cutoff frequency, using zeros, poles, and gain (‘zpk’) representation.

z : array_like
Zeros of the analog IIR filter transfer function.
p : array_like
Poles of the analog IIR filter transfer function.
k : float
System gain of the analog IIR filter transfer function.
wo : float
Desired stopband center, as angular frequency (e.g. rad/s). Defaults to no change.
bw : float
Desired stopband width, as angular frequency (e.g. rad/s). Defaults to 1.
z : ndarray
Zeros of the transformed band-stop filter transfer function.
p : ndarray
Poles of the transformed band-stop filter transfer function.
k : float
System gain of the transformed band-stop filter.

This is derived from the s-plane substitution

This is the “wideband” transformation, producing a stopband with geometric (log frequency) symmetry about wo.

butter(N, Wn, btype="low", analog=False, output="ba")

Butterworth digital and analog filter design.

Design an Nth-order digital or analog Butterworth filter and return the filter coefficients.

N : int
The order of the filter.
Wn : array_like
A scalar or length-2 sequence giving the critical frequencies. For a Butterworth filter, this is the point at which the gain drops to 1/sqrt(2) that of the passband (the “-3 dB point”). For digital filters, Wn is normalized from 0 to 1, where 1 is the Nyquist frequency, pi radians/sample. (Wn is thus in half-cycles / sample.) For analog filters, Wn is an angular frequency (e.g. rad/s).
btype : {‘lowpass’, ‘highpass’, ‘bandpass’, ‘bandstop’}, optional
The type of filter. Default is ‘lowpass’.
analog : bool, optional
When True, return an analog filter, otherwise a digital filter is returned.
output : {‘ba’, ‘zpk’, ‘sos’}, optional
Type of output: numerator/denominator (‘ba’), pole-zero (‘zpk’), or second-order sections (‘sos’). Default is ‘ba’.
b, a : ndarray, ndarray
Numerator (b) and denominator (a) polynomials of the IIR filter. Only returned if output='ba'.
z, p, k : ndarray, ndarray, float
Zeros, poles, and system gain of the IIR filter transfer function. Only returned if output='zpk'.
sos : ndarray
Second-order sections representation of the IIR filter. Only returned if output=='sos'.

buttord, buttap

The Butterworth filter has maximally flat frequency response in the passband.

The 'sos' output parameter was added in 0.16.0.

Plot the filter’s frequency response, showing the critical points:

>>> from scipy import signal
>>> import matplotlib.pyplot as plt
>>> b, a = signal.butter(4, 100, 'low', analog=True)
>>> w, h = signal.freqs(b, a)
>>> plt.semilogx(w, 20 * np.log10(abs(h)))
>>> plt.title('Butterworth filter frequency response')
>>> plt.xlabel('Frequency [radians / second]')
>>> plt.ylabel('Amplitude [dB]')
>>> plt.margins(0, 0.1)
>>> plt.grid(which='both', axis='both')
>>> plt.axvline(100, color='green') # cutoff frequency
>>> plt.show()
cheby1(N, rp, Wn, btype="low", analog=False, output="ba")

Chebyshev type I digital and analog filter design.

Design an Nth-order digital or analog Chebyshev type I filter and return the filter coefficients.

N : int
The order of the filter.
rp : float
The maximum ripple allowed below unity gain in the passband. Specified in decibels, as a positive number.
Wn : array_like
A scalar or length-2 sequence giving the critical frequencies. For Type I filters, this is the point in the transition band at which the gain first drops below -rp. For digital filters, Wn is normalized from 0 to 1, where 1 is the Nyquist frequency, pi radians/sample. (Wn is thus in half-cycles / sample.) For analog filters, Wn is an angular frequency (e.g. rad/s).
btype : {‘lowpass’, ‘highpass’, ‘bandpass’, ‘bandstop’}, optional
The type of filter. Default is ‘lowpass’.
analog : bool, optional
When True, return an analog filter, otherwise a digital filter is returned.
output : {‘ba’, ‘zpk’, ‘sos’}, optional
Type of output: numerator/denominator (‘ba’), pole-zero (‘zpk’), or second-order sections (‘sos’). Default is ‘ba’.
b, a : ndarray, ndarray
Numerator (b) and denominator (a) polynomials of the IIR filter. Only returned if output='ba'.
z, p, k : ndarray, ndarray, float
Zeros, poles, and system gain of the IIR filter transfer function. Only returned if output='zpk'.
sos : ndarray
Second-order sections representation of the IIR filter. Only returned if output=='sos'.

cheb1ord, cheb1ap

The Chebyshev type I filter maximizes the rate of cutoff between the frequency response’s passband and stopband, at the expense of ripple in the passband and increased ringing in the step response.

Type I filters roll off faster than Type II (cheby2), but Type II filters do not have any ripple in the passband.

The equiripple passband has N maxima or minima (for example, a 5th-order filter has 3 maxima and 2 minima). Consequently, the DC gain is unity for odd-order filters, or -rp dB for even-order filters.

The 'sos' output parameter was added in 0.16.0.

Plot the filter’s frequency response, showing the critical points:

>>> from scipy import signal
>>> import matplotlib.pyplot as plt
>>> b, a = signal.cheby1(4, 5, 100, 'low', analog=True)
>>> w, h = signal.freqs(b, a)
>>> plt.semilogx(w, 20 * np.log10(abs(h)))
>>> plt.title('Chebyshev Type I frequency response (rp=5)')
>>> plt.xlabel('Frequency [radians / second]')
>>> plt.ylabel('Amplitude [dB]')
>>> plt.margins(0, 0.1)
>>> plt.grid(which='both', axis='both')
>>> plt.axvline(100, color='green') # cutoff frequency
>>> plt.axhline(-5, color='green') # rp
>>> plt.show()
cheby2(N, rs, Wn, btype="low", analog=False, output="ba")

Chebyshev type II digital and analog filter design.

Design an Nth-order digital or analog Chebyshev type II filter and return the filter coefficients.

N : int
The order of the filter.
rs : float
The minimum attenuation required in the stop band. Specified in decibels, as a positive number.
Wn : array_like
A scalar or length-2 sequence giving the critical frequencies. For Type II filters, this is the point in the transition band at which the gain first reaches -rs. For digital filters, Wn is normalized from 0 to 1, where 1 is the Nyquist frequency, pi radians/sample. (Wn is thus in half-cycles / sample.) For analog filters, Wn is an angular frequency (e.g. rad/s).
btype : {‘lowpass’, ‘highpass’, ‘bandpass’, ‘bandstop’}, optional
The type of filter. Default is ‘lowpass’.
analog : bool, optional
When True, return an analog filter, otherwise a digital filter is returned.
output : {‘ba’, ‘zpk’, ‘sos’}, optional
Type of output: numerator/denominator (‘ba’), pole-zero (‘zpk’), or second-order sections (‘sos’). Default is ‘ba’.
b, a : ndarray, ndarray
Numerator (b) and denominator (a) polynomials of the IIR filter. Only returned if output='ba'.
z, p, k : ndarray, ndarray, float
Zeros, poles, and system gain of the IIR filter transfer function. Only returned if output='zpk'.
sos : ndarray
Second-order sections representation of the IIR filter. Only returned if output=='sos'.

cheb2ord, cheb2ap

The Chebyshev type II filter maximizes the rate of cutoff between the frequency response’s passband and stopband, at the expense of ripple in the stopband and increased ringing in the step response.

Type II filters do not roll off as fast as Type I (cheby1).

The 'sos' output parameter was added in 0.16.0.

Plot the filter’s frequency response, showing the critical points:

>>> from scipy import signal
>>> import matplotlib.pyplot as plt
>>> b, a = signal.cheby2(4, 40, 100, 'low', analog=True)
>>> w, h = signal.freqs(b, a)
>>> plt.semilogx(w, 20 * np.log10(abs(h)))
>>> plt.title('Chebyshev Type II frequency response (rs=40)')
>>> plt.xlabel('Frequency [radians / second]')
>>> plt.ylabel('Amplitude [dB]')
>>> plt.margins(0, 0.1)
>>> plt.grid(which='both', axis='both')
>>> plt.axvline(100, color='green') # cutoff frequency
>>> plt.axhline(-40, color='green') # rs
>>> plt.show()
ellip(N, rp, rs, Wn, btype="low", analog=False, output="ba")

Elliptic (Cauer) digital and analog filter design.

Design an Nth-order digital or analog elliptic filter and return the filter coefficients.

N : int
The order of the filter.
rp : float
The maximum ripple allowed below unity gain in the passband. Specified in decibels, as a positive number.
rs : float
The minimum attenuation required in the stop band. Specified in decibels, as a positive number.
Wn : array_like
A scalar or length-2 sequence giving the critical frequencies. For elliptic filters, this is the point in the transition band at which the gain first drops below -rp. For digital filters, Wn is normalized from 0 to 1, where 1 is the Nyquist frequency, pi radians/sample. (Wn is thus in half-cycles / sample.) For analog filters, Wn is an angular frequency (e.g. rad/s).
btype : {‘lowpass’, ‘highpass’, ‘bandpass’, ‘bandstop’}, optional
The type of filter. Default is ‘lowpass’.
analog : bool, optional
When True, return an analog filter, otherwise a digital filter is returned.
output : {‘ba’, ‘zpk’, ‘sos’}, optional
Type of output: numerator/denominator (‘ba’), pole-zero (‘zpk’), or second-order sections (‘sos’). Default is ‘ba’.
b, a : ndarray, ndarray
Numerator (b) and denominator (a) polynomials of the IIR filter. Only returned if output='ba'.
z, p, k : ndarray, ndarray, float
Zeros, poles, and system gain of the IIR filter transfer function. Only returned if output='zpk'.
sos : ndarray
Second-order sections representation of the IIR filter. Only returned if output=='sos'.

ellipord, ellipap

Also known as Cauer or Zolotarev filters, the elliptical filter maximizes the rate of transition between the frequency response’s passband and stopband, at the expense of ripple in both, and increased ringing in the step response.

As rp approaches 0, the elliptical filter becomes a Chebyshev type II filter (cheby2). As rs approaches 0, it becomes a Chebyshev type I filter (cheby1). As both approach 0, it becomes a Butterworth filter (butter).

The equiripple passband has N maxima or minima (for example, a 5th-order filter has 3 maxima and 2 minima). Consequently, the DC gain is unity for odd-order filters, or -rp dB for even-order filters.

The 'sos' output parameter was added in 0.16.0.

Plot the filter’s frequency response, showing the critical points:

>>> from scipy import signal
>>> import matplotlib.pyplot as plt
>>> b, a = signal.ellip(4, 5, 40, 100, 'low', analog=True)
>>> w, h = signal.freqs(b, a)
>>> plt.semilogx(w, 20 * np.log10(abs(h)))
>>> plt.title('Elliptic filter frequency response (rp=5, rs=40)')
>>> plt.xlabel('Frequency [radians / second]')
>>> plt.ylabel('Amplitude [dB]')
>>> plt.margins(0, 0.1)
>>> plt.grid(which='both', axis='both')
>>> plt.axvline(100, color='green') # cutoff frequency
>>> plt.axhline(-40, color='green') # rs
>>> plt.axhline(-5, color='green') # rp
>>> plt.show()
bessel(N, Wn, btype="low", analog=False, output="ba", norm="phase")

Bessel/Thomson digital and analog filter design.

Design an Nth-order digital or analog Bessel filter and return the filter coefficients.

N : int
The order of the filter.
Wn : array_like
A scalar or length-2 sequence giving the critical frequencies (defined by the norm parameter). For analog filters, Wn is an angular frequency (e.g. rad/s). For digital filters, Wn is normalized from 0 to 1, where 1 is the Nyquist frequency, pi radians/sample. (Wn is thus in half-cycles / sample.)
btype : {‘lowpass’, ‘highpass’, ‘bandpass’, ‘bandstop’}, optional
The type of filter. Default is ‘lowpass’.
analog : bool, optional
When True, return an analog filter, otherwise a digital filter is returned. (See Notes.)
output : {‘ba’, ‘zpk’, ‘sos’}, optional
Type of output: numerator/denominator (‘ba’), pole-zero (‘zpk’), or second-order sections (‘sos’). Default is ‘ba’.
norm : {‘phase’, ‘delay’, ‘mag’}, optional

Critical frequency normalization:

phase

The filter is normalized such that the phase response reaches its midpoint at angular (e.g. rad/s) frequency Wn. This happens for both low-pass and high-pass filters, so this is the “phase-matched” case.

The magnitude response asymptotes are the same as a Butterworth filter of the same order with a cutoff of Wn.

This is the default, and matches MATLAB’s implementation.

delay
The filter is normalized such that the group delay in the passband is 1/Wn (e.g. seconds). This is the “natural” type obtained by solving Bessel polynomials.
mag
The filter is normalized such that the gain magnitude is -3 dB at angular frequency Wn.

New in version 0.18.0.

b, a : ndarray, ndarray
Numerator (b) and denominator (a) polynomials of the IIR filter. Only returned if output='ba'.
z, p, k : ndarray, ndarray, float
Zeros, poles, and system gain of the IIR filter transfer function. Only returned if output='zpk'.
sos : ndarray
Second-order sections representation of the IIR filter. Only returned if output=='sos'.

Also known as a Thomson filter, the analog Bessel filter has maximally flat group delay and maximally linear phase response, with very little ringing in the step response. [1]_

The Bessel is inherently an analog filter. This function generates digital Bessel filters using the bilinear transform, which does not preserve the phase response of the analog filter. As such, it is only approximately correct at frequencies below about fs/4. To get maximally-flat group delay at higher frequencies, the analog Bessel filter must be transformed using phase-preserving techniques.

See besselap for implementation details and references.

The 'sos' output parameter was added in 0.16.0.

Plot the phase-normalized frequency response, showing the relationship to the Butterworth’s cutoff frequency (green):

>>> from scipy import signal
>>> import matplotlib.pyplot as plt
>>> b, a = signal.butter(4, 100, 'low', analog=True)
>>> w, h = signal.freqs(b, a)
>>> plt.semilogx(w, 20 * np.log10(np.abs(h)), color='silver', ls='dashed')
>>> b, a = signal.bessel(4, 100, 'low', analog=True, norm='phase')
>>> w, h = signal.freqs(b, a)
>>> plt.semilogx(w, 20 * np.log10(np.abs(h)))
>>> plt.title('Bessel filter magnitude response (with Butterworth)')
>>> plt.xlabel('Frequency [radians / second]')
>>> plt.ylabel('Amplitude [dB]')
>>> plt.margins(0, 0.1)
>>> plt.grid(which='both', axis='both')
>>> plt.axvline(100, color='green')  # cutoff frequency
>>> plt.show()

and the phase midpoint:

>>> plt.figure()
>>> plt.semilogx(w, np.unwrap(np.angle(h)))
>>> plt.axvline(100, color='green')  # cutoff frequency
>>> plt.axhline(-np.pi, color='red')  # phase midpoint
>>> plt.title('Bessel filter phase response')
>>> plt.xlabel('Frequency [radians / second]')
>>> plt.ylabel('Phase [radians]')
>>> plt.margins(0, 0.1)
>>> plt.grid(which='both', axis='both')
>>> plt.show()

Plot the magnitude-normalized frequency response, showing the -3 dB cutoff:

>>> b, a = signal.bessel(3, 10, 'low', analog=True, norm='mag')
>>> w, h = signal.freqs(b, a)
>>> plt.semilogx(w, 20 * np.log10(np.abs(h)))
>>> plt.axhline(-3, color='red')  # -3 dB magnitude
>>> plt.axvline(10, color='green')  # cutoff frequency
>>> plt.title('Magnitude-normalized Bessel filter frequency response')
>>> plt.xlabel('Frequency [radians / second]')
>>> plt.ylabel('Amplitude [dB]')
>>> plt.margins(0, 0.1)
>>> plt.grid(which='both', axis='both')
>>> plt.show()

Plot the delay-normalized filter, showing the maximally-flat group delay at 0.1 seconds:

>>> b, a = signal.bessel(5, 1/0.1, 'low', analog=True, norm='delay')
>>> w, h = signal.freqs(b, a)
>>> plt.figure()
>>> plt.semilogx(w[1:], -np.diff(np.unwrap(np.angle(h)))/np.diff(w))
>>> plt.axhline(0.1, color='red')  # 0.1 seconds group delay
>>> plt.title('Bessel filter group delay')
>>> plt.xlabel('Frequency [radians / second]')
>>> plt.ylabel('Group delay [seconds]')
>>> plt.margins(0, 0.1)
>>> plt.grid(which='both', axis='both')
>>> plt.show()
[1]Thomson, W.E., “Delay Networks having Maximally Flat Frequency Characteristics”, Proceedings of the Institution of Electrical Engineers, Part III, November 1949, Vol. 96, No. 44, pp. 487-490.
maxflat()
yulewalk()
band_stop_obj(wp, ind, passb, stopb, gpass, gstop, type)

Band Stop Objective Function for order minimization.

Returns the non-integer order for an analog band stop filter.

wp : scalar
Edge of passband passb.
ind : int, {0, 1}
Index specifying which passb edge to vary (0 or 1).
passb : ndarray
Two element sequence of fixed passband edges.
stopb : ndarray
Two element sequence of fixed stopband edges.
gstop : float
Amount of attenuation in stopband in dB.
gpass : float
Amount of ripple in the passband in dB.
type : {‘butter’, ‘cheby’, ‘ellip’}
Type of filter.
n : scalar
Filter order (possibly non-integer).
buttord(wp, ws, gpass, gstop, analog=False)

Butterworth filter order selection.

Return the order of the lowest order digital or analog Butterworth filter that loses no more than gpass dB in the passband and has at least gstop dB attenuation in the stopband.

wp, ws : float

Passband and stopband edge frequencies. For digital filters, these are normalized from 0 to 1, where 1 is the Nyquist frequency, pi radians/sample. (wp and ws are thus in half-cycles / sample.) For example:

  • Lowpass: wp = 0.2, ws = 0.3
  • Highpass: wp = 0.3, ws = 0.2
  • Bandpass: wp = [0.2, 0.5], ws = [0.1, 0.6]
  • Bandstop: wp = [0.1, 0.6], ws = [0.2, 0.5]

For analog filters, wp and ws are angular frequencies (e.g. rad/s).

gpass : float
The maximum loss in the passband (dB).
gstop : float
The minimum attenuation in the stopband (dB).
analog : bool, optional
When True, return an analog filter, otherwise a digital filter is returned.
ord : int
The lowest order for a Butterworth filter which meets specs.
wn : ndarray or float
The Butterworth natural frequency (i.e. the “3dB frequency”). Should be used with butter to give filter results.

butter : Filter design using order and critical points cheb1ord : Find order and critical points from passband and stopband spec cheb2ord, ellipord iirfilter : General filter design using order and critical frequencies iirdesign : General filter design using passband and stopband spec

Design an analog bandpass filter with passband within 3 dB from 20 to 50 rad/s, while rejecting at least -40 dB below 14 and above 60 rad/s. Plot its frequency response, showing the passband and stopband constraints in gray.

>>> from scipy import signal
>>> import matplotlib.pyplot as plt
>>> N, Wn = signal.buttord([20, 50], [14, 60], 3, 40, True)
>>> b, a = signal.butter(N, Wn, 'band', True)
>>> w, h = signal.freqs(b, a, np.logspace(1, 2, 500))
>>> plt.semilogx(w, 20 * np.log10(abs(h)))
>>> plt.title('Butterworth bandpass filter fit to constraints')
>>> plt.xlabel('Frequency [radians / second]')
>>> plt.ylabel('Amplitude [dB]')
>>> plt.grid(which='both', axis='both')
>>> plt.fill([1,  14,  14,   1], [-40, -40, 99, 99], '0.9', lw=0) # stop
>>> plt.fill([20, 20,  50,  50], [-99, -3, -3, -99], '0.9', lw=0) # pass
>>> plt.fill([60, 60, 1e9, 1e9], [99, -40, -40, 99], '0.9', lw=0) # stop
>>> plt.axis([10, 100, -60, 3])
>>> plt.show()
cheb1ord(wp, ws, gpass, gstop, analog=False)

Chebyshev type I filter order selection.

Return the order of the lowest order digital or analog Chebyshev Type I filter that loses no more than gpass dB in the passband and has at least gstop dB attenuation in the stopband.

wp, ws : float

Passband and stopband edge frequencies. For digital filters, these are normalized from 0 to 1, where 1 is the Nyquist frequency, pi radians/sample. (wp and ws are thus in half-cycles / sample.) For example:

  • Lowpass: wp = 0.2, ws = 0.3
  • Highpass: wp = 0.3, ws = 0.2
  • Bandpass: wp = [0.2, 0.5], ws = [0.1, 0.6]
  • Bandstop: wp = [0.1, 0.6], ws = [0.2, 0.5]

For analog filters, wp and ws are angular frequencies (e.g. rad/s).

gpass : float
The maximum loss in the passband (dB).
gstop : float
The minimum attenuation in the stopband (dB).
analog : bool, optional
When True, return an analog filter, otherwise a digital filter is returned.
ord : int
The lowest order for a Chebyshev type I filter that meets specs.
wn : ndarray or float
The Chebyshev natural frequency (the “3dB frequency”) for use with cheby1 to give filter results.

cheby1 : Filter design using order and critical points buttord : Find order and critical points from passband and stopband spec cheb2ord, ellipord iirfilter : General filter design using order and critical frequencies iirdesign : General filter design using passband and stopband spec

Design a digital lowpass filter such that the passband is within 3 dB up to 0.2*(fs/2), while rejecting at least -40 dB above 0.3*(fs/2). Plot its frequency response, showing the passband and stopband constraints in gray.

>>> from scipy import signal
>>> import matplotlib.pyplot as plt
>>> N, Wn = signal.cheb1ord(0.2, 0.3, 3, 40)
>>> b, a = signal.cheby1(N, 3, Wn, 'low')
>>> w, h = signal.freqz(b, a)
>>> plt.semilogx(w / np.pi, 20 * np.log10(abs(h)))
>>> plt.title('Chebyshev I lowpass filter fit to constraints')
>>> plt.xlabel('Normalized frequency')
>>> plt.ylabel('Amplitude [dB]')
>>> plt.grid(which='both', axis='both')
>>> plt.fill([.01, 0.2, 0.2, .01], [-3, -3, -99, -99], '0.9', lw=0) # stop
>>> plt.fill([0.3, 0.3,   2,   2], [ 9, -40, -40,  9], '0.9', lw=0) # pass
>>> plt.axis([0.08, 1, -60, 3])
>>> plt.show()
cheb2ord(wp, ws, gpass, gstop, analog=False)

Chebyshev type II filter order selection.

Return the order of the lowest order digital or analog Chebyshev Type II filter that loses no more than gpass dB in the passband and has at least gstop dB attenuation in the stopband.

wp, ws : float

Passband and stopband edge frequencies. For digital filters, these are normalized from 0 to 1, where 1 is the Nyquist frequency, pi radians/sample. (wp and ws are thus in half-cycles / sample.) For example:

  • Lowpass: wp = 0.2, ws = 0.3
  • Highpass: wp = 0.3, ws = 0.2
  • Bandpass: wp = [0.2, 0.5], ws = [0.1, 0.6]
  • Bandstop: wp = [0.1, 0.6], ws = [0.2, 0.5]

For analog filters, wp and ws are angular frequencies (e.g. rad/s).

gpass : float
The maximum loss in the passband (dB).
gstop : float
The minimum attenuation in the stopband (dB).
analog : bool, optional
When True, return an analog filter, otherwise a digital filter is returned.
ord : int
The lowest order for a Chebyshev type II filter that meets specs.
wn : ndarray or float
The Chebyshev natural frequency (the “3dB frequency”) for use with cheby2 to give filter results.

cheby2 : Filter design using order and critical points buttord : Find order and critical points from passband and stopband spec cheb1ord, ellipord iirfilter : General filter design using order and critical frequencies iirdesign : General filter design using passband and stopband spec

Design a digital bandstop filter which rejects -60 dB from 0.2*(fs/2) to 0.5*(fs/2), while staying within 3 dB below 0.1*(fs/2) or above 0.6*(fs/2). Plot its frequency response, showing the passband and stopband constraints in gray.

>>> from scipy import signal
>>> import matplotlib.pyplot as plt
>>> N, Wn = signal.cheb2ord([0.1, 0.6], [0.2, 0.5], 3, 60)
>>> b, a = signal.cheby2(N, 60, Wn, 'stop')
>>> w, h = signal.freqz(b, a)
>>> plt.semilogx(w / np.pi, 20 * np.log10(abs(h)))
>>> plt.title('Chebyshev II bandstop filter fit to constraints')
>>> plt.xlabel('Normalized frequency')
>>> plt.ylabel('Amplitude [dB]')
>>> plt.grid(which='both', axis='both')
>>> plt.fill([.01, .1, .1, .01], [-3,  -3, -99, -99], '0.9', lw=0) # stop
>>> plt.fill([.2,  .2, .5,  .5], [ 9, -60, -60,   9], '0.9', lw=0) # pass
>>> plt.fill([.6,  .6,  2,   2], [-99, -3,  -3, -99], '0.9', lw=0) # stop
>>> plt.axis([0.06, 1, -80, 3])
>>> plt.show()
ellipord(wp, ws, gpass, gstop, analog=False)

Elliptic (Cauer) filter order selection.

Return the order of the lowest order digital or analog elliptic filter that loses no more than gpass dB in the passband and has at least gstop dB attenuation in the stopband.

wp, ws : float

Passband and stopband edge frequencies. For digital filters, these are normalized from 0 to 1, where 1 is the Nyquist frequency, pi radians/sample. (wp and ws are thus in half-cycles / sample.) For example:

  • Lowpass: wp = 0.2, ws = 0.3
  • Highpass: wp = 0.3, ws = 0.2
  • Bandpass: wp = [0.2, 0.5], ws = [0.1, 0.6]
  • Bandstop: wp = [0.1, 0.6], ws = [0.2, 0.5]

For analog filters, wp and ws are angular frequencies (e.g. rad/s).

gpass : float
The maximum loss in the passband (dB).
gstop : float
The minimum attenuation in the stopband (dB).
analog : bool, optional
When True, return an analog filter, otherwise a digital filter is returned.
ord : int
The lowest order for an Elliptic (Cauer) filter that meets specs.
wn : ndarray or float
The Chebyshev natural frequency (the “3dB frequency”) for use with ellip to give filter results.

ellip : Filter design using order and critical points buttord : Find order and critical points from passband and stopband spec cheb1ord, cheb2ord iirfilter : General filter design using order and critical frequencies iirdesign : General filter design using passband and stopband spec

Design an analog highpass filter such that the passband is within 3 dB above 30 rad/s, while rejecting -60 dB at 10 rad/s. Plot its frequency response, showing the passband and stopband constraints in gray.

>>> from scipy import signal
>>> import matplotlib.pyplot as plt
>>> N, Wn = signal.ellipord(30, 10, 3, 60, True)
>>> b, a = signal.ellip(N, 3, 60, Wn, 'high', True)
>>> w, h = signal.freqs(b, a, np.logspace(0, 3, 500))
>>> plt.semilogx(w, 20 * np.log10(abs(h)))
>>> plt.title('Elliptical highpass filter fit to constraints')
>>> plt.xlabel('Frequency [radians / second]')
>>> plt.ylabel('Amplitude [dB]')
>>> plt.grid(which='both', axis='both')
>>> plt.fill([.1, 10,  10,  .1], [1e4, 1e4, -60, -60], '0.9', lw=0) # stop
>>> plt.fill([30, 30, 1e9, 1e9], [-99,  -3,  -3, -99], '0.9', lw=0) # pass
>>> plt.axis([1, 300, -80, 3])
>>> plt.show()
buttap(N)

Return (z,p,k) for analog prototype of Nth-order Butterworth filter.

The filter will have an angular (e.g. rad/s) cutoff frequency of 1.

butter : Filter design function using this prototype

cheb1ap(N, rp)

Return (z,p,k) for Nth-order Chebyshev type I analog lowpass filter.

The returned filter prototype has rp decibels of ripple in the passband.

The filter’s angular (e.g. rad/s) cutoff frequency is normalized to 1, defined as the point at which the gain first drops below -rp.

cheby1 : Filter design function using this prototype

cheb2ap(N, rs)

Return (z,p,k) for Nth-order Chebyshev type I analog lowpass filter.

The returned filter prototype has rs decibels of ripple in the stopband.

The filter’s angular (e.g. rad/s) cutoff frequency is normalized to 1, defined as the point at which the gain first reaches -rs.

cheby2 : Filter design function using this prototype

_vratio(u, ineps, mp)
_kratio(m, k_ratio)
ellipap(N, rp, rs)

Return (z,p,k) of Nth-order elliptic analog lowpass filter.

The filter is a normalized prototype that has rp decibels of ripple in the passband and a stopband rs decibels down.

The filter’s angular (e.g. rad/s) cutoff frequency is normalized to 1, defined as the point at which the gain first drops below -rp.

ellip : Filter design function using this prototype

[1]Lutova, Tosic, and Evans, “Filter Design for Signal Processing”, Chapters 5 and 12.
_falling_factorial(x, n)

r Return the factorial of x to the n falling.

This is defined as:

This can more efficiently calculate ratios of factorials, since:

n!/m! == falling_factorial(n, n-m)

where n >= m

skipping the factors that cancel out

the usual factorial n! == ff(n, n)

_bessel_poly(n, reverse=False)

Return the coefficients of Bessel polynomial of degree n

If reverse is true, a reverse Bessel polynomial is output.

Output is a list of coefficients: [1] = 1 [1, 1] = 1*s + 1 [1, 3, 3] = 1*s^2 + 3*s + 3 [1, 6, 15, 15] = 1*s^3 + 6*s^2 + 15*s + 15 [1, 10, 45, 105, 105] = 1*s^4 + 10*s^3 + 45*s^2 + 105*s + 105 etc.

Output is a Python list of arbitrary precision long ints, so n is only limited by your hardware’s memory.

Sequence is http://oeis.org/A001498 , and output can be confirmed to match http://oeis.org/A001498/b001498.txt :

>>> i = 0
>>> for n in range(51):
...     for x in _bessel_poly(n, reverse=True):
...         print(i, x)
...         i += 1
_campos_zeros(n)

Return approximate zero locations of Bessel polynomials y_n(x) for order n using polynomial fit (Campos-Calderon 2011)

_aberth(f, fp, x0, tol=1e-15, maxiter=50)

Given a function f, its first derivative fp, and a set of initial guesses x0, simultaneously find the roots of the polynomial using the Aberth-Ehrlich method.

len(x0) should equal the number of roots of f.

(This is not a complete implementation of Bini’s algorithm.)

_bessel_zeros(N)

Find zeros of ordinary Bessel polynomial of order N, by root-finding of modified Bessel function of the second kind

_norm_factor(p, k)

Numerically find frequency shift to apply to delay-normalized filter such that -3 dB point is at 1 rad/sec.

p is an array_like of polynomial poles k is a float gain

First 10 values are listed in “Bessel Scale Factors” table, “Bessel Filters Polynomials, Poles and Circuit Elements 2003, C. Bond.”

besselap(N, norm="phase")

Return (z,p,k) for analog prototype of an Nth-order Bessel filter.

N : int
The order of the filter.
norm : {‘phase’, ‘delay’, ‘mag’}, optional

Frequency normalization:

phase

The filter is normalized such that the phase response reaches its midpoint at an angular (e.g. rad/s) cutoff frequency of 1. This happens for both low-pass and high-pass filters, so this is the “phase-matched” case. [6]

The magnitude response asymptotes are the same as a Butterworth filter of the same order with a cutoff of Wn.

This is the default, and matches MATLAB’s implementation.

delay
The filter is normalized such that the group delay in the passband is 1 (e.g. 1 second). This is the “natural” type obtained by solving Bessel polynomials
mag
The filter is normalized such that the gain magnitude is -3 dB at angular frequency 1. This is called “frequency normalization” by Bond. [1]_

New in version 0.18.0.

z : ndarray
Zeros of the transfer function. Is always an empty array.
p : ndarray
Poles of the transfer function.
k : scalar
Gain of the transfer function. For phase-normalized, this is always 1.

bessel : Filter design function using this prototype

To find the pole locations, approximate starting points are generated [2] for the zeros of the ordinary Bessel polynomial [3], then the Aberth-Ehrlich method [4] [5] is used on the Kv(x) Bessel function to calculate more accurate zeros, and these locations are then inverted about the unit circle.

[1]C.R. Bond, “Bessel Filter Constants”, http://www.crbond.com/papers/bsf.pdf
[2]Campos and Calderon, “Approximate closed-form formulas for the zeros of the Bessel Polynomials”, :arXiv:`1105.0957`.
[3]Thomson, W.E., “Delay Networks having Maximally Flat Frequency Characteristics”, Proceedings of the Institution of Electrical Engineers, Part III, November 1949, Vol. 96, No. 44, pp. 487-490.
[4]Aberth, “Iteration Methods for Finding all Zeros of a Polynomial Simultaneously”, Mathematics of Computation, Vol. 27, No. 122, April 1973
[5]Ehrlich, “A modified Newton method for polynomials”, Communications of the ACM, Vol. 10, Issue 2, pp. 107-108, Feb. 1967, :DOI:`10.1145/363067.363115`
[6]Miller and Bohn, “A Bessel Filter Crossover, and Its Relation to Others”, RaneNote 147, 1998, http://www.rane.com/note147.html
iirnotch(w0, Q)

Design second-order IIR notch digital filter.

A notch filter is a band-stop filter with a narrow bandwidth (high quality factor). It rejects a narrow frequency band and leaves the rest of the spectrum little changed.

w0 : float
Normalized frequency to remove from a signal. It is a scalar that must satisfy 0 < w0 < 1, with w0 = 1 corresponding to half of the sampling frequency.
Q : float
Quality factor. Dimensionless parameter that characterizes notch filter -3 dB bandwidth bw relative to its center frequency, Q = w0/bw.
b, a : ndarray, ndarray
Numerator (b) and denominator (a) polynomials of the IIR filter.

iirpeak

[1]Sophocles J. Orfanidis, “Introduction To Signal Processing”, Prentice-Hall, 1996

Design and plot filter to remove the 60Hz component from a signal sampled at 200Hz, using a quality factor Q = 30

>>> from scipy import signal
>>> import numpy as np
>>> import matplotlib.pyplot as plt
>>> fs = 200.0  # Sample frequency (Hz)
>>> f0 = 60.0  # Frequency to be removed from signal (Hz)
>>> Q = 30.0  # Quality factor
>>> w0 = f0/(fs/2)  # Normalized Frequency
>>> # Design notch filter
>>> b, a = signal.iirnotch(w0, Q)
>>> # Frequency response
>>> w, h = signal.freqz(b, a)
>>> # Generate frequency axis
>>> freq = w*fs/(2*np.pi)
>>> # Plot
>>> fig, ax = plt.subplots(2, 1, figsize=(8, 6))
>>> ax[0].plot(freq, 20*np.log10(abs(h)), color='blue')
>>> ax[0].set_title("Frequency Response")
>>> ax[0].set_ylabel("Amplitude (dB)", color='blue')
>>> ax[0].set_xlim([0, 100])
>>> ax[0].set_ylim([-25, 10])
>>> ax[0].grid()
>>> ax[1].plot(freq, np.unwrap(np.angle(h))*180/np.pi, color='green')
>>> ax[1].set_ylabel("Angle (degrees)", color='green')
>>> ax[1].set_xlabel("Frequency (Hz)")
>>> ax[1].set_xlim([0, 100])
>>> ax[1].set_yticks([-90, -60, -30, 0, 30, 60, 90])
>>> ax[1].set_ylim([-90, 90])
>>> ax[1].grid()
>>> plt.show()
iirpeak(w0, Q)

Design second-order IIR peak (resonant) digital filter.

A peak filter is a band-pass filter with a narrow bandwidth (high quality factor). It rejects components outside a narrow frequency band.

w0 : float
Normalized frequency to be retained in a signal. It is a scalar that must satisfy 0 < w0 < 1, with w0 = 1 corresponding to half of the sampling frequency.
Q : float
Quality factor. Dimensionless parameter that characterizes peak filter -3 dB bandwidth bw relative to its center frequency, Q = w0/bw.
b, a : ndarray, ndarray
Numerator (b) and denominator (a) polynomials of the IIR filter.

iirnotch

[1]Sophocles J. Orfanidis, “Introduction To Signal Processing”, Prentice-Hall, 1996

Design and plot filter to remove the frequencies other than the 300Hz component from a signal sampled at 1000Hz, using a quality factor Q = 30

>>> from scipy import signal
>>> import numpy as np
>>> import matplotlib.pyplot as plt
>>> fs = 1000.0  # Sample frequency (Hz)
>>> f0 = 300.0  # Frequency to be retained (Hz)
>>> Q = 30.0  # Quality factor
>>> w0 = f0/(fs/2)  # Normalized Frequency
>>> # Design peak filter
>>> b, a = signal.iirpeak(w0, Q)
>>> # Frequency response
>>> w, h = signal.freqz(b, a)
>>> # Generate frequency axis
>>> freq = w*fs/(2*np.pi)
>>> # Plot
>>> fig, ax = plt.subplots(2, 1, figsize=(8, 6))
>>> ax[0].plot(freq, 20*np.log10(abs(h)), color='blue')
>>> ax[0].set_title("Frequency Response")
>>> ax[0].set_ylabel("Amplitude (dB)", color='blue')
>>> ax[0].set_xlim([0, 500])
>>> ax[0].set_ylim([-50, 10])
>>> ax[0].grid()
>>> ax[1].plot(freq, np.unwrap(np.angle(h))*180/np.pi, color='green')
>>> ax[1].set_ylabel("Angle (degrees)", color='green')
>>> ax[1].set_xlabel("Frequency (Hz)")
>>> ax[1].set_xlim([0, 500])
>>> ax[1].set_yticks([-90, -60, -30, 0, 30, 60, 90])
>>> ax[1].set_ylim([-90, 90])
>>> ax[1].grid()
>>> plt.show()
_design_notch_peak_filter(w0, Q, ftype)

Design notch or peak digital filter.

w0 : float
Normalized frequency to remove from a signal. It is a scalar that must satisfy 0 < w0 < 1, with w0 = 1 corresponding to half of the sampling frequency.
Q : float
Quality factor. Dimensionless parameter that characterizes notch filter -3 dB bandwidth bw relative to its center frequency, Q = w0/bw.
ftype : str

The type of IIR filter to design:

  • notch filter : notch
  • peak filter : peak
b, a : ndarray, ndarray
Numerator (b) and denominator (a) polynomials of the IIR filter.