Hide keyboard shortcuts

Hot-keys on this page

r m x p   toggle line displays

j k   next/prev highlighted chunk

0   (zero) top of page

1   (one) first highlighted chunk

1"""Filter design. 

2""" 

3import math 

4import operator 

5import warnings 

6 

7import numpy 

8import numpy as np 

9from numpy import (atleast_1d, poly, polyval, roots, real, asarray, 

10 resize, pi, absolute, logspace, r_, sqrt, tan, log10, 

11 arctan, arcsinh, sin, exp, cosh, arccosh, ceil, conjugate, 

12 zeros, sinh, append, concatenate, prod, ones, full, array, 

13 mintypecode) 

14from numpy.polynomial.polynomial import polyval as npp_polyval 

15from numpy.polynomial.polynomial import polyvalfromroots 

16 

17from scipy import special, optimize, fft as sp_fft 

18from scipy.special import comb, factorial 

19 

20 

21__all__ = ['findfreqs', 'freqs', 'freqz', 'tf2zpk', 'zpk2tf', 'normalize', 

22 'lp2lp', 'lp2hp', 'lp2bp', 'lp2bs', 'bilinear', 'iirdesign', 

23 'iirfilter', 'butter', 'cheby1', 'cheby2', 'ellip', 'bessel', 

24 'band_stop_obj', 'buttord', 'cheb1ord', 'cheb2ord', 'ellipord', 

25 'buttap', 'cheb1ap', 'cheb2ap', 'ellipap', 'besselap', 

26 'BadCoefficients', 'freqs_zpk', 'freqz_zpk', 

27 'tf2sos', 'sos2tf', 'zpk2sos', 'sos2zpk', 'group_delay', 

28 'sosfreqz', 'iirnotch', 'iirpeak', 'bilinear_zpk', 

29 'lp2lp_zpk', 'lp2hp_zpk', 'lp2bp_zpk', 'lp2bs_zpk'] 

30 

31 

32class BadCoefficients(UserWarning): 

33 """Warning about badly conditioned filter coefficients""" 

34 pass 

35 

36 

37abs = absolute 

38 

39 

40def _is_int_type(x): 

41 """ 

42 Check if input is of a scalar integer type (so ``5`` and ``array(5)`` will 

43 pass, while ``5.0`` and ``array([5])`` will fail. 

44 """ 

45 if np.ndim(x) != 0: 

46 # Older versions of NumPy did not raise for np.array([1]).__index__() 

47 # This is safe to remove when support for those versions is dropped 

48 return False 

49 try: 

50 operator.index(x) 

51 except TypeError: 

52 return False 

53 else: 

54 return True 

55 

56 

57def findfreqs(num, den, N, kind='ba'): 

58 """ 

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

60 

61 Parameters 

62 ---------- 

63 num, den : array_like, 1-D 

64 The polynomial coefficients of the numerator and denominator of the 

65 transfer function of the filter or LTI system, where the coefficients 

66 are ordered from highest to lowest degree. Or, the roots of the 

67 transfer function numerator and denominator (i.e., zeroes and poles). 

68 N : int 

69 The length of the array to be computed. 

70 kind : str {'ba', 'zp'}, optional 

71 Specifies whether the numerator and denominator are specified by their 

72 polynomial coefficients ('ba'), or their roots ('zp'). 

73 

74 Returns 

75 ------- 

76 w : (N,) ndarray 

77 A 1-D array of frequencies, logarithmically spaced. 

78 

79 Examples 

80 -------- 

81 Find a set of nine frequencies that span the "interesting part" of the 

82 frequency response for the filter with the transfer function 

83 

84 H(s) = s / (s^2 + 8s + 25) 

85 

86 >>> from scipy import signal 

87 >>> signal.findfreqs([1, 0], [1, 8, 25], N=9) 

88 array([ 1.00000000e-02, 3.16227766e-02, 1.00000000e-01, 

89 3.16227766e-01, 1.00000000e+00, 3.16227766e+00, 

90 1.00000000e+01, 3.16227766e+01, 1.00000000e+02]) 

91 """ 

92 if kind == 'ba': 

93 ep = atleast_1d(roots(den)) + 0j 

94 tz = atleast_1d(roots(num)) + 0j 

95 elif kind == 'zp': 

96 ep = atleast_1d(den) + 0j 

97 tz = atleast_1d(num) + 0j 

98 else: 

99 raise ValueError("input must be one of {'ba', 'zp'}") 

100 

101 if len(ep) == 0: 

102 ep = atleast_1d(-1000) + 0j 

103 

104 ez = r_['-1', 

105 numpy.compress(ep.imag >= 0, ep, axis=-1), 

106 numpy.compress((abs(tz) < 1e5) & (tz.imag >= 0), tz, axis=-1)] 

107 

108 integ = abs(ez) < 1e-10 

109 hfreq = numpy.around(numpy.log10(numpy.max(3 * abs(ez.real + integ) + 

110 1.5 * ez.imag)) + 0.5) 

111 lfreq = numpy.around(numpy.log10(0.1 * numpy.min(abs(real(ez + integ)) + 

112 2 * ez.imag)) - 0.5) 

113 

114 w = logspace(lfreq, hfreq, N) 

115 return w 

116 

117 

118def freqs(b, a, worN=200, plot=None): 

119 """ 

120 Compute frequency response of analog filter. 

121 

122 Given the M-order numerator `b` and N-order denominator `a` of an analog 

123 filter, compute its frequency response:: 

124 

125 b[0]*(jw)**M + b[1]*(jw)**(M-1) + ... + b[M] 

126 H(w) = ---------------------------------------------- 

127 a[0]*(jw)**N + a[1]*(jw)**(N-1) + ... + a[N] 

128 

129 Parameters 

130 ---------- 

131 b : array_like 

132 Numerator of a linear filter. 

133 a : array_like 

134 Denominator of a linear filter. 

135 worN : {None, int, array_like}, optional 

136 If None, then compute at 200 frequencies around the interesting parts 

137 of the response curve (determined by pole-zero locations). If a single 

138 integer, then compute at that many frequencies. Otherwise, compute the 

139 response at the angular frequencies (e.g., rad/s) given in `worN`. 

140 plot : callable, optional 

141 A callable that takes two arguments. If given, the return parameters 

142 `w` and `h` are passed to plot. Useful for plotting the frequency 

143 response inside `freqs`. 

144 

145 Returns 

146 ------- 

147 w : ndarray 

148 The angular frequencies at which `h` was computed. 

149 h : ndarray 

150 The frequency response. 

151 

152 See Also 

153 -------- 

154 freqz : Compute the frequency response of a digital filter. 

155 

156 Notes 

157 ----- 

158 Using Matplotlib's "plot" function as the callable for `plot` produces 

159 unexpected results, this plots the real part of the complex transfer 

160 function, not the magnitude. Try ``lambda w, h: plot(w, abs(h))``. 

161 

162 Examples 

163 -------- 

164 >>> from scipy.signal import freqs, iirfilter 

165 

166 >>> b, a = iirfilter(4, [1, 10], 1, 60, analog=True, ftype='cheby1') 

167 

168 >>> w, h = freqs(b, a, worN=np.logspace(-1, 2, 1000)) 

169 

170 >>> import matplotlib.pyplot as plt 

171 >>> plt.semilogx(w, 20 * np.log10(abs(h))) 

172 >>> plt.xlabel('Frequency') 

173 >>> plt.ylabel('Amplitude response [dB]') 

174 >>> plt.grid() 

175 >>> plt.show() 

176 

177 """ 

178 if worN is None: 

179 # For backwards compatibility 

180 w = findfreqs(b, a, 200) 

181 elif _is_int_type(worN): 

182 w = findfreqs(b, a, worN) 

183 else: 

184 w = atleast_1d(worN) 

185 

186 s = 1j * w 

187 h = polyval(b, s) / polyval(a, s) 

188 if plot is not None: 

189 plot(w, h) 

190 

191 return w, h 

192 

193 

194def freqs_zpk(z, p, k, worN=200): 

195 """ 

196 Compute frequency response of analog filter. 

197 

198 Given the zeros `z`, poles `p`, and gain `k` of a filter, compute its 

199 frequency response:: 

200 

201 (jw-z[0]) * (jw-z[1]) * ... * (jw-z[-1]) 

202 H(w) = k * ---------------------------------------- 

203 (jw-p[0]) * (jw-p[1]) * ... * (jw-p[-1]) 

204 

205 Parameters 

206 ---------- 

207 z : array_like 

208 Zeroes of a linear filter 

209 p : array_like 

210 Poles of a linear filter 

211 k : scalar 

212 Gain of a linear filter 

213 worN : {None, int, array_like}, optional 

214 If None, then compute at 200 frequencies around the interesting parts 

215 of the response curve (determined by pole-zero locations). If a single 

216 integer, then compute at that many frequencies. Otherwise, compute the 

217 response at the angular frequencies (e.g., rad/s) given in `worN`. 

218 

219 Returns 

220 ------- 

221 w : ndarray 

222 The angular frequencies at which `h` was computed. 

223 h : ndarray 

224 The frequency response. 

225 

226 See Also 

227 -------- 

228 freqs : Compute the frequency response of an analog filter in TF form 

229 freqz : Compute the frequency response of a digital filter in TF form 

230 freqz_zpk : Compute the frequency response of a digital filter in ZPK form 

231 

232 Notes 

233 ----- 

234 .. versionadded:: 0.19.0 

235 

236 Examples 

237 -------- 

238 >>> from scipy.signal import freqs_zpk, iirfilter 

239 

240 >>> z, p, k = iirfilter(4, [1, 10], 1, 60, analog=True, ftype='cheby1', 

241 ... output='zpk') 

242 

243 >>> w, h = freqs_zpk(z, p, k, worN=np.logspace(-1, 2, 1000)) 

244 

245 >>> import matplotlib.pyplot as plt 

246 >>> plt.semilogx(w, 20 * np.log10(abs(h))) 

247 >>> plt.xlabel('Frequency') 

248 >>> plt.ylabel('Amplitude response [dB]') 

249 >>> plt.grid() 

250 >>> plt.show() 

251 

252 """ 

253 k = np.asarray(k) 

254 if k.size > 1: 

255 raise ValueError('k must be a single scalar gain') 

256 

257 if worN is None: 

258 # For backwards compatibility 

259 w = findfreqs(z, p, 200, kind='zp') 

260 elif _is_int_type(worN): 

261 w = findfreqs(z, p, worN, kind='zp') 

262 else: 

263 w = worN 

264 

265 w = atleast_1d(w) 

266 s = 1j * w 

267 num = polyvalfromroots(s, z) 

268 den = polyvalfromroots(s, p) 

269 h = k * num/den 

270 return w, h 

271 

272 

273def freqz(b, a=1, worN=512, whole=False, plot=None, fs=2*pi, include_nyquist=False): 

274 """ 

275 Compute the frequency response of a digital filter. 

276 

277 Given the M-order numerator `b` and N-order denominator `a` of a digital 

278 filter, compute its frequency response:: 

279 

280 jw -jw -jwM 

281 jw B(e ) b[0] + b[1]e + ... + b[M]e 

282 H(e ) = ------ = ----------------------------------- 

283 jw -jw -jwN 

284 A(e ) a[0] + a[1]e + ... + a[N]e 

285 

286 Parameters 

287 ---------- 

288 b : array_like 

289 Numerator of a linear filter. If `b` has dimension greater than 1, 

290 it is assumed that the coefficients are stored in the first dimension, 

291 and ``b.shape[1:]``, ``a.shape[1:]``, and the shape of the frequencies 

292 array must be compatible for broadcasting. 

293 a : array_like 

294 Denominator of a linear filter. If `b` has dimension greater than 1, 

295 it is assumed that the coefficients are stored in the first dimension, 

296 and ``b.shape[1:]``, ``a.shape[1:]``, and the shape of the frequencies 

297 array must be compatible for broadcasting. 

298 worN : {None, int, array_like}, optional 

299 If a single integer, then compute at that many frequencies (default is 

300 N=512). This is a convenient alternative to:: 

301 

302 np.linspace(0, fs if whole else fs/2, N, endpoint=include_nyquist) 

303 

304 Using a number that is fast for FFT computations can result in 

305 faster computations (see Notes). 

306 

307 If an array_like, compute the response at the frequencies given. 

308 These are in the same units as `fs`. 

309 whole : bool, optional 

310 Normally, frequencies are computed from 0 to the Nyquist frequency, 

311 fs/2 (upper-half of unit-circle). If `whole` is True, compute 

312 frequencies from 0 to fs. Ignored if w is array_like. 

313 plot : callable 

314 A callable that takes two arguments. If given, the return parameters 

315 `w` and `h` are passed to plot. Useful for plotting the frequency 

316 response inside `freqz`. 

317 fs : float, optional 

318 The sampling frequency of the digital system. Defaults to 2*pi 

319 radians/sample (so w is from 0 to pi). 

320 

321 .. versionadded:: 1.2.0 

322 include_nyquist : bool, optional 

323 If `whole` is False and `worN` is an integer, setting `include_nyquist` to True 

324 will include the last frequency (Nyquist frequency) and is otherwise ignored. 

325 

326 .. versionadded:: 1.5.0 

327 

328 Returns 

329 ------- 

330 w : ndarray 

331 The frequencies at which `h` was computed, in the same units as `fs`. 

332 By default, `w` is normalized to the range [0, pi) (radians/sample). 

333 h : ndarray 

334 The frequency response, as complex numbers. 

335 

336 See Also 

337 -------- 

338 freqz_zpk 

339 sosfreqz 

340 

341 Notes 

342 ----- 

343 Using Matplotlib's :func:`matplotlib.pyplot.plot` function as the callable 

344 for `plot` produces unexpected results, as this plots the real part of the 

345 complex transfer function, not the magnitude. 

346 Try ``lambda w, h: plot(w, np.abs(h))``. 

347 

348 A direct computation via (R)FFT is used to compute the frequency response 

349 when the following conditions are met: 

350 

351 1. An integer value is given for `worN`. 

352 2. `worN` is fast to compute via FFT (i.e., 

353 `next_fast_len(worN) <scipy.fft.next_fast_len>` equals `worN`). 

354 3. The denominator coefficients are a single value (``a.shape[0] == 1``). 

355 4. `worN` is at least as long as the numerator coefficients 

356 (``worN >= b.shape[0]``). 

357 5. If ``b.ndim > 1``, then ``b.shape[-1] == 1``. 

358 

359 For long FIR filters, the FFT approach can have lower error and be much 

360 faster than the equivalent direct polynomial calculation. 

361 

362 Examples 

363 -------- 

364 >>> from scipy import signal 

365 >>> b = signal.firwin(80, 0.5, window=('kaiser', 8)) 

366 >>> w, h = signal.freqz(b) 

367 

368 >>> import matplotlib.pyplot as plt 

369 >>> fig, ax1 = plt.subplots() 

370 >>> ax1.set_title('Digital filter frequency response') 

371 

372 >>> ax1.plot(w, 20 * np.log10(abs(h)), 'b') 

373 >>> ax1.set_ylabel('Amplitude [dB]', color='b') 

374 >>> ax1.set_xlabel('Frequency [rad/sample]') 

375 

376 >>> ax2 = ax1.twinx() 

377 >>> angles = np.unwrap(np.angle(h)) 

378 >>> ax2.plot(w, angles, 'g') 

379 >>> ax2.set_ylabel('Angle (radians)', color='g') 

380 >>> ax2.grid() 

381 >>> ax2.axis('tight') 

382 >>> plt.show() 

383 

384 Broadcasting Examples 

385 

386 Suppose we have two FIR filters whose coefficients are stored in the 

387 rows of an array with shape (2, 25). For this demonstration, we'll 

388 use random data: 

389 

390 >>> np.random.seed(42) 

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

392 

393 To compute the frequency response for these two filters with one call 

394 to `freqz`, we must pass in ``b.T``, because `freqz` expects the first 

395 axis to hold the coefficients. We must then extend the shape with a 

396 trivial dimension of length 1 to allow broadcasting with the array 

397 of frequencies. That is, we pass in ``b.T[..., np.newaxis]``, which has 

398 shape (25, 2, 1): 

399 

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

401 >>> w.shape 

402 (1024,) 

403 >>> h.shape 

404 (2, 1024) 

405 

406 Now, suppose we have two transfer functions, with the same numerator 

407 coefficients ``b = [0.5, 0.5]``. The coefficients for the two denominators 

408 are stored in the first dimension of the 2-D array `a`:: 

409 

410 a = [ 1 1 ] 

411 [ -0.25, -0.5 ] 

412 

413 >>> b = np.array([0.5, 0.5]) 

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

415 

416 Only `a` is more than 1-D. To make it compatible for 

417 broadcasting with the frequencies, we extend it with a trivial dimension 

418 in the call to `freqz`: 

419 

420 >>> w, h = signal.freqz(b, a[..., np.newaxis], worN=1024) 

421 >>> w.shape 

422 (1024,) 

423 >>> h.shape 

424 (2, 1024) 

425 

426 """ 

427 b = atleast_1d(b) 

428 a = atleast_1d(a) 

429 

430 if worN is None: 

431 # For backwards compatibility 

432 worN = 512 

433 

434 h = None 

435 

436 if _is_int_type(worN): 

437 N = operator.index(worN) 

438 del worN 

439 if N < 0: 

440 raise ValueError('worN must be nonnegative, got %s' % (N,)) 

441 lastpoint = 2 * pi if whole else pi 

442 # if include_nyquist is true and whole is false, w should include end point 

443 w = np.linspace(0, lastpoint, N, endpoint=include_nyquist and not whole) 

444 if (a.size == 1 and N >= b.shape[0] and 

445 sp_fft.next_fast_len(N) == N and 

446 (b.ndim == 1 or (b.shape[-1] == 1))): 

447 # if N is fast, 2 * N will be fast, too, so no need to check 

448 n_fft = N if whole else N * 2 

449 if np.isrealobj(b) and np.isrealobj(a): 

450 fft_func = sp_fft.rfft 

451 else: 

452 fft_func = sp_fft.fft 

453 h = fft_func(b, n=n_fft, axis=0)[:N] 

454 h /= a 

455 if fft_func is sp_fft.rfft and whole: 

456 # exclude DC and maybe Nyquist (no need to use axis_reverse 

457 # here because we can build reversal with the truncation) 

458 stop = -1 if n_fft % 2 == 1 else -2 

459 h_flip = slice(stop, 0, -1) 

460 h = np.concatenate((h, h[h_flip].conj())) 

461 if b.ndim > 1: 

462 # Last axis of h has length 1, so drop it. 

463 h = h[..., 0] 

464 # Rotate the first axis of h to the end. 

465 h = np.rollaxis(h, 0, h.ndim) 

466 else: 

467 w = atleast_1d(worN) 

468 del worN 

469 w = 2*pi*w/fs 

470 

471 if h is None: # still need to compute using freqs w 

472 zm1 = exp(-1j * w) 

473 h = (npp_polyval(zm1, b, tensor=False) / 

474 npp_polyval(zm1, a, tensor=False)) 

475 

476 w = w*fs/(2*pi) 

477 

478 if plot is not None: 

479 plot(w, h) 

480 

481 return w, h 

482 

483 

484def freqz_zpk(z, p, k, worN=512, whole=False, fs=2*pi): 

485 r""" 

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

487 

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

489 response: 

490 

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

492 

493 where :math:`k` is the `gain`, :math:`Z` are the `zeros` and :math:`P` are 

494 the `poles`. 

495 

496 Parameters 

497 ---------- 

498 z : array_like 

499 Zeroes of a linear filter 

500 p : array_like 

501 Poles of a linear filter 

502 k : scalar 

503 Gain of a linear filter 

504 worN : {None, int, array_like}, optional 

505 If a single integer, then compute at that many frequencies (default is 

506 N=512). 

507 

508 If an array_like, compute the response at the frequencies given. 

509 These are in the same units as `fs`. 

510 whole : bool, optional 

511 Normally, frequencies are computed from 0 to the Nyquist frequency, 

512 fs/2 (upper-half of unit-circle). If `whole` is True, compute 

513 frequencies from 0 to fs. Ignored if w is array_like. 

514 fs : float, optional 

515 The sampling frequency of the digital system. Defaults to 2*pi 

516 radians/sample (so w is from 0 to pi). 

517 

518 .. versionadded:: 1.2.0 

519 

520 Returns 

521 ------- 

522 w : ndarray 

523 The frequencies at which `h` was computed, in the same units as `fs`. 

524 By default, `w` is normalized to the range [0, pi) (radians/sample). 

525 h : ndarray 

526 The frequency response, as complex numbers. 

527 

528 See Also 

529 -------- 

530 freqs : Compute the frequency response of an analog filter in TF form 

531 freqs_zpk : Compute the frequency response of an analog filter in ZPK form 

532 freqz : Compute the frequency response of a digital filter in TF form 

533 

534 Notes 

535 ----- 

536 .. versionadded:: 0.19.0 

537 

538 Examples 

539 -------- 

540 Design a 4th-order digital Butterworth filter with cut-off of 100 Hz in a 

541 system with sample rate of 1000 Hz, and plot the frequency response: 

542 

543 >>> from scipy import signal 

544 >>> z, p, k = signal.butter(4, 100, output='zpk', fs=1000) 

545 >>> w, h = signal.freqz_zpk(z, p, k, fs=1000) 

546 

547 >>> import matplotlib.pyplot as plt 

548 >>> fig = plt.figure() 

549 >>> ax1 = fig.add_subplot(1, 1, 1) 

550 >>> ax1.set_title('Digital filter frequency response') 

551 

552 >>> ax1.plot(w, 20 * np.log10(abs(h)), 'b') 

553 >>> ax1.set_ylabel('Amplitude [dB]', color='b') 

554 >>> ax1.set_xlabel('Frequency [Hz]') 

555 >>> ax1.grid() 

556 

557 >>> ax2 = ax1.twinx() 

558 >>> angles = np.unwrap(np.angle(h)) 

559 >>> ax2.plot(w, angles, 'g') 

560 >>> ax2.set_ylabel('Angle [radians]', color='g') 

561 

562 >>> plt.axis('tight') 

563 >>> plt.show() 

564 

565 """ 

566 z, p = map(atleast_1d, (z, p)) 

567 

568 if whole: 

569 lastpoint = 2 * pi 

570 else: 

571 lastpoint = pi 

572 

573 if worN is None: 

574 # For backwards compatibility 

575 w = numpy.linspace(0, lastpoint, 512, endpoint=False) 

576 elif _is_int_type(worN): 

577 w = numpy.linspace(0, lastpoint, worN, endpoint=False) 

578 else: 

579 w = atleast_1d(worN) 

580 w = 2*pi*w/fs 

581 

582 zm1 = exp(1j * w) 

583 h = k * polyvalfromroots(zm1, z) / polyvalfromroots(zm1, p) 

584 

585 w = w*fs/(2*pi) 

586 

587 return w, h 

588 

589 

590def group_delay(system, w=512, whole=False, fs=2*pi): 

591 r"""Compute the group delay of a digital filter. 

592 

593 The group delay measures by how many samples amplitude envelopes of 

594 various spectral components of a signal are delayed by a filter. 

595 It is formally defined as the derivative of continuous (unwrapped) phase:: 

596 

597 d jw 

598 D(w) = - -- arg H(e) 

599 dw 

600 

601 Parameters 

602 ---------- 

603 system : tuple of array_like (b, a) 

604 Numerator and denominator coefficients of a filter transfer function. 

605 w : {None, int, array_like}, optional 

606 If a single integer, then compute at that many frequencies (default is 

607 N=512). 

608 

609 If an array_like, compute the delay at the frequencies given. These 

610 are in the same units as `fs`. 

611 whole : bool, optional 

612 Normally, frequencies are computed from 0 to the Nyquist frequency, 

613 fs/2 (upper-half of unit-circle). If `whole` is True, compute 

614 frequencies from 0 to fs. Ignored if w is array_like. 

615 fs : float, optional 

616 The sampling frequency of the digital system. Defaults to 2*pi 

617 radians/sample (so w is from 0 to pi). 

618 

619 .. versionadded:: 1.2.0 

620 

621 Returns 

622 ------- 

623 w : ndarray 

624 The frequencies at which group delay was computed, in the same units 

625 as `fs`. By default, `w` is normalized to the range [0, pi) 

626 (radians/sample). 

627 gd : ndarray 

628 The group delay. 

629 

630 Notes 

631 ----- 

632 The similar function in MATLAB is called `grpdelay`. 

633 

634 If the transfer function :math:`H(z)` has zeros or poles on the unit 

635 circle, the group delay at corresponding frequencies is undefined. 

636 When such a case arises the warning is raised and the group delay 

637 is set to 0 at those frequencies. 

638 

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

640 

641 .. versionadded:: 0.16.0 

642 

643 See Also 

644 -------- 

645 freqz : Frequency response of a digital filter 

646 

647 References 

648 ---------- 

649 .. [1] Richard G. Lyons, "Understanding Digital Signal Processing, 

650 3rd edition", p. 830. 

651 

652 Examples 

653 -------- 

654 >>> from scipy import signal 

655 >>> b, a = signal.iirdesign(0.1, 0.3, 5, 50, ftype='cheby1') 

656 >>> w, gd = signal.group_delay((b, a)) 

657 

658 >>> import matplotlib.pyplot as plt 

659 >>> plt.title('Digital filter group delay') 

660 >>> plt.plot(w, gd) 

661 >>> plt.ylabel('Group delay [samples]') 

662 >>> plt.xlabel('Frequency [rad/sample]') 

663 >>> plt.show() 

664 

665 """ 

666 if w is None: 

667 # For backwards compatibility 

668 w = 512 

669 

670 if _is_int_type(w): 

671 if whole: 

672 w = np.linspace(0, 2 * pi, w, endpoint=False) 

673 else: 

674 w = np.linspace(0, pi, w, endpoint=False) 

675 else: 

676 w = np.atleast_1d(w) 

677 w = 2*pi*w/fs 

678 

679 b, a = map(np.atleast_1d, system) 

680 c = np.convolve(b, a[::-1]) 

681 cr = c * np.arange(c.size) 

682 z = np.exp(-1j * w) 

683 num = np.polyval(cr[::-1], z) 

684 den = np.polyval(c[::-1], z) 

685 singular = np.absolute(den) < 10 * EPSILON 

686 if np.any(singular): 

687 warnings.warn( 

688 "The group delay is singular at frequencies [{0}], setting to 0". 

689 format(", ".join("{0:.3f}".format(ws) for ws in w[singular])) 

690 ) 

691 

692 gd = np.zeros_like(w) 

693 gd[~singular] = np.real(num[~singular] / den[~singular]) - a.size + 1 

694 

695 w = w*fs/(2*pi) 

696 

697 return w, gd 

698 

699 

700def _validate_sos(sos): 

701 """Helper to validate a SOS input""" 

702 sos = np.atleast_2d(sos) 

703 if sos.ndim != 2: 

704 raise ValueError('sos array must be 2D') 

705 n_sections, m = sos.shape 

706 if m != 6: 

707 raise ValueError('sos array must be shape (n_sections, 6)') 

708 if not (sos[:, 3] == 1).all(): 

709 raise ValueError('sos[:, 3] should be all ones') 

710 return sos, n_sections 

711 

712 

713def sosfreqz(sos, worN=512, whole=False, fs=2*pi): 

714 r""" 

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

716 

717 Given `sos`, an array with shape (n, 6) of second order sections of 

718 a digital filter, compute the frequency response of the system function:: 

719 

720 B0(z) B1(z) B{n-1}(z) 

721 H(z) = ----- * ----- * ... * --------- 

722 A0(z) A1(z) A{n-1}(z) 

723 

724 for z = exp(omega*1j), where B{k}(z) and A{k}(z) are numerator and 

725 denominator of the transfer function of the k-th second order section. 

726 

727 Parameters 

728 ---------- 

729 sos : array_like 

730 Array of second-order filter coefficients, must have shape 

731 ``(n_sections, 6)``. Each row corresponds to a second-order 

732 section, with the first three columns providing the numerator 

733 coefficients and the last three providing the denominator 

734 coefficients. 

735 worN : {None, int, array_like}, optional 

736 If a single integer, then compute at that many frequencies (default is 

737 N=512). Using a number that is fast for FFT computations can result 

738 in faster computations (see Notes of `freqz`). 

739 

740 If an array_like, compute the response at the frequencies given (must 

741 be 1-D). These are in the same units as `fs`. 

742 whole : bool, optional 

743 Normally, frequencies are computed from 0 to the Nyquist frequency, 

744 fs/2 (upper-half of unit-circle). If `whole` is True, compute 

745 frequencies from 0 to fs. 

746 fs : float, optional 

747 The sampling frequency of the digital system. Defaults to 2*pi 

748 radians/sample (so w is from 0 to pi). 

749 

750 .. versionadded:: 1.2.0 

751 

752 Returns 

753 ------- 

754 w : ndarray 

755 The frequencies at which `h` was computed, in the same units as `fs`. 

756 By default, `w` is normalized to the range [0, pi) (radians/sample). 

757 h : ndarray 

758 The frequency response, as complex numbers. 

759 

760 See Also 

761 -------- 

762 freqz, sosfilt 

763 

764 Notes 

765 ----- 

766 .. versionadded:: 0.19.0 

767 

768 Examples 

769 -------- 

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

771 

772 >>> from scipy import signal 

773 >>> sos = signal.ellip(15, 0.5, 60, (0.2, 0.4), btype='bandpass', 

774 ... output='sos') 

775 

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

777 

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

779 

780 Plot the response. 

781 

782 >>> import matplotlib.pyplot as plt 

783 >>> plt.subplot(2, 1, 1) 

784 >>> db = 20*np.log10(np.maximum(np.abs(h), 1e-5)) 

785 >>> plt.plot(w/np.pi, db) 

786 >>> plt.ylim(-75, 5) 

787 >>> plt.grid(True) 

788 >>> plt.yticks([0, -20, -40, -60]) 

789 >>> plt.ylabel('Gain [dB]') 

790 >>> plt.title('Frequency Response') 

791 >>> plt.subplot(2, 1, 2) 

792 >>> plt.plot(w/np.pi, np.angle(h)) 

793 >>> plt.grid(True) 

794 >>> plt.yticks([-np.pi, -0.5*np.pi, 0, 0.5*np.pi, np.pi], 

795 ... [r'$-\pi$', r'$-\pi/2$', '0', r'$\pi/2$', r'$\pi$']) 

796 >>> plt.ylabel('Phase [rad]') 

797 >>> plt.xlabel('Normalized frequency (1.0 = Nyquist)') 

798 >>> plt.show() 

799 

800 If the same filter is implemented as a single transfer function, 

801 numerical error corrupts the frequency response: 

802 

803 >>> b, a = signal.ellip(15, 0.5, 60, (0.2, 0.4), btype='bandpass', 

804 ... output='ba') 

805 >>> w, h = signal.freqz(b, a, worN=1500) 

806 >>> plt.subplot(2, 1, 1) 

807 >>> db = 20*np.log10(np.maximum(np.abs(h), 1e-5)) 

808 >>> plt.plot(w/np.pi, db) 

809 >>> plt.ylim(-75, 5) 

810 >>> plt.grid(True) 

811 >>> plt.yticks([0, -20, -40, -60]) 

812 >>> plt.ylabel('Gain [dB]') 

813 >>> plt.title('Frequency Response') 

814 >>> plt.subplot(2, 1, 2) 

815 >>> plt.plot(w/np.pi, np.angle(h)) 

816 >>> plt.grid(True) 

817 >>> plt.yticks([-np.pi, -0.5*np.pi, 0, 0.5*np.pi, np.pi], 

818 ... [r'$-\pi$', r'$-\pi/2$', '0', r'$\pi/2$', r'$\pi$']) 

819 >>> plt.ylabel('Phase [rad]') 

820 >>> plt.xlabel('Normalized frequency (1.0 = Nyquist)') 

821 >>> plt.show() 

822 

823 """ 

824 

825 sos, n_sections = _validate_sos(sos) 

826 if n_sections == 0: 

827 raise ValueError('Cannot compute frequencies with no sections') 

828 h = 1. 

829 for row in sos: 

830 w, rowh = freqz(row[:3], row[3:], worN=worN, whole=whole, fs=fs) 

831 h *= rowh 

832 return w, h 

833 

834 

835def _cplxreal(z, tol=None): 

836 """ 

837 Split into complex and real parts, combining conjugate pairs. 

838 

839 The 1-D input vector `z` is split up into its complex (`zc`) and real (`zr`) 

840 elements. Every complex element must be part of a complex-conjugate pair, 

841 which are combined into a single number (with positive imaginary part) in 

842 the output. Two complex numbers are considered a conjugate pair if their 

843 real and imaginary parts differ in magnitude by less than ``tol * abs(z)``. 

844 

845 Parameters 

846 ---------- 

847 z : array_like 

848 Vector of complex numbers to be sorted and split 

849 tol : float, optional 

850 Relative tolerance for testing realness and conjugate equality. 

851 Default is ``100 * spacing(1)`` of `z`'s data type (i.e., 2e-14 for 

852 float64) 

853 

854 Returns 

855 ------- 

856 zc : ndarray 

857 Complex elements of `z`, with each pair represented by a single value 

858 having positive imaginary part, sorted first by real part, and then 

859 by magnitude of imaginary part. The pairs are averaged when combined 

860 to reduce error. 

861 zr : ndarray 

862 Real elements of `z` (those having imaginary part less than 

863 `tol` times their magnitude), sorted by value. 

864 

865 Raises 

866 ------ 

867 ValueError 

868 If there are any complex numbers in `z` for which a conjugate 

869 cannot be found. 

870 

871 See Also 

872 -------- 

873 _cplxpair 

874 

875 Examples 

876 -------- 

877 >>> a = [4, 3, 1, 2-2j, 2+2j, 2-1j, 2+1j, 2-1j, 2+1j, 1+1j, 1-1j] 

878 >>> zc, zr = _cplxreal(a) 

879 >>> print(zc) 

880 [ 1.+1.j 2.+1.j 2.+1.j 2.+2.j] 

881 >>> print(zr) 

882 [ 1. 3. 4.] 

883 """ 

884 

885 z = atleast_1d(z) 

886 if z.size == 0: 

887 return z, z 

888 elif z.ndim != 1: 

889 raise ValueError('_cplxreal only accepts 1-D input') 

890 

891 if tol is None: 

892 # Get tolerance from dtype of input 

893 tol = 100 * np.finfo((1.0 * z).dtype).eps 

894 

895 # Sort by real part, magnitude of imaginary part (speed up further sorting) 

896 z = z[np.lexsort((abs(z.imag), z.real))] 

897 

898 # Split reals from conjugate pairs 

899 real_indices = abs(z.imag) <= tol * abs(z) 

900 zr = z[real_indices].real 

901 

902 if len(zr) == len(z): 

903 # Input is entirely real 

904 return array([]), zr 

905 

906 # Split positive and negative halves of conjugates 

907 z = z[~real_indices] 

908 zp = z[z.imag > 0] 

909 zn = z[z.imag < 0] 

910 

911 if len(zp) != len(zn): 

912 raise ValueError('Array contains complex value with no matching ' 

913 'conjugate.') 

914 

915 # Find runs of (approximately) the same real part 

916 same_real = np.diff(zp.real) <= tol * abs(zp[:-1]) 

917 diffs = numpy.diff(concatenate(([0], same_real, [0]))) 

918 run_starts = numpy.nonzero(diffs > 0)[0] 

919 run_stops = numpy.nonzero(diffs < 0)[0] 

920 

921 # Sort each run by their imaginary parts 

922 for i in range(len(run_starts)): 

923 start = run_starts[i] 

924 stop = run_stops[i] + 1 

925 for chunk in (zp[start:stop], zn[start:stop]): 

926 chunk[...] = chunk[np.lexsort([abs(chunk.imag)])] 

927 

928 # Check that negatives match positives 

929 if any(abs(zp - zn.conj()) > tol * abs(zn)): 

930 raise ValueError('Array contains complex value with no matching ' 

931 'conjugate.') 

932 

933 # Average out numerical inaccuracy in real vs imag parts of pairs 

934 zc = (zp + zn.conj()) / 2 

935 

936 return zc, zr 

937 

938 

939def _cplxpair(z, tol=None): 

940 """ 

941 Sort into pairs of complex conjugates. 

942 

943 Complex conjugates in `z` are sorted by increasing real part. In each 

944 pair, the number with negative imaginary part appears first. 

945 

946 If pairs have identical real parts, they are sorted by increasing 

947 imaginary magnitude. 

948 

949 Two complex numbers are considered a conjugate pair if their real and 

950 imaginary parts differ in magnitude by less than ``tol * abs(z)``. The 

951 pairs are forced to be exact complex conjugates by averaging the positive 

952 and negative values. 

953 

954 Purely real numbers are also sorted, but placed after the complex 

955 conjugate pairs. A number is considered real if its imaginary part is 

956 smaller than `tol` times the magnitude of the number. 

957 

958 Parameters 

959 ---------- 

960 z : array_like 

961 1-D input array to be sorted. 

962 tol : float, optional 

963 Relative tolerance for testing realness and conjugate equality. 

964 Default is ``100 * spacing(1)`` of `z`'s data type (i.e., 2e-14 for 

965 float64) 

966 

967 Returns 

968 ------- 

969 y : ndarray 

970 Complex conjugate pairs followed by real numbers. 

971 

972 Raises 

973 ------ 

974 ValueError 

975 If there are any complex numbers in `z` for which a conjugate 

976 cannot be found. 

977 

978 See Also 

979 -------- 

980 _cplxreal 

981 

982 Examples 

983 -------- 

984 >>> a = [4, 3, 1, 2-2j, 2+2j, 2-1j, 2+1j, 2-1j, 2+1j, 1+1j, 1-1j] 

985 >>> z = _cplxpair(a) 

986 >>> print(z) 

987 [ 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 

988 3.+0.j 4.+0.j] 

989 """ 

990 

991 z = atleast_1d(z) 

992 if z.size == 0 or np.isrealobj(z): 

993 return np.sort(z) 

994 

995 if z.ndim != 1: 

996 raise ValueError('z must be 1-D') 

997 

998 zc, zr = _cplxreal(z, tol) 

999 

1000 # Interleave complex values and their conjugates, with negative imaginary 

1001 # parts first in each pair 

1002 zc = np.dstack((zc.conj(), zc)).flatten() 

1003 z = np.append(zc, zr) 

1004 return z 

1005 

1006 

1007def tf2zpk(b, a): 

1008 r"""Return zero, pole, gain (z, p, k) representation from a numerator, 

1009 denominator representation of a linear filter. 

1010 

1011 Parameters 

1012 ---------- 

1013 b : array_like 

1014 Numerator polynomial coefficients. 

1015 a : array_like 

1016 Denominator polynomial coefficients. 

1017 

1018 Returns 

1019 ------- 

1020 z : ndarray 

1021 Zeros of the transfer function. 

1022 p : ndarray 

1023 Poles of the transfer function. 

1024 k : float 

1025 System gain. 

1026 

1027 Notes 

1028 ----- 

1029 If some values of `b` are too close to 0, they are removed. In that case, 

1030 a BadCoefficients warning is emitted. 

1031 

1032 The `b` and `a` arrays are interpreted as coefficients for positive, 

1033 descending powers of the transfer function variable. So the inputs 

1034 :math:`b = [b_0, b_1, ..., b_M]` and :math:`a =[a_0, a_1, ..., a_N]` 

1035 can represent an analog filter of the form: 

1036 

1037 .. math:: 

1038 

1039 H(s) = \frac 

1040 {b_0 s^M + b_1 s^{(M-1)} + \cdots + b_M} 

1041 {a_0 s^N + a_1 s^{(N-1)} + \cdots + a_N} 

1042 

1043 or a discrete-time filter of the form: 

1044 

1045 .. math:: 

1046 

1047 H(z) = \frac 

1048 {b_0 z^M + b_1 z^{(M-1)} + \cdots + b_M} 

1049 {a_0 z^N + a_1 z^{(N-1)} + \cdots + a_N} 

1050 

1051 This "positive powers" form is found more commonly in controls 

1052 engineering. If `M` and `N` are equal (which is true for all filters 

1053 generated by the bilinear transform), then this happens to be equivalent 

1054 to the "negative powers" discrete-time form preferred in DSP: 

1055 

1056 .. math:: 

1057 

1058 H(z) = \frac 

1059 {b_0 + b_1 z^{-1} + \cdots + b_M z^{-M}} 

1060 {a_0 + a_1 z^{-1} + \cdots + a_N z^{-N}} 

1061 

1062 Although this is true for common filters, remember that this is not true 

1063 in the general case. If `M` and `N` are not equal, the discrete-time 

1064 transfer function coefficients must first be converted to the "positive 

1065 powers" form before finding the poles and zeros. 

1066 

1067 """ 

1068 b, a = normalize(b, a) 

1069 b = (b + 0.0) / a[0] 

1070 a = (a + 0.0) / a[0] 

1071 k = b[0] 

1072 b /= b[0] 

1073 z = roots(b) 

1074 p = roots(a) 

1075 return z, p, k 

1076 

1077 

1078def zpk2tf(z, p, k): 

1079 """ 

1080 Return polynomial transfer function representation from zeros and poles 

1081 

1082 Parameters 

1083 ---------- 

1084 z : array_like 

1085 Zeros of the transfer function. 

1086 p : array_like 

1087 Poles of the transfer function. 

1088 k : float 

1089 System gain. 

1090 

1091 Returns 

1092 ------- 

1093 b : ndarray 

1094 Numerator polynomial coefficients. 

1095 a : ndarray 

1096 Denominator polynomial coefficients. 

1097 

1098 """ 

1099 z = atleast_1d(z) 

1100 k = atleast_1d(k) 

1101 if len(z.shape) > 1: 

1102 temp = poly(z[0]) 

1103 b = zeros((z.shape[0], z.shape[1] + 1), temp.dtype.char) 

1104 if len(k) == 1: 

1105 k = [k[0]] * z.shape[0] 

1106 for i in range(z.shape[0]): 

1107 b[i] = k[i] * poly(z[i]) 

1108 else: 

1109 b = k * poly(z) 

1110 a = atleast_1d(poly(p)) 

1111 

1112 # Use real output if possible. Copied from numpy.poly, since 

1113 # we can't depend on a specific version of numpy. 

1114 if issubclass(b.dtype.type, numpy.complexfloating): 

1115 # if complex roots are all complex conjugates, the roots are real. 

1116 roots = numpy.asarray(z, complex) 

1117 pos_roots = numpy.compress(roots.imag > 0, roots) 

1118 neg_roots = numpy.conjugate(numpy.compress(roots.imag < 0, roots)) 

1119 if len(pos_roots) == len(neg_roots): 

1120 if numpy.all(numpy.sort_complex(neg_roots) == 

1121 numpy.sort_complex(pos_roots)): 

1122 b = b.real.copy() 

1123 

1124 if issubclass(a.dtype.type, numpy.complexfloating): 

1125 # if complex roots are all complex conjugates, the roots are real. 

1126 roots = numpy.asarray(p, complex) 

1127 pos_roots = numpy.compress(roots.imag > 0, roots) 

1128 neg_roots = numpy.conjugate(numpy.compress(roots.imag < 0, roots)) 

1129 if len(pos_roots) == len(neg_roots): 

1130 if numpy.all(numpy.sort_complex(neg_roots) == 

1131 numpy.sort_complex(pos_roots)): 

1132 a = a.real.copy() 

1133 

1134 return b, a 

1135 

1136 

1137def tf2sos(b, a, pairing='nearest'): 

1138 """ 

1139 Return second-order sections from transfer function representation 

1140 

1141 Parameters 

1142 ---------- 

1143 b : array_like 

1144 Numerator polynomial coefficients. 

1145 a : array_like 

1146 Denominator polynomial coefficients. 

1147 pairing : {'nearest', 'keep_odd'}, optional 

1148 The method to use to combine pairs of poles and zeros into sections. 

1149 See `zpk2sos`. 

1150 

1151 Returns 

1152 ------- 

1153 sos : ndarray 

1154 Array of second-order filter coefficients, with shape 

1155 ``(n_sections, 6)``. See `sosfilt` for the SOS filter format 

1156 specification. 

1157 

1158 See Also 

1159 -------- 

1160 zpk2sos, sosfilt 

1161 

1162 Notes 

1163 ----- 

1164 It is generally discouraged to convert from TF to SOS format, since doing 

1165 so usually will not improve numerical precision errors. Instead, consider 

1166 designing filters in ZPK format and converting directly to SOS. TF is 

1167 converted to SOS by first converting to ZPK format, then converting 

1168 ZPK to SOS. 

1169 

1170 .. versionadded:: 0.16.0 

1171 """ 

1172 return zpk2sos(*tf2zpk(b, a), pairing=pairing) 

1173 

1174 

1175def sos2tf(sos): 

1176 """ 

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

1178 

1179 Parameters 

1180 ---------- 

1181 sos : array_like 

1182 Array of second-order filter coefficients, must have shape 

1183 ``(n_sections, 6)``. See `sosfilt` for the SOS filter format 

1184 specification. 

1185 

1186 Returns 

1187 ------- 

1188 b : ndarray 

1189 Numerator polynomial coefficients. 

1190 a : ndarray 

1191 Denominator polynomial coefficients. 

1192 

1193 Notes 

1194 ----- 

1195 .. versionadded:: 0.16.0 

1196 """ 

1197 sos = np.asarray(sos) 

1198 b = [1.] 

1199 a = [1.] 

1200 n_sections = sos.shape[0] 

1201 for section in range(n_sections): 

1202 b = np.polymul(b, sos[section, :3]) 

1203 a = np.polymul(a, sos[section, 3:]) 

1204 return b, a 

1205 

1206 

1207def sos2zpk(sos): 

1208 """ 

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

1210 

1211 Parameters 

1212 ---------- 

1213 sos : array_like 

1214 Array of second-order filter coefficients, must have shape 

1215 ``(n_sections, 6)``. See `sosfilt` for the SOS filter format 

1216 specification. 

1217 

1218 Returns 

1219 ------- 

1220 z : ndarray 

1221 Zeros of the transfer function. 

1222 p : ndarray 

1223 Poles of the transfer function. 

1224 k : float 

1225 System gain. 

1226 

1227 Notes 

1228 ----- 

1229 The number of zeros and poles returned will be ``n_sections * 2`` 

1230 even if some of these are (effectively) zero. 

1231 

1232 .. versionadded:: 0.16.0 

1233 """ 

1234 sos = np.asarray(sos) 

1235 n_sections = sos.shape[0] 

1236 z = np.zeros(n_sections*2, np.complex128) 

1237 p = np.zeros(n_sections*2, np.complex128) 

1238 k = 1. 

1239 for section in range(n_sections): 

1240 zpk = tf2zpk(sos[section, :3], sos[section, 3:]) 

1241 z[2*section:2*section+len(zpk[0])] = zpk[0] 

1242 p[2*section:2*section+len(zpk[1])] = zpk[1] 

1243 k *= zpk[2] 

1244 return z, p, k 

1245 

1246 

1247def _nearest_real_complex_idx(fro, to, which): 

1248 """Get the next closest real or complex element based on distance""" 

1249 assert which in ('real', 'complex') 

1250 order = np.argsort(np.abs(fro - to)) 

1251 mask = np.isreal(fro[order]) 

1252 if which == 'complex': 

1253 mask = ~mask 

1254 return order[np.nonzero(mask)[0][0]] 

1255 

1256 

1257def zpk2sos(z, p, k, pairing='nearest'): 

1258 """ 

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

1260 

1261 Parameters 

1262 ---------- 

1263 z : array_like 

1264 Zeros of the transfer function. 

1265 p : array_like 

1266 Poles of the transfer function. 

1267 k : float 

1268 System gain. 

1269 pairing : {'nearest', 'keep_odd'}, optional 

1270 The method to use to combine pairs of poles and zeros into sections. 

1271 See Notes below. 

1272 

1273 Returns 

1274 ------- 

1275 sos : ndarray 

1276 Array of second-order filter coefficients, with shape 

1277 ``(n_sections, 6)``. See `sosfilt` for the SOS filter format 

1278 specification. 

1279 

1280 See Also 

1281 -------- 

1282 sosfilt 

1283 

1284 Notes 

1285 ----- 

1286 The algorithm used to convert ZPK to SOS format is designed to 

1287 minimize errors due to numerical precision issues. The pairing 

1288 algorithm attempts to minimize the peak gain of each biquadratic 

1289 section. This is done by pairing poles with the nearest zeros, starting 

1290 with the poles closest to the unit circle. 

1291 

1292 *Algorithms* 

1293 

1294 The current algorithms are designed specifically for use with digital 

1295 filters. (The output coefficients are not correct for analog filters.) 

1296 

1297 The steps in the ``pairing='nearest'`` and ``pairing='keep_odd'`` 

1298 algorithms are mostly shared. The ``nearest`` algorithm attempts to 

1299 minimize the peak gain, while ``'keep_odd'`` minimizes peak gain under 

1300 the constraint that odd-order systems should retain one section 

1301 as first order. The algorithm steps and are as follows: 

1302 

1303 As a pre-processing step, add poles or zeros to the origin as 

1304 necessary to obtain the same number of poles and zeros for pairing. 

1305 If ``pairing == 'nearest'`` and there are an odd number of poles, 

1306 add an additional pole and a zero at the origin. 

1307 

1308 The following steps are then iterated over until no more poles or 

1309 zeros remain: 

1310 

1311 1. Take the (next remaining) pole (complex or real) closest to the 

1312 unit circle to begin a new filter section. 

1313 

1314 2. If the pole is real and there are no other remaining real poles [#]_, 

1315 add the closest real zero to the section and leave it as a first 

1316 order section. Note that after this step we are guaranteed to be 

1317 left with an even number of real poles, complex poles, real zeros, 

1318 and complex zeros for subsequent pairing iterations. 

1319 

1320 3. Else: 

1321 

1322 1. If the pole is complex and the zero is the only remaining real 

1323 zero*, then pair the pole with the *next* closest zero 

1324 (guaranteed to be complex). This is necessary to ensure that 

1325 there will be a real zero remaining to eventually create a 

1326 first-order section (thus keeping the odd order). 

1327 

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

1329 real). 

1330 

1331 3. Proceed to complete the second-order section by adding another 

1332 pole and zero to the current pole and zero in the section: 

1333 

1334 1. If the current pole and zero are both complex, add their 

1335 conjugates. 

1336 

1337 2. Else if the pole is complex and the zero is real, add the 

1338 conjugate pole and the next closest real zero. 

1339 

1340 3. Else if the pole is real and the zero is complex, add the 

1341 conjugate zero and the real pole closest to those zeros. 

1342 

1343 4. Else (we must have a real pole and real zero) add the next 

1344 real pole closest to the unit circle, and then add the real 

1345 zero closest to that pole. 

1346 

1347 .. [#] This conditional can only be met for specific odd-order inputs 

1348 with the ``pairing == 'keep_odd'`` method. 

1349 

1350 .. versionadded:: 0.16.0 

1351 

1352 Examples 

1353 -------- 

1354 

1355 Design a 6th order low-pass elliptic digital filter for a system with a 

1356 sampling rate of 8000 Hz that has a pass-band corner frequency of 

1357 1000 Hz. The ripple in the pass-band should not exceed 0.087 dB, and 

1358 the attenuation in the stop-band should be at least 90 dB. 

1359 

1360 In the following call to `signal.ellip`, we could use ``output='sos'``, 

1361 but for this example, we'll use ``output='zpk'``, and then convert to SOS 

1362 format with `zpk2sos`: 

1363 

1364 >>> from scipy import signal 

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

1366 

1367 Now convert to SOS format. 

1368 

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

1370 

1371 The coefficients of the numerators of the sections: 

1372 

1373 >>> sos[:, :3] 

1374 array([[ 0.0014154 , 0.00248707, 0.0014154 ], 

1375 [ 1. , 0.72965193, 1. ], 

1376 [ 1. , 0.17594966, 1. ]]) 

1377 

1378 The symmetry in the coefficients occurs because all the zeros are on the 

1379 unit circle. 

1380 

1381 The coefficients of the denominators of the sections: 

1382 

1383 >>> sos[:, 3:] 

1384 array([[ 1. , -1.32543251, 0.46989499], 

1385 [ 1. , -1.26117915, 0.6262586 ], 

1386 [ 1. , -1.25707217, 0.86199667]]) 

1387 

1388 The next example shows the effect of the `pairing` option. We have a 

1389 system with three poles and three zeros, so the SOS array will have 

1390 shape (2, 6). The means there is, in effect, an extra pole and an extra 

1391 zero at the origin in the SOS representation. 

1392 

1393 >>> z1 = np.array([-1, -0.5-0.5j, -0.5+0.5j]) 

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

1395 

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

1397 

1398 >>> signal.zpk2sos(z1, p1, 1) 

1399 array([[ 1. , 1. , 0.5 , 1. , -0.75, 0. ], 

1400 [ 1. , 1. , 0. , 1. , -1.6 , 0.65]]) 

1401 

1402 The first section has the zeros {-0.5-0.05j, -0.5+0.5j} and the poles 

1403 {0, 0.75}, and the second section has the zeros {-1, 0} and poles 

1404 {0.8+0.1j, 0.8-0.1j}. Note that the extra pole and zero at the origin 

1405 have been assigned to different sections. 

1406 

1407 With ``pairing='keep_odd'``, we obtain: 

1408 

1409 >>> signal.zpk2sos(z1, p1, 1, pairing='keep_odd') 

1410 array([[ 1. , 1. , 0. , 1. , -0.75, 0. ], 

1411 [ 1. , 1. , 0.5 , 1. , -1.6 , 0.65]]) 

1412 

1413 The extra pole and zero at the origin are in the same section. 

1414 The first section is, in effect, a first-order section. 

1415 

1416 """ 

1417 # TODO in the near future: 

1418 # 1. Add SOS capability to `filtfilt`, `freqz`, etc. somehow (#3259). 

1419 # 2. Make `decimate` use `sosfilt` instead of `lfilter`. 

1420 # 3. Make sosfilt automatically simplify sections to first order 

1421 # when possible. Note this might make `sosfiltfilt` a bit harder (ICs). 

1422 # 4. Further optimizations of the section ordering / pole-zero pairing. 

1423 # See the wiki for other potential issues. 

1424 

1425 valid_pairings = ['nearest', 'keep_odd'] 

1426 if pairing not in valid_pairings: 

1427 raise ValueError('pairing must be one of %s, not %s' 

1428 % (valid_pairings, pairing)) 

1429 if len(z) == len(p) == 0: 

1430 return array([[k, 0., 0., 1., 0., 0.]]) 

1431 

1432 # ensure we have the same number of poles and zeros, and make copies 

1433 p = np.concatenate((p, np.zeros(max(len(z) - len(p), 0)))) 

1434 z = np.concatenate((z, np.zeros(max(len(p) - len(z), 0)))) 

1435 n_sections = (max(len(p), len(z)) + 1) // 2 

1436 sos = zeros((n_sections, 6)) 

1437 

1438 if len(p) % 2 == 1 and pairing == 'nearest': 

1439 p = np.concatenate((p, [0.])) 

1440 z = np.concatenate((z, [0.])) 

1441 assert len(p) == len(z) 

1442 

1443 # Ensure we have complex conjugate pairs 

1444 # (note that _cplxreal only gives us one element of each complex pair): 

1445 z = np.concatenate(_cplxreal(z)) 

1446 p = np.concatenate(_cplxreal(p)) 

1447 

1448 p_sos = np.zeros((n_sections, 2), np.complex128) 

1449 z_sos = np.zeros_like(p_sos) 

1450 for si in range(n_sections): 

1451 # Select the next "worst" pole 

1452 p1_idx = np.argmin(np.abs(1 - np.abs(p))) 

1453 p1 = p[p1_idx] 

1454 p = np.delete(p, p1_idx) 

1455 

1456 # Pair that pole with a zero 

1457 

1458 if np.isreal(p1) and np.isreal(p).sum() == 0: 

1459 # Special case to set a first-order section 

1460 z1_idx = _nearest_real_complex_idx(z, p1, 'real') 

1461 z1 = z[z1_idx] 

1462 z = np.delete(z, z1_idx) 

1463 p2 = z2 = 0 

1464 else: 

1465 if not np.isreal(p1) and np.isreal(z).sum() == 1: 

1466 # Special case to ensure we choose a complex zero to pair 

1467 # with so later (setting up a first-order section) 

1468 z1_idx = _nearest_real_complex_idx(z, p1, 'complex') 

1469 assert not np.isreal(z[z1_idx]) 

1470 else: 

1471 # Pair the pole with the closest zero (real or complex) 

1472 z1_idx = np.argmin(np.abs(p1 - z)) 

1473 z1 = z[z1_idx] 

1474 z = np.delete(z, z1_idx) 

1475 

1476 # Now that we have p1 and z1, figure out what p2 and z2 need to be 

1477 if not np.isreal(p1): 

1478 if not np.isreal(z1): # complex pole, complex zero 

1479 p2 = p1.conj() 

1480 z2 = z1.conj() 

1481 else: # complex pole, real zero 

1482 p2 = p1.conj() 

1483 z2_idx = _nearest_real_complex_idx(z, p1, 'real') 

1484 z2 = z[z2_idx] 

1485 assert np.isreal(z2) 

1486 z = np.delete(z, z2_idx) 

1487 else: 

1488 if not np.isreal(z1): # real pole, complex zero 

1489 z2 = z1.conj() 

1490 p2_idx = _nearest_real_complex_idx(p, z1, 'real') 

1491 p2 = p[p2_idx] 

1492 assert np.isreal(p2) 

1493 else: # real pole, real zero 

1494 # pick the next "worst" pole to use 

1495 idx = np.nonzero(np.isreal(p))[0] 

1496 assert len(idx) > 0 

1497 p2_idx = idx[np.argmin(np.abs(np.abs(p[idx]) - 1))] 

1498 p2 = p[p2_idx] 

1499 # find a real zero to match the added pole 

1500 assert np.isreal(p2) 

1501 z2_idx = _nearest_real_complex_idx(z, p2, 'real') 

1502 z2 = z[z2_idx] 

1503 assert np.isreal(z2) 

1504 z = np.delete(z, z2_idx) 

1505 p = np.delete(p, p2_idx) 

1506 p_sos[si] = [p1, p2] 

1507 z_sos[si] = [z1, z2] 

1508 assert len(p) == len(z) == 0 # we've consumed all poles and zeros 

1509 del p, z 

1510 

1511 # Construct the system, reversing order so the "worst" are last 

1512 p_sos = np.reshape(p_sos[::-1], (n_sections, 2)) 

1513 z_sos = np.reshape(z_sos[::-1], (n_sections, 2)) 

1514 gains = np.ones(n_sections, np.array(k).dtype) 

1515 gains[0] = k 

1516 for si in range(n_sections): 

1517 x = zpk2tf(z_sos[si], p_sos[si], gains[si]) 

1518 sos[si] = np.concatenate(x) 

1519 return sos 

1520 

1521 

1522def _align_nums(nums): 

1523 """Aligns the shapes of multiple numerators. 

1524 

1525 Given an array of numerator coefficient arrays [[a_1, a_2,..., 

1526 a_n],..., [b_1, b_2,..., b_m]], this function pads shorter numerator 

1527 arrays with zero's so that all numerators have the same length. Such 

1528 alignment is necessary for functions like 'tf2ss', which needs the 

1529 alignment when dealing with SIMO transfer functions. 

1530 

1531 Parameters 

1532 ---------- 

1533 nums: array_like 

1534 Numerator or list of numerators. Not necessarily with same length. 

1535 

1536 Returns 

1537 ------- 

1538 nums: array 

1539 The numerator. If `nums` input was a list of numerators then a 2-D 

1540 array with padded zeros for shorter numerators is returned. Otherwise 

1541 returns ``np.asarray(nums)``. 

1542 """ 

1543 try: 

1544 # The statement can throw a ValueError if one 

1545 # of the numerators is a single digit and another 

1546 # is array-like e.g. if nums = [5, [1, 2, 3]] 

1547 nums = asarray(nums) 

1548 

1549 if not np.issubdtype(nums.dtype, np.number): 

1550 raise ValueError("dtype of numerator is non-numeric") 

1551 

1552 return nums 

1553 

1554 except ValueError: 

1555 nums = [np.atleast_1d(num) for num in nums] 

1556 max_width = max(num.size for num in nums) 

1557 

1558 # pre-allocate 

1559 aligned_nums = np.zeros((len(nums), max_width)) 

1560 

1561 # Create numerators with padded zeros 

1562 for index, num in enumerate(nums): 

1563 aligned_nums[index, -num.size:] = num 

1564 

1565 return aligned_nums 

1566 

1567 

1568def normalize(b, a): 

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

1570 

1571 If values of `b` are too close to 0, they are removed. In that case, a 

1572 BadCoefficients warning is emitted. 

1573 

1574 Parameters 

1575 ---------- 

1576 b: array_like 

1577 Numerator of the transfer function. Can be a 2-D array to normalize 

1578 multiple transfer functions. 

1579 a: array_like 

1580 Denominator of the transfer function. At most 1-D. 

1581 

1582 Returns 

1583 ------- 

1584 num: array 

1585 The numerator of the normalized transfer function. At least a 1-D 

1586 array. A 2-D array if the input `num` is a 2-D array. 

1587 den: 1-D array 

1588 The denominator of the normalized transfer function. 

1589 

1590 Notes 

1591 ----- 

1592 Coefficients for both the numerator and denominator should be specified in 

1593 descending exponent order (e.g., ``s^2 + 3s + 5`` would be represented as 

1594 ``[1, 3, 5]``). 

1595 """ 

1596 num, den = b, a 

1597 

1598 den = np.atleast_1d(den) 

1599 num = np.atleast_2d(_align_nums(num)) 

1600 

1601 if den.ndim != 1: 

1602 raise ValueError("Denominator polynomial must be rank-1 array.") 

1603 if num.ndim > 2: 

1604 raise ValueError("Numerator polynomial must be rank-1 or" 

1605 " rank-2 array.") 

1606 if np.all(den == 0): 

1607 raise ValueError("Denominator must have at least on nonzero element.") 

1608 

1609 # Trim leading zeros in denominator, leave at least one. 

1610 den = np.trim_zeros(den, 'f') 

1611 

1612 # Normalize transfer function 

1613 num, den = num / den[0], den / den[0] 

1614 

1615 # Count numerator columns that are all zero 

1616 leading_zeros = 0 

1617 for col in num.T: 

1618 if np.allclose(col, 0, atol=1e-14): 

1619 leading_zeros += 1 

1620 else: 

1621 break 

1622 

1623 # Trim leading zeros of numerator 

1624 if leading_zeros > 0: 

1625 warnings.warn("Badly conditioned filter coefficients (numerator): the " 

1626 "results may be meaningless", BadCoefficients) 

1627 # Make sure at least one column remains 

1628 if leading_zeros == num.shape[1]: 

1629 leading_zeros -= 1 

1630 num = num[:, leading_zeros:] 

1631 

1632 # Squeeze first dimension if singular 

1633 if num.shape[0] == 1: 

1634 num = num[0, :] 

1635 

1636 return num, den 

1637 

1638 

1639def lp2lp(b, a, wo=1.0): 

1640 r""" 

1641 Transform a lowpass filter prototype to a different frequency. 

1642 

1643 Return an analog low-pass filter with cutoff frequency `wo` 

1644 from an analog low-pass filter prototype with unity cutoff frequency, in 

1645 transfer function ('ba') representation. 

1646 

1647 Parameters 

1648 ---------- 

1649 b : array_like 

1650 Numerator polynomial coefficients. 

1651 a : array_like 

1652 Denominator polynomial coefficients. 

1653 wo : float 

1654 Desired cutoff, as angular frequency (e.g. rad/s). 

1655 Defaults to no change. 

1656 

1657 Returns 

1658 ------- 

1659 b : array_like 

1660 Numerator polynomial coefficients of the transformed low-pass filter. 

1661 a : array_like 

1662 Denominator polynomial coefficients of the transformed low-pass filter. 

1663 

1664 See Also 

1665 -------- 

1666 lp2hp, lp2bp, lp2bs, bilinear 

1667 lp2lp_zpk 

1668 

1669 Notes 

1670 ----- 

1671 This is derived from the s-plane substitution 

1672 

1673 .. math:: s \rightarrow \frac{s}{\omega_0} 

1674 

1675 Examples 

1676 -------- 

1677 

1678 >>> from scipy import signal 

1679 >>> import matplotlib.pyplot as plt 

1680 

1681 >>> lp = signal.lti([1.0], [1.0, 1.0]) 

1682 >>> lp2 = signal.lti(*signal.lp2lp(lp.num, lp.den, 2)) 

1683 >>> w, mag_lp, p_lp = lp.bode() 

1684 >>> w, mag_lp2, p_lp2 = lp2.bode(w) 

1685 

1686 >>> plt.plot(w, mag_lp, label='Lowpass') 

1687 >>> plt.plot(w, mag_lp2, label='Transformed Lowpass') 

1688 >>> plt.semilogx() 

1689 >>> plt.grid() 

1690 >>> plt.xlabel('Frequency [rad/s]') 

1691 >>> plt.ylabel('Magnitude [dB]') 

1692 >>> plt.legend() 

1693 

1694 """ 

1695 a, b = map(atleast_1d, (a, b)) 

1696 try: 

1697 wo = float(wo) 

1698 except TypeError: 

1699 wo = float(wo[0]) 

1700 d = len(a) 

1701 n = len(b) 

1702 M = max((d, n)) 

1703 pwo = pow(wo, numpy.arange(M - 1, -1, -1)) 

1704 start1 = max((n - d, 0)) 

1705 start2 = max((d - n, 0)) 

1706 b = b * pwo[start1] / pwo[start2:] 

1707 a = a * pwo[start1] / pwo[start1:] 

1708 return normalize(b, a) 

1709 

1710 

1711def lp2hp(b, a, wo=1.0): 

1712 r""" 

1713 Transform a lowpass filter prototype to a highpass filter. 

1714 

1715 Return an analog high-pass filter with cutoff frequency `wo` 

1716 from an analog low-pass filter prototype with unity cutoff frequency, in 

1717 transfer function ('ba') representation. 

1718 

1719 Parameters 

1720 ---------- 

1721 b : array_like 

1722 Numerator polynomial coefficients. 

1723 a : array_like 

1724 Denominator polynomial coefficients. 

1725 wo : float 

1726 Desired cutoff, as angular frequency (e.g., rad/s). 

1727 Defaults to no change. 

1728 

1729 Returns 

1730 ------- 

1731 b : array_like 

1732 Numerator polynomial coefficients of the transformed high-pass filter. 

1733 a : array_like 

1734 Denominator polynomial coefficients of the transformed high-pass filter. 

1735 

1736 See Also 

1737 -------- 

1738 lp2lp, lp2bp, lp2bs, bilinear 

1739 lp2hp_zpk 

1740 

1741 Notes 

1742 ----- 

1743 This is derived from the s-plane substitution 

1744 

1745 .. math:: s \rightarrow \frac{\omega_0}{s} 

1746 

1747 This maintains symmetry of the lowpass and highpass responses on a 

1748 logarithmic scale. 

1749 

1750 Examples 

1751 -------- 

1752 >>> from scipy import signal 

1753 >>> import matplotlib.pyplot as plt 

1754 

1755 >>> lp = signal.lti([1.0], [1.0, 1.0]) 

1756 >>> hp = signal.lti(*signal.lp2hp(lp.num, lp.den)) 

1757 >>> w, mag_lp, p_lp = lp.bode() 

1758 >>> w, mag_hp, p_hp = hp.bode(w) 

1759 

1760 >>> plt.plot(w, mag_lp, label='Lowpass') 

1761 >>> plt.plot(w, mag_hp, label='Highpass') 

1762 >>> plt.semilogx() 

1763 >>> plt.grid() 

1764 >>> plt.xlabel('Frequency [rad/s]') 

1765 >>> plt.ylabel('Magnitude [dB]') 

1766 >>> plt.legend() 

1767 

1768 """ 

1769 a, b = map(atleast_1d, (a, b)) 

1770 try: 

1771 wo = float(wo) 

1772 except TypeError: 

1773 wo = float(wo[0]) 

1774 d = len(a) 

1775 n = len(b) 

1776 if wo != 1: 

1777 pwo = pow(wo, numpy.arange(max((d, n)))) 

1778 else: 

1779 pwo = numpy.ones(max((d, n)), b.dtype.char) 

1780 if d >= n: 

1781 outa = a[::-1] * pwo 

1782 outb = resize(b, (d,)) 

1783 outb[n:] = 0.0 

1784 outb[:n] = b[::-1] * pwo[:n] 

1785 else: 

1786 outb = b[::-1] * pwo 

1787 outa = resize(a, (n,)) 

1788 outa[d:] = 0.0 

1789 outa[:d] = a[::-1] * pwo[:d] 

1790 

1791 return normalize(outb, outa) 

1792 

1793 

1794def lp2bp(b, a, wo=1.0, bw=1.0): 

1795 r""" 

1796 Transform a lowpass filter prototype to a bandpass filter. 

1797 

1798 Return an analog band-pass filter with center frequency `wo` and 

1799 bandwidth `bw` from an analog low-pass filter prototype with unity 

1800 cutoff frequency, in transfer function ('ba') representation. 

1801 

1802 Parameters 

1803 ---------- 

1804 b : array_like 

1805 Numerator polynomial coefficients. 

1806 a : array_like 

1807 Denominator polynomial coefficients. 

1808 wo : float 

1809 Desired passband center, as angular frequency (e.g., rad/s). 

1810 Defaults to no change. 

1811 bw : float 

1812 Desired passband width, as angular frequency (e.g., rad/s). 

1813 Defaults to 1. 

1814 

1815 Returns 

1816 ------- 

1817 b : array_like 

1818 Numerator polynomial coefficients of the transformed band-pass filter. 

1819 a : array_like 

1820 Denominator polynomial coefficients of the transformed band-pass filter. 

1821 

1822 See Also 

1823 -------- 

1824 lp2lp, lp2hp, lp2bs, bilinear 

1825 lp2bp_zpk 

1826 

1827 Notes 

1828 ----- 

1829 This is derived from the s-plane substitution 

1830 

1831 .. math:: s \rightarrow \frac{s^2 + {\omega_0}^2}{s \cdot \mathrm{BW}} 

1832 

1833 This is the "wideband" transformation, producing a passband with 

1834 geometric (log frequency) symmetry about `wo`. 

1835 

1836 Examples 

1837 -------- 

1838 >>> from scipy import signal 

1839 >>> import matplotlib.pyplot as plt 

1840 

1841 >>> lp = signal.lti([1.0], [1.0, 1.0]) 

1842 >>> bp = signal.lti(*signal.lp2bp(lp.num, lp.den)) 

1843 >>> w, mag_lp, p_lp = lp.bode() 

1844 >>> w, mag_bp, p_bp = bp.bode(w) 

1845 

1846 >>> plt.plot(w, mag_lp, label='Lowpass') 

1847 >>> plt.plot(w, mag_bp, label='Bandpass') 

1848 >>> plt.semilogx() 

1849 >>> plt.grid() 

1850 >>> plt.xlabel('Frequency [rad/s]') 

1851 >>> plt.ylabel('Magnitude [dB]') 

1852 >>> plt.legend() 

1853 """ 

1854 

1855 a, b = map(atleast_1d, (a, b)) 

1856 D = len(a) - 1 

1857 N = len(b) - 1 

1858 artype = mintypecode((a, b)) 

1859 ma = max([N, D]) 

1860 Np = N + ma 

1861 Dp = D + ma 

1862 bprime = numpy.zeros(Np + 1, artype) 

1863 aprime = numpy.zeros(Dp + 1, artype) 

1864 wosq = wo * wo 

1865 for j in range(Np + 1): 

1866 val = 0.0 

1867 for i in range(0, N + 1): 

1868 for k in range(0, i + 1): 

1869 if ma - i + 2 * k == j: 

1870 val += comb(i, k) * b[N - i] * (wosq) ** (i - k) / bw ** i 

1871 bprime[Np - j] = val 

1872 for j in range(Dp + 1): 

1873 val = 0.0 

1874 for i in range(0, D + 1): 

1875 for k in range(0, i + 1): 

1876 if ma - i + 2 * k == j: 

1877 val += comb(i, k) * a[D - i] * (wosq) ** (i - k) / bw ** i 

1878 aprime[Dp - j] = val 

1879 

1880 return normalize(bprime, aprime) 

1881 

1882 

1883def lp2bs(b, a, wo=1.0, bw=1.0): 

1884 r""" 

1885 Transform a lowpass filter prototype to a bandstop filter. 

1886 

1887 Return an analog band-stop filter with center frequency `wo` and 

1888 bandwidth `bw` from an analog low-pass filter prototype with unity 

1889 cutoff frequency, in transfer function ('ba') representation. 

1890 

1891 Parameters 

1892 ---------- 

1893 b : array_like 

1894 Numerator polynomial coefficients. 

1895 a : array_like 

1896 Denominator polynomial coefficients. 

1897 wo : float 

1898 Desired stopband center, as angular frequency (e.g., rad/s). 

1899 Defaults to no change. 

1900 bw : float 

1901 Desired stopband width, as angular frequency (e.g., rad/s). 

1902 Defaults to 1. 

1903 

1904 Returns 

1905 ------- 

1906 b : array_like 

1907 Numerator polynomial coefficients of the transformed band-stop filter. 

1908 a : array_like 

1909 Denominator polynomial coefficients of the transformed band-stop filter. 

1910 

1911 See Also 

1912 -------- 

1913 lp2lp, lp2hp, lp2bp, bilinear 

1914 lp2bs_zpk 

1915 

1916 Notes 

1917 ----- 

1918 This is derived from the s-plane substitution 

1919 

1920 .. math:: s \rightarrow \frac{s \cdot \mathrm{BW}}{s^2 + {\omega_0}^2} 

1921 

1922 This is the "wideband" transformation, producing a stopband with 

1923 geometric (log frequency) symmetry about `wo`. 

1924 

1925 Examples 

1926 -------- 

1927 >>> from scipy import signal 

1928 >>> import matplotlib.pyplot as plt 

1929 

1930 >>> lp = signal.lti([1.0], [1.0, 1.5]) 

1931 >>> bs = signal.lti(*signal.lp2bs(lp.num, lp.den)) 

1932 >>> w, mag_lp, p_lp = lp.bode() 

1933 >>> w, mag_bs, p_bs = bs.bode(w) 

1934 >>> plt.plot(w, mag_lp, label='Lowpass') 

1935 >>> plt.plot(w, mag_bs, label='Bandstop') 

1936 >>> plt.semilogx() 

1937 >>> plt.grid() 

1938 >>> plt.xlabel('Frequency [rad/s]') 

1939 >>> plt.ylabel('Magnitude [dB]') 

1940 >>> plt.legend() 

1941 """ 

1942 a, b = map(atleast_1d, (a, b)) 

1943 D = len(a) - 1 

1944 N = len(b) - 1 

1945 artype = mintypecode((a, b)) 

1946 M = max([N, D]) 

1947 Np = M + M 

1948 Dp = M + M 

1949 bprime = numpy.zeros(Np + 1, artype) 

1950 aprime = numpy.zeros(Dp + 1, artype) 

1951 wosq = wo * wo 

1952 for j in range(Np + 1): 

1953 val = 0.0 

1954 for i in range(0, N + 1): 

1955 for k in range(0, M - i + 1): 

1956 if i + 2 * k == j: 

1957 val += (comb(M - i, k) * b[N - i] * 

1958 (wosq) ** (M - i - k) * bw ** i) 

1959 bprime[Np - j] = val 

1960 for j in range(Dp + 1): 

1961 val = 0.0 

1962 for i in range(0, D + 1): 

1963 for k in range(0, M - i + 1): 

1964 if i + 2 * k == j: 

1965 val += (comb(M - i, k) * a[D - i] * 

1966 (wosq) ** (M - i - k) * bw ** i) 

1967 aprime[Dp - j] = val 

1968 

1969 return normalize(bprime, aprime) 

1970 

1971 

1972def bilinear(b, a, fs=1.0): 

1973 r""" 

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

1975 

1976 Transform a set of poles and zeros from the analog s-plane to the digital 

1977 z-plane using Tustin's method, which substitutes ``(z-1) / (z+1)`` for 

1978 ``s``, maintaining the shape of the frequency response. 

1979 

1980 Parameters 

1981 ---------- 

1982 b : array_like 

1983 Numerator of the analog filter transfer function. 

1984 a : array_like 

1985 Denominator of the analog filter transfer function. 

1986 fs : float 

1987 Sample rate, as ordinary frequency (e.g., hertz). No prewarping is 

1988 done in this function. 

1989 

1990 Returns 

1991 ------- 

1992 z : ndarray 

1993 Numerator of the transformed digital filter transfer function. 

1994 p : ndarray 

1995 Denominator of the transformed digital filter transfer function. 

1996 

1997 See Also 

1998 -------- 

1999 lp2lp, lp2hp, lp2bp, lp2bs 

2000 bilinear_zpk 

2001 

2002 Examples 

2003 -------- 

2004 >>> from scipy import signal 

2005 >>> import matplotlib.pyplot as plt 

2006 

2007 >>> fs = 100 

2008 >>> bf = 2 * np.pi * np.array([7, 13]) 

2009 >>> filts = signal.lti(*signal.butter(4, bf, btype='bandpass', analog=True)) 

2010 >>> filtz = signal.lti(*signal.bilinear(filts.num, filts.den, fs)) 

2011 >>> wz, hz = signal.freqz(filtz.num, filtz.den) 

2012 >>> ws, hs = signal.freqs(filts.num, filts.den, worN=fs*wz) 

2013 

2014 >>> plt.semilogx(wz*fs/(2*np.pi), 20*np.log10(np.abs(hz).clip(1e-15)), label=r'$|H(j \omega)|$') 

2015 >>> plt.semilogx(wz*fs/(2*np.pi), 20*np.log10(np.abs(hs).clip(1e-15)), label=r'$|H_z(e^{j \omega})|$') 

2016 >>> plt.legend() 

2017 >>> plt.xlabel('Frequency [Hz]') 

2018 >>> plt.ylabel('Magnitude [dB]') 

2019 >>> plt.grid() 

2020 """ 

2021 fs = float(fs) 

2022 a, b = map(atleast_1d, (a, b)) 

2023 D = len(a) - 1 

2024 N = len(b) - 1 

2025 artype = float 

2026 M = max([N, D]) 

2027 Np = M 

2028 Dp = M 

2029 bprime = numpy.zeros(Np + 1, artype) 

2030 aprime = numpy.zeros(Dp + 1, artype) 

2031 for j in range(Np + 1): 

2032 val = 0.0 

2033 for i in range(N + 1): 

2034 for k in range(i + 1): 

2035 for l in range(M - i + 1): 

2036 if k + l == j: 

2037 val += (comb(i, k) * comb(M - i, l) * b[N - i] * 

2038 pow(2 * fs, i) * (-1) ** k) 

2039 bprime[j] = real(val) 

2040 for j in range(Dp + 1): 

2041 val = 0.0 

2042 for i in range(D + 1): 

2043 for k in range(i + 1): 

2044 for l in range(M - i + 1): 

2045 if k + l == j: 

2046 val += (comb(i, k) * comb(M - i, l) * a[D - i] * 

2047 pow(2 * fs, i) * (-1) ** k) 

2048 aprime[j] = real(val) 

2049 

2050 return normalize(bprime, aprime) 

2051 

2052 

2053def _validate_gpass_gstop(gpass, gstop): 

2054 

2055 if gpass <= 0.0: 

2056 raise ValueError("gpass should be larger than 0.0") 

2057 elif gstop <= 0.0: 

2058 raise ValueError("gstop should be larger than 0.0") 

2059 elif gpass > gstop: 

2060 raise ValueError("gpass should be smaller than gstop") 

2061 

2062 

2063def iirdesign(wp, ws, gpass, gstop, analog=False, ftype='ellip', output='ba', 

2064 fs=None): 

2065 """Complete IIR digital and analog filter design. 

2066 

2067 Given passband and stopband frequencies and gains, construct an analog or 

2068 digital IIR filter of minimum order for a given basic type. Return the 

2069 output in numerator, denominator ('ba'), pole-zero ('zpk') or second order 

2070 sections ('sos') form. 

2071 

2072 Parameters 

2073 ---------- 

2074 wp, ws : float 

2075 Passband and stopband edge frequencies. 

2076 For digital filters, these are in the same units as `fs`. By default, 

2077 `fs` is 2 half-cycles/sample, so these are normalized from 0 to 1, 

2078 where 1 is the Nyquist frequency. For example: 

2079 

2080 - Lowpass: wp = 0.2, ws = 0.3 

2081 - Highpass: wp = 0.3, ws = 0.2 

2082 - Bandpass: wp = [0.2, 0.5], ws = [0.1, 0.6] 

2083 - Bandstop: wp = [0.1, 0.6], ws = [0.2, 0.5] 

2084 

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

2086 gpass : float 

2087 The maximum loss in the passband (dB). 

2088 gstop : float 

2089 The minimum attenuation in the stopband (dB). 

2090 analog : bool, optional 

2091 When True, return an analog filter, otherwise a digital filter is 

2092 returned. 

2093 ftype : str, optional 

2094 The type of IIR filter to design: 

2095 

2096 - Butterworth : 'butter' 

2097 - Chebyshev I : 'cheby1' 

2098 - Chebyshev II : 'cheby2' 

2099 - Cauer/elliptic: 'ellip' 

2100 - Bessel/Thomson: 'bessel' 

2101 

2102 output : {'ba', 'zpk', 'sos'}, optional 

2103 Type of output: numerator/denominator ('ba'), pole-zero ('zpk'), or 

2104 second-order sections ('sos'). Default is 'ba' for backwards 

2105 compatibility, but 'sos' should be used for general-purpose filtering. 

2106 fs : float, optional 

2107 The sampling frequency of the digital system. 

2108 

2109 .. versionadded:: 1.2.0 

2110 

2111 Returns 

2112 ------- 

2113 b, a : ndarray, ndarray 

2114 Numerator (`b`) and denominator (`a`) polynomials of the IIR filter. 

2115 Only returned if ``output='ba'``. 

2116 z, p, k : ndarray, ndarray, float 

2117 Zeros, poles, and system gain of the IIR filter transfer 

2118 function. Only returned if ``output='zpk'``. 

2119 sos : ndarray 

2120 Second-order sections representation of the IIR filter. 

2121 Only returned if ``output=='sos'``. 

2122 

2123 See Also 

2124 -------- 

2125 butter : Filter design using order and critical points 

2126 cheby1, cheby2, ellip, bessel 

2127 buttord : Find order and critical points from passband and stopband spec 

2128 cheb1ord, cheb2ord, ellipord 

2129 iirfilter : General filter design using order and critical frequencies 

2130 

2131 Notes 

2132 ----- 

2133 The ``'sos'`` output parameter was added in 0.16.0. 

2134 

2135 Examples 

2136 -------- 

2137 

2138 >>> from scipy import signal 

2139 >>> import matplotlib.pyplot as plt 

2140 >>> import matplotlib.ticker 

2141 

2142 >>> wp = 0.2 

2143 >>> ws = 0.3 

2144 >>> gpass = 1 

2145 >>> gstop = 40 

2146 

2147 >>> system = signal.iirdesign(wp, ws, gpass, gstop) 

2148 >>> w, h = signal.freqz(*system) 

2149 

2150 >>> fig, ax1 = plt.subplots() 

2151 >>> ax1.set_title('Digital filter frequency response') 

2152 >>> ax1.plot(w, 20 * np.log10(abs(h)), 'b') 

2153 >>> ax1.set_ylabel('Amplitude [dB]', color='b') 

2154 >>> ax1.set_xlabel('Frequency [rad/sample]') 

2155 >>> ax1.grid() 

2156 >>> ax1.set_ylim([-120, 20]) 

2157 >>> ax2 = ax1.twinx() 

2158 >>> angles = np.unwrap(np.angle(h)) 

2159 >>> ax2.plot(w, angles, 'g') 

2160 >>> ax2.set_ylabel('Angle (radians)', color='g') 

2161 >>> ax2.grid() 

2162 >>> ax2.axis('tight') 

2163 >>> ax2.set_ylim([-6, 1]) 

2164 >>> nticks = 8 

2165 >>> ax1.yaxis.set_major_locator(matplotlib.ticker.LinearLocator(nticks)) 

2166 >>> ax2.yaxis.set_major_locator(matplotlib.ticker.LinearLocator(nticks)) 

2167 

2168 """ 

2169 try: 

2170 ordfunc = filter_dict[ftype][1] 

2171 except KeyError: 

2172 raise ValueError("Invalid IIR filter type: %s" % ftype) 

2173 except IndexError: 

2174 raise ValueError(("%s does not have order selection. Use " 

2175 "iirfilter function.") % ftype) 

2176 

2177 _validate_gpass_gstop(gpass, gstop) 

2178 

2179 wp = atleast_1d(wp) 

2180 ws = atleast_1d(ws) 

2181 

2182 band_type = 2 * (len(wp) - 1) 

2183 band_type += 1 

2184 if wp[0] >= ws[0]: 

2185 band_type += 1 

2186 

2187 btype = {1: 'lowpass', 2: 'highpass', 

2188 3: 'bandstop', 4: 'bandpass'}[band_type] 

2189 

2190 N, Wn = ordfunc(wp, ws, gpass, gstop, analog=analog, fs=fs) 

2191 return iirfilter(N, Wn, rp=gpass, rs=gstop, analog=analog, btype=btype, 

2192 ftype=ftype, output=output, fs=fs) 

2193 

2194 

2195def iirfilter(N, Wn, rp=None, rs=None, btype='band', analog=False, 

2196 ftype='butter', output='ba', fs=None): 

2197 """ 

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

2199 

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

2201 coefficients. 

2202 

2203 Parameters 

2204 ---------- 

2205 N : int 

2206 The order of the filter. 

2207 Wn : array_like 

2208 A scalar or length-2 sequence giving the critical frequencies. 

2209 

2210 For digital filters, `Wn` are in the same units as `fs`. By default, 

2211 `fs` is 2 half-cycles/sample, so these are normalized from 0 to 1, 

2212 where 1 is the Nyquist frequency. (`Wn` is thus in 

2213 half-cycles / sample.) 

2214 

2215 For analog filters, `Wn` is an angular frequency (e.g., rad/s). 

2216 rp : float, optional 

2217 For Chebyshev and elliptic filters, provides the maximum ripple 

2218 in the passband. (dB) 

2219 rs : float, optional 

2220 For Chebyshev and elliptic filters, provides the minimum attenuation 

2221 in the stop band. (dB) 

2222 btype : {'bandpass', 'lowpass', 'highpass', 'bandstop'}, optional 

2223 The type of filter. Default is 'bandpass'. 

2224 analog : bool, optional 

2225 When True, return an analog filter, otherwise a digital filter is 

2226 returned. 

2227 ftype : str, optional 

2228 The type of IIR filter to design: 

2229 

2230 - Butterworth : 'butter' 

2231 - Chebyshev I : 'cheby1' 

2232 - Chebyshev II : 'cheby2' 

2233 - Cauer/elliptic: 'ellip' 

2234 - Bessel/Thomson: 'bessel' 

2235 

2236 output : {'ba', 'zpk', 'sos'}, optional 

2237 Type of output: numerator/denominator ('ba'), pole-zero ('zpk'), or 

2238 second-order sections ('sos'). Default is 'ba' for backwards 

2239 compatibility, but 'sos' should be used for general-purpose filtering. 

2240 fs : float, optional 

2241 The sampling frequency of the digital system. 

2242 

2243 .. versionadded:: 1.2.0 

2244 

2245 Returns 

2246 ------- 

2247 b, a : ndarray, ndarray 

2248 Numerator (`b`) and denominator (`a`) polynomials of the IIR filter. 

2249 Only returned if ``output='ba'``. 

2250 z, p, k : ndarray, ndarray, float 

2251 Zeros, poles, and system gain of the IIR filter transfer 

2252 function. Only returned if ``output='zpk'``. 

2253 sos : ndarray 

2254 Second-order sections representation of the IIR filter. 

2255 Only returned if ``output=='sos'``. 

2256 

2257 See Also 

2258 -------- 

2259 butter : Filter design using order and critical points 

2260 cheby1, cheby2, ellip, bessel 

2261 buttord : Find order and critical points from passband and stopband spec 

2262 cheb1ord, cheb2ord, ellipord 

2263 iirdesign : General filter design using passband and stopband spec 

2264 

2265 Notes 

2266 ----- 

2267 The ``'sos'`` output parameter was added in 0.16.0. 

2268 

2269 Examples 

2270 -------- 

2271 Generate a 17th-order Chebyshev II analog bandpass filter from 50 Hz to 

2272 200 Hz and plot the frequency response: 

2273 

2274 >>> from scipy import signal 

2275 >>> import matplotlib.pyplot as plt 

2276 

2277 >>> b, a = signal.iirfilter(17, [2*np.pi*50, 2*np.pi*200], rs=60, 

2278 ... btype='band', analog=True, ftype='cheby2') 

2279 >>> w, h = signal.freqs(b, a, 1000) 

2280 >>> fig = plt.figure() 

2281 >>> ax = fig.add_subplot(1, 1, 1) 

2282 >>> ax.semilogx(w / (2*np.pi), 20 * np.log10(np.maximum(abs(h), 1e-5))) 

2283 >>> ax.set_title('Chebyshev Type II bandpass frequency response') 

2284 >>> ax.set_xlabel('Frequency [Hz]') 

2285 >>> ax.set_ylabel('Amplitude [dB]') 

2286 >>> ax.axis((10, 1000, -100, 10)) 

2287 >>> ax.grid(which='both', axis='both') 

2288 >>> plt.show() 

2289 

2290 Create a digital filter with the same properties, in a system with 

2291 sampling rate of 2000 Hz, and plot the frequency response. (Second-order 

2292 sections implementation is required to ensure stability of a filter of 

2293 this order): 

2294 

2295 >>> sos = signal.iirfilter(17, [50, 200], rs=60, btype='band', 

2296 ... analog=False, ftype='cheby2', fs=2000, 

2297 ... output='sos') 

2298 >>> w, h = signal.sosfreqz(sos, 2000, fs=2000) 

2299 >>> fig = plt.figure() 

2300 >>> ax = fig.add_subplot(1, 1, 1) 

2301 >>> ax.semilogx(w, 20 * np.log10(np.maximum(abs(h), 1e-5))) 

2302 >>> ax.set_title('Chebyshev Type II bandpass frequency response') 

2303 >>> ax.set_xlabel('Frequency [Hz]') 

2304 >>> ax.set_ylabel('Amplitude [dB]') 

2305 >>> ax.axis((10, 1000, -100, 10)) 

2306 >>> ax.grid(which='both', axis='both') 

2307 >>> plt.show() 

2308 

2309 """ 

2310 ftype, btype, output = [x.lower() for x in (ftype, btype, output)] 

2311 Wn = asarray(Wn) 

2312 if fs is not None: 

2313 if analog: 

2314 raise ValueError("fs cannot be specified for an analog filter") 

2315 Wn = 2*Wn/fs 

2316 

2317 try: 

2318 btype = band_dict[btype] 

2319 except KeyError: 

2320 raise ValueError("'%s' is an invalid bandtype for filter." % btype) 

2321 

2322 try: 

2323 typefunc = filter_dict[ftype][0] 

2324 except KeyError: 

2325 raise ValueError("'%s' is not a valid basic IIR filter." % ftype) 

2326 

2327 if output not in ['ba', 'zpk', 'sos']: 

2328 raise ValueError("'%s' is not a valid output form." % output) 

2329 

2330 if rp is not None and rp < 0: 

2331 raise ValueError("passband ripple (rp) must be positive") 

2332 

2333 if rs is not None and rs < 0: 

2334 raise ValueError("stopband attenuation (rs) must be positive") 

2335 

2336 # Get analog lowpass prototype 

2337 if typefunc == buttap: 

2338 z, p, k = typefunc(N) 

2339 elif typefunc == besselap: 

2340 z, p, k = typefunc(N, norm=bessel_norms[ftype]) 

2341 elif typefunc == cheb1ap: 

2342 if rp is None: 

2343 raise ValueError("passband ripple (rp) must be provided to " 

2344 "design a Chebyshev I filter.") 

2345 z, p, k = typefunc(N, rp) 

2346 elif typefunc == cheb2ap: 

2347 if rs is None: 

2348 raise ValueError("stopband attenuation (rs) must be provided to " 

2349 "design an Chebyshev II filter.") 

2350 z, p, k = typefunc(N, rs) 

2351 elif typefunc == ellipap: 

2352 if rs is None or rp is None: 

2353 raise ValueError("Both rp and rs must be provided to design an " 

2354 "elliptic filter.") 

2355 z, p, k = typefunc(N, rp, rs) 

2356 else: 

2357 raise NotImplementedError("'%s' not implemented in iirfilter." % ftype) 

2358 

2359 # Pre-warp frequencies for digital filter design 

2360 if not analog: 

2361 if numpy.any(Wn <= 0) or numpy.any(Wn >= 1): 

2362 if fs is not None: 

2363 raise ValueError("Digital filter critical frequencies " 

2364 "must be 0 < Wn < fs/2 (fs={} -> fs/2={})".format(fs, fs/2)) 

2365 raise ValueError("Digital filter critical frequencies " 

2366 "must be 0 < Wn < 1") 

2367 fs = 2.0 

2368 warped = 2 * fs * tan(pi * Wn / fs) 

2369 else: 

2370 warped = Wn 

2371 

2372 # transform to lowpass, bandpass, highpass, or bandstop 

2373 if btype in ('lowpass', 'highpass'): 

2374 if numpy.size(Wn) != 1: 

2375 raise ValueError('Must specify a single critical frequency Wn for lowpass or highpass filter') 

2376 

2377 if btype == 'lowpass': 

2378 z, p, k = lp2lp_zpk(z, p, k, wo=warped) 

2379 elif btype == 'highpass': 

2380 z, p, k = lp2hp_zpk(z, p, k, wo=warped) 

2381 elif btype in ('bandpass', 'bandstop'): 

2382 try: 

2383 bw = warped[1] - warped[0] 

2384 wo = sqrt(warped[0] * warped[1]) 

2385 except IndexError: 

2386 raise ValueError('Wn must specify start and stop frequencies for bandpass or bandstop filter') 

2387 

2388 if btype == 'bandpass': 

2389 z, p, k = lp2bp_zpk(z, p, k, wo=wo, bw=bw) 

2390 elif btype == 'bandstop': 

2391 z, p, k = lp2bs_zpk(z, p, k, wo=wo, bw=bw) 

2392 else: 

2393 raise NotImplementedError("'%s' not implemented in iirfilter." % btype) 

2394 

2395 # Find discrete equivalent if necessary 

2396 if not analog: 

2397 z, p, k = bilinear_zpk(z, p, k, fs=fs) 

2398 

2399 # Transform to proper out type (pole-zero, state-space, numer-denom) 

2400 if output == 'zpk': 

2401 return z, p, k 

2402 elif output == 'ba': 

2403 return zpk2tf(z, p, k) 

2404 elif output == 'sos': 

2405 return zpk2sos(z, p, k) 

2406 

2407 

2408def _relative_degree(z, p): 

2409 """ 

2410 Return relative degree of transfer function from zeros and poles 

2411 """ 

2412 degree = len(p) - len(z) 

2413 if degree < 0: 

2414 raise ValueError("Improper transfer function. " 

2415 "Must have at least as many poles as zeros.") 

2416 else: 

2417 return degree 

2418 

2419 

2420def bilinear_zpk(z, p, k, fs): 

2421 r""" 

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

2423 

2424 Transform a set of poles and zeros from the analog s-plane to the digital 

2425 z-plane using Tustin's method, which substitutes ``(z-1) / (z+1)`` for 

2426 ``s``, maintaining the shape of the frequency response. 

2427 

2428 Parameters 

2429 ---------- 

2430 z : array_like 

2431 Zeros of the analog filter transfer function. 

2432 p : array_like 

2433 Poles of the analog filter transfer function. 

2434 k : float 

2435 System gain of the analog filter transfer function. 

2436 fs : float 

2437 Sample rate, as ordinary frequency (e.g., hertz). No prewarping is 

2438 done in this function. 

2439 

2440 Returns 

2441 ------- 

2442 z : ndarray 

2443 Zeros of the transformed digital filter transfer function. 

2444 p : ndarray 

2445 Poles of the transformed digital filter transfer function. 

2446 k : float 

2447 System gain of the transformed digital filter. 

2448 

2449 See Also 

2450 -------- 

2451 lp2lp_zpk, lp2hp_zpk, lp2bp_zpk, lp2bs_zpk 

2452 bilinear 

2453 

2454 Notes 

2455 ----- 

2456 .. versionadded:: 1.1.0 

2457 

2458 Examples 

2459 -------- 

2460 >>> from scipy import signal 

2461 >>> import matplotlib.pyplot as plt 

2462 

2463 >>> fs = 100 

2464 >>> bf = 2 * np.pi * np.array([7, 13]) 

2465 >>> filts = signal.lti(*signal.butter(4, bf, btype='bandpass', analog=True, output='zpk')) 

2466 >>> filtz = signal.lti(*signal.bilinear_zpk(filts.zeros, filts.poles, filts.gain, fs)) 

2467 >>> wz, hz = signal.freqz_zpk(filtz.zeros, filtz.poles, filtz.gain) 

2468 >>> ws, hs = signal.freqs_zpk(filts.zeros, filts.poles, filts.gain, worN=fs*wz) 

2469 >>> plt.semilogx(wz*fs/(2*np.pi), 20*np.log10(np.abs(hz).clip(1e-15)), label=r'$|H(j \omega)|$') 

2470 >>> plt.semilogx(wz*fs/(2*np.pi), 20*np.log10(np.abs(hs).clip(1e-15)), label=r'$|H_z(e^{j \omega})|$') 

2471 >>> plt.legend() 

2472 >>> plt.xlabel('Frequency [Hz]') 

2473 >>> plt.ylabel('Magnitude [dB]') 

2474 >>> plt.grid() 

2475 """ 

2476 z = atleast_1d(z) 

2477 p = atleast_1d(p) 

2478 

2479 degree = _relative_degree(z, p) 

2480 

2481 fs2 = 2.0*fs 

2482 

2483 # Bilinear transform the poles and zeros 

2484 z_z = (fs2 + z) / (fs2 - z) 

2485 p_z = (fs2 + p) / (fs2 - p) 

2486 

2487 # Any zeros that were at infinity get moved to the Nyquist frequency 

2488 z_z = append(z_z, -ones(degree)) 

2489 

2490 # Compensate for gain change 

2491 k_z = k * real(prod(fs2 - z) / prod(fs2 - p)) 

2492 

2493 return z_z, p_z, k_z 

2494 

2495 

2496def lp2lp_zpk(z, p, k, wo=1.0): 

2497 r""" 

2498 Transform a lowpass filter prototype to a different frequency. 

2499 

2500 Return an analog low-pass filter with cutoff frequency `wo` 

2501 from an analog low-pass filter prototype with unity cutoff frequency, 

2502 using zeros, poles, and gain ('zpk') representation. 

2503 

2504 Parameters 

2505 ---------- 

2506 z : array_like 

2507 Zeros of the analog filter transfer function. 

2508 p : array_like 

2509 Poles of the analog filter transfer function. 

2510 k : float 

2511 System gain of the analog filter transfer function. 

2512 wo : float 

2513 Desired cutoff, as angular frequency (e.g., rad/s). 

2514 Defaults to no change. 

2515 

2516 Returns 

2517 ------- 

2518 z : ndarray 

2519 Zeros of the transformed low-pass filter transfer function. 

2520 p : ndarray 

2521 Poles of the transformed low-pass filter transfer function. 

2522 k : float 

2523 System gain of the transformed low-pass filter. 

2524 

2525 See Also 

2526 -------- 

2527 lp2hp_zpk, lp2bp_zpk, lp2bs_zpk, bilinear 

2528 lp2lp 

2529 

2530 Notes 

2531 ----- 

2532 This is derived from the s-plane substitution 

2533 

2534 .. math:: s \rightarrow \frac{s}{\omega_0} 

2535 

2536 .. versionadded:: 1.1.0 

2537 

2538 """ 

2539 z = atleast_1d(z) 

2540 p = atleast_1d(p) 

2541 wo = float(wo) # Avoid int wraparound 

2542 

2543 degree = _relative_degree(z, p) 

2544 

2545 # Scale all points radially from origin to shift cutoff frequency 

2546 z_lp = wo * z 

2547 p_lp = wo * p 

2548 

2549 # Each shifted pole decreases gain by wo, each shifted zero increases it. 

2550 # Cancel out the net change to keep overall gain the same 

2551 k_lp = k * wo**degree 

2552 

2553 return z_lp, p_lp, k_lp 

2554 

2555 

2556def lp2hp_zpk(z, p, k, wo=1.0): 

2557 r""" 

2558 Transform a lowpass filter prototype to a highpass filter. 

2559 

2560 Return an analog high-pass filter with cutoff frequency `wo` 

2561 from an analog low-pass filter prototype with unity cutoff frequency, 

2562 using zeros, poles, and gain ('zpk') representation. 

2563 

2564 Parameters 

2565 ---------- 

2566 z : array_like 

2567 Zeros of the analog filter transfer function. 

2568 p : array_like 

2569 Poles of the analog filter transfer function. 

2570 k : float 

2571 System gain of the analog filter transfer function. 

2572 wo : float 

2573 Desired cutoff, as angular frequency (e.g., rad/s). 

2574 Defaults to no change. 

2575 

2576 Returns 

2577 ------- 

2578 z : ndarray 

2579 Zeros of the transformed high-pass filter transfer function. 

2580 p : ndarray 

2581 Poles of the transformed high-pass filter transfer function. 

2582 k : float 

2583 System gain of the transformed high-pass filter. 

2584 

2585 See Also 

2586 -------- 

2587 lp2lp_zpk, lp2bp_zpk, lp2bs_zpk, bilinear 

2588 lp2hp 

2589 

2590 Notes 

2591 ----- 

2592 This is derived from the s-plane substitution 

2593 

2594 .. math:: s \rightarrow \frac{\omega_0}{s} 

2595 

2596 This maintains symmetry of the lowpass and highpass responses on a 

2597 logarithmic scale. 

2598 

2599 .. versionadded:: 1.1.0 

2600 

2601 """ 

2602 z = atleast_1d(z) 

2603 p = atleast_1d(p) 

2604 wo = float(wo) 

2605 

2606 degree = _relative_degree(z, p) 

2607 

2608 # Invert positions radially about unit circle to convert LPF to HPF 

2609 # Scale all points radially from origin to shift cutoff frequency 

2610 z_hp = wo / z 

2611 p_hp = wo / p 

2612 

2613 # If lowpass had zeros at infinity, inverting moves them to origin. 

2614 z_hp = append(z_hp, zeros(degree)) 

2615 

2616 # Cancel out gain change caused by inversion 

2617 k_hp = k * real(prod(-z) / prod(-p)) 

2618 

2619 return z_hp, p_hp, k_hp 

2620 

2621 

2622def lp2bp_zpk(z, p, k, wo=1.0, bw=1.0): 

2623 r""" 

2624 Transform a lowpass filter prototype to a bandpass filter. 

2625 

2626 Return an analog band-pass filter with center frequency `wo` and 

2627 bandwidth `bw` from an analog low-pass filter prototype with unity 

2628 cutoff frequency, using zeros, poles, and gain ('zpk') representation. 

2629 

2630 Parameters 

2631 ---------- 

2632 z : array_like 

2633 Zeros of the analog filter transfer function. 

2634 p : array_like 

2635 Poles of the analog filter transfer function. 

2636 k : float 

2637 System gain of the analog filter transfer function. 

2638 wo : float 

2639 Desired passband center, as angular frequency (e.g., rad/s). 

2640 Defaults to no change. 

2641 bw : float 

2642 Desired passband width, as angular frequency (e.g., rad/s). 

2643 Defaults to 1. 

2644 

2645 Returns 

2646 ------- 

2647 z : ndarray 

2648 Zeros of the transformed band-pass filter transfer function. 

2649 p : ndarray 

2650 Poles of the transformed band-pass filter transfer function. 

2651 k : float 

2652 System gain of the transformed band-pass filter. 

2653 

2654 See Also 

2655 -------- 

2656 lp2lp_zpk, lp2hp_zpk, lp2bs_zpk, bilinear 

2657 lp2bp 

2658 

2659 Notes 

2660 ----- 

2661 This is derived from the s-plane substitution 

2662 

2663 .. math:: s \rightarrow \frac{s^2 + {\omega_0}^2}{s \cdot \mathrm{BW}} 

2664 

2665 This is the "wideband" transformation, producing a passband with 

2666 geometric (log frequency) symmetry about `wo`. 

2667 

2668 .. versionadded:: 1.1.0 

2669 

2670 """ 

2671 z = atleast_1d(z) 

2672 p = atleast_1d(p) 

2673 wo = float(wo) 

2674 bw = float(bw) 

2675 

2676 degree = _relative_degree(z, p) 

2677 

2678 # Scale poles and zeros to desired bandwidth 

2679 z_lp = z * bw/2 

2680 p_lp = p * bw/2 

2681 

2682 # Square root needs to produce complex result, not NaN 

2683 z_lp = z_lp.astype(complex) 

2684 p_lp = p_lp.astype(complex) 

2685 

2686 # Duplicate poles and zeros and shift from baseband to +wo and -wo 

2687 z_bp = concatenate((z_lp + sqrt(z_lp**2 - wo**2), 

2688 z_lp - sqrt(z_lp**2 - wo**2))) 

2689 p_bp = concatenate((p_lp + sqrt(p_lp**2 - wo**2), 

2690 p_lp - sqrt(p_lp**2 - wo**2))) 

2691 

2692 # Move degree zeros to origin, leaving degree zeros at infinity for BPF 

2693 z_bp = append(z_bp, zeros(degree)) 

2694 

2695 # Cancel out gain change from frequency scaling 

2696 k_bp = k * bw**degree 

2697 

2698 return z_bp, p_bp, k_bp 

2699 

2700 

2701def lp2bs_zpk(z, p, k, wo=1.0, bw=1.0): 

2702 r""" 

2703 Transform a lowpass filter prototype to a bandstop filter. 

2704 

2705 Return an analog band-stop filter with center frequency `wo` and 

2706 stopband width `bw` from an analog low-pass filter prototype with unity 

2707 cutoff frequency, using zeros, poles, and gain ('zpk') representation. 

2708 

2709 Parameters 

2710 ---------- 

2711 z : array_like 

2712 Zeros of the analog filter transfer function. 

2713 p : array_like 

2714 Poles of the analog filter transfer function. 

2715 k : float 

2716 System gain of the analog filter transfer function. 

2717 wo : float 

2718 Desired stopband center, as angular frequency (e.g., rad/s). 

2719 Defaults to no change. 

2720 bw : float 

2721 Desired stopband width, as angular frequency (e.g., rad/s). 

2722 Defaults to 1. 

2723 

2724 Returns 

2725 ------- 

2726 z : ndarray 

2727 Zeros of the transformed band-stop filter transfer function. 

2728 p : ndarray 

2729 Poles of the transformed band-stop filter transfer function. 

2730 k : float 

2731 System gain of the transformed band-stop filter. 

2732 

2733 See Also 

2734 -------- 

2735 lp2lp_zpk, lp2hp_zpk, lp2bp_zpk, bilinear 

2736 lp2bs 

2737 

2738 Notes 

2739 ----- 

2740 This is derived from the s-plane substitution 

2741 

2742 .. math:: s \rightarrow \frac{s \cdot \mathrm{BW}}{s^2 + {\omega_0}^2} 

2743 

2744 This is the "wideband" transformation, producing a stopband with 

2745 geometric (log frequency) symmetry about `wo`. 

2746 

2747 .. versionadded:: 1.1.0 

2748 

2749 """ 

2750 z = atleast_1d(z) 

2751 p = atleast_1d(p) 

2752 wo = float(wo) 

2753 bw = float(bw) 

2754 

2755 degree = _relative_degree(z, p) 

2756 

2757 # Invert to a highpass filter with desired bandwidth 

2758 z_hp = (bw/2) / z 

2759 p_hp = (bw/2) / p 

2760 

2761 # Square root needs to produce complex result, not NaN 

2762 z_hp = z_hp.astype(complex) 

2763 p_hp = p_hp.astype(complex) 

2764 

2765 # Duplicate poles and zeros and shift from baseband to +wo and -wo 

2766 z_bs = concatenate((z_hp + sqrt(z_hp**2 - wo**2), 

2767 z_hp - sqrt(z_hp**2 - wo**2))) 

2768 p_bs = concatenate((p_hp + sqrt(p_hp**2 - wo**2), 

2769 p_hp - sqrt(p_hp**2 - wo**2))) 

2770 

2771 # Move any zeros that were at infinity to the center of the stopband 

2772 z_bs = append(z_bs, full(degree, +1j*wo)) 

2773 z_bs = append(z_bs, full(degree, -1j*wo)) 

2774 

2775 # Cancel out gain change caused by inversion 

2776 k_bs = k * real(prod(-z) / prod(-p)) 

2777 

2778 return z_bs, p_bs, k_bs 

2779 

2780 

2781def butter(N, Wn, btype='low', analog=False, output='ba', fs=None): 

2782 """ 

2783 Butterworth digital and analog filter design. 

2784 

2785 Design an Nth-order digital or analog Butterworth filter and return 

2786 the filter coefficients. 

2787 

2788 Parameters 

2789 ---------- 

2790 N : int 

2791 The order of the filter. 

2792 Wn : array_like 

2793 The critical frequency or frequencies. For lowpass and highpass 

2794 filters, Wn is a scalar; for bandpass and bandstop filters, 

2795 Wn is a length-2 sequence. 

2796 

2797 For a Butterworth filter, this is the point at which the gain 

2798 drops to 1/sqrt(2) that of the passband (the "-3 dB point"). 

2799 

2800 For digital filters, `Wn` are in the same units as `fs`. By default, 

2801 `fs` is 2 half-cycles/sample, so these are normalized from 0 to 1, 

2802 where 1 is the Nyquist frequency. (`Wn` is thus in 

2803 half-cycles / sample.) 

2804 

2805 For analog filters, `Wn` is an angular frequency (e.g. rad/s). 

2806 btype : {'lowpass', 'highpass', 'bandpass', 'bandstop'}, optional 

2807 The type of filter. Default is 'lowpass'. 

2808 analog : bool, optional 

2809 When True, return an analog filter, otherwise a digital filter is 

2810 returned. 

2811 output : {'ba', 'zpk', 'sos'}, optional 

2812 Type of output: numerator/denominator ('ba'), pole-zero ('zpk'), or 

2813 second-order sections ('sos'). Default is 'ba' for backwards 

2814 compatibility, but 'sos' should be used for general-purpose filtering. 

2815 fs : float, optional 

2816 The sampling frequency of the digital system. 

2817 

2818 .. versionadded:: 1.2.0 

2819 

2820 Returns 

2821 ------- 

2822 b, a : ndarray, ndarray 

2823 Numerator (`b`) and denominator (`a`) polynomials of the IIR filter. 

2824 Only returned if ``output='ba'``. 

2825 z, p, k : ndarray, ndarray, float 

2826 Zeros, poles, and system gain of the IIR filter transfer 

2827 function. Only returned if ``output='zpk'``. 

2828 sos : ndarray 

2829 Second-order sections representation of the IIR filter. 

2830 Only returned if ``output=='sos'``. 

2831 

2832 See Also 

2833 -------- 

2834 buttord, buttap 

2835 

2836 Notes 

2837 ----- 

2838 The Butterworth filter has maximally flat frequency response in the 

2839 passband. 

2840 

2841 The ``'sos'`` output parameter was added in 0.16.0. 

2842 

2843 If the transfer function form ``[b, a]`` is requested, numerical 

2844 problems can occur since the conversion between roots and 

2845 the polynomial coefficients is a numerically sensitive operation, 

2846 even for N >= 4. It is recommended to work with the SOS 

2847 representation. 

2848 

2849 Examples 

2850 -------- 

2851 Design an analog filter and plot its frequency response, showing the 

2852 critical points: 

2853 

2854 >>> from scipy import signal 

2855 >>> import matplotlib.pyplot as plt 

2856 

2857 >>> b, a = signal.butter(4, 100, 'low', analog=True) 

2858 >>> w, h = signal.freqs(b, a) 

2859 >>> plt.semilogx(w, 20 * np.log10(abs(h))) 

2860 >>> plt.title('Butterworth filter frequency response') 

2861 >>> plt.xlabel('Frequency [radians / second]') 

2862 >>> plt.ylabel('Amplitude [dB]') 

2863 >>> plt.margins(0, 0.1) 

2864 >>> plt.grid(which='both', axis='both') 

2865 >>> plt.axvline(100, color='green') # cutoff frequency 

2866 >>> plt.show() 

2867 

2868 Generate a signal made up of 10 Hz and 20 Hz, sampled at 1 kHz 

2869 

2870 >>> t = np.linspace(0, 1, 1000, False) # 1 second 

2871 >>> sig = np.sin(2*np.pi*10*t) + np.sin(2*np.pi*20*t) 

2872 >>> fig, (ax1, ax2) = plt.subplots(2, 1, sharex=True) 

2873 >>> ax1.plot(t, sig) 

2874 >>> ax1.set_title('10 Hz and 20 Hz sinusoids') 

2875 >>> ax1.axis([0, 1, -2, 2]) 

2876 

2877 Design a digital high-pass filter at 15 Hz to remove the 10 Hz tone, and 

2878 apply it to the signal. (It's recommended to use second-order sections 

2879 format when filtering, to avoid numerical error with transfer function 

2880 (``ba``) format): 

2881 

2882 >>> sos = signal.butter(10, 15, 'hp', fs=1000, output='sos') 

2883 >>> filtered = signal.sosfilt(sos, sig) 

2884 >>> ax2.plot(t, filtered) 

2885 >>> ax2.set_title('After 15 Hz high-pass filter') 

2886 >>> ax2.axis([0, 1, -2, 2]) 

2887 >>> ax2.set_xlabel('Time [seconds]') 

2888 >>> plt.tight_layout() 

2889 >>> plt.show() 

2890 """ 

2891 return iirfilter(N, Wn, btype=btype, analog=analog, 

2892 output=output, ftype='butter', fs=fs) 

2893 

2894 

2895def cheby1(N, rp, Wn, btype='low', analog=False, output='ba', fs=None): 

2896 """ 

2897 Chebyshev type I digital and analog filter design. 

2898 

2899 Design an Nth-order digital or analog Chebyshev type I filter and 

2900 return the filter coefficients. 

2901 

2902 Parameters 

2903 ---------- 

2904 N : int 

2905 The order of the filter. 

2906 rp : float 

2907 The maximum ripple allowed below unity gain in the passband. 

2908 Specified in decibels, as a positive number. 

2909 Wn : array_like 

2910 A scalar or length-2 sequence giving the critical frequencies. 

2911 For Type I filters, this is the point in the transition band at which 

2912 the gain first drops below -`rp`. 

2913 

2914 For digital filters, `Wn` are in the same units as `fs`. By default, 

2915 `fs` is 2 half-cycles/sample, so these are normalized from 0 to 1, 

2916 where 1 is the Nyquist frequency. (`Wn` is thus in 

2917 half-cycles / sample.) 

2918 

2919 For analog filters, `Wn` is an angular frequency (e.g., rad/s). 

2920 btype : {'lowpass', 'highpass', 'bandpass', 'bandstop'}, optional 

2921 The type of filter. Default is 'lowpass'. 

2922 analog : bool, optional 

2923 When True, return an analog filter, otherwise a digital filter is 

2924 returned. 

2925 output : {'ba', 'zpk', 'sos'}, optional 

2926 Type of output: numerator/denominator ('ba'), pole-zero ('zpk'), or 

2927 second-order sections ('sos'). Default is 'ba' for backwards 

2928 compatibility, but 'sos' should be used for general-purpose filtering. 

2929 fs : float, optional 

2930 The sampling frequency of the digital system. 

2931 

2932 .. versionadded:: 1.2.0 

2933 

2934 Returns 

2935 ------- 

2936 b, a : ndarray, ndarray 

2937 Numerator (`b`) and denominator (`a`) polynomials of the IIR filter. 

2938 Only returned if ``output='ba'``. 

2939 z, p, k : ndarray, ndarray, float 

2940 Zeros, poles, and system gain of the IIR filter transfer 

2941 function. Only returned if ``output='zpk'``. 

2942 sos : ndarray 

2943 Second-order sections representation of the IIR filter. 

2944 Only returned if ``output=='sos'``. 

2945 

2946 See Also 

2947 -------- 

2948 cheb1ord, cheb1ap 

2949 

2950 Notes 

2951 ----- 

2952 The Chebyshev type I filter maximizes the rate of cutoff between the 

2953 frequency response's passband and stopband, at the expense of ripple in 

2954 the passband and increased ringing in the step response. 

2955 

2956 Type I filters roll off faster than Type II (`cheby2`), but Type II 

2957 filters do not have any ripple in the passband. 

2958 

2959 The equiripple passband has N maxima or minima (for example, a 

2960 5th-order filter has 3 maxima and 2 minima). Consequently, the DC gain is 

2961 unity for odd-order filters, or -rp dB for even-order filters. 

2962 

2963 The ``'sos'`` output parameter was added in 0.16.0. 

2964 

2965 Examples 

2966 -------- 

2967 Design an analog filter and plot its frequency response, showing the 

2968 critical points: 

2969 

2970 >>> from scipy import signal 

2971 >>> import matplotlib.pyplot as plt 

2972 

2973 >>> b, a = signal.cheby1(4, 5, 100, 'low', analog=True) 

2974 >>> w, h = signal.freqs(b, a) 

2975 >>> plt.semilogx(w, 20 * np.log10(abs(h))) 

2976 >>> plt.title('Chebyshev Type I frequency response (rp=5)') 

2977 >>> plt.xlabel('Frequency [radians / second]') 

2978 >>> plt.ylabel('Amplitude [dB]') 

2979 >>> plt.margins(0, 0.1) 

2980 >>> plt.grid(which='both', axis='both') 

2981 >>> plt.axvline(100, color='green') # cutoff frequency 

2982 >>> plt.axhline(-5, color='green') # rp 

2983 >>> plt.show() 

2984 

2985 Generate a signal made up of 10 Hz and 20 Hz, sampled at 1 kHz 

2986 

2987 >>> t = np.linspace(0, 1, 1000, False) # 1 second 

2988 >>> sig = np.sin(2*np.pi*10*t) + np.sin(2*np.pi*20*t) 

2989 >>> fig, (ax1, ax2) = plt.subplots(2, 1, sharex=True) 

2990 >>> ax1.plot(t, sig) 

2991 >>> ax1.set_title('10 Hz and 20 Hz sinusoids') 

2992 >>> ax1.axis([0, 1, -2, 2]) 

2993 

2994 Design a digital high-pass filter at 15 Hz to remove the 10 Hz tone, and 

2995 apply it to the signal. (It's recommended to use second-order sections 

2996 format when filtering, to avoid numerical error with transfer function 

2997 (``ba``) format): 

2998 

2999 >>> sos = signal.cheby1(10, 1, 15, 'hp', fs=1000, output='sos') 

3000 >>> filtered = signal.sosfilt(sos, sig) 

3001 >>> ax2.plot(t, filtered) 

3002 >>> ax2.set_title('After 15 Hz high-pass filter') 

3003 >>> ax2.axis([0, 1, -2, 2]) 

3004 >>> ax2.set_xlabel('Time [seconds]') 

3005 >>> plt.tight_layout() 

3006 >>> plt.show() 

3007 """ 

3008 return iirfilter(N, Wn, rp=rp, btype=btype, analog=analog, 

3009 output=output, ftype='cheby1', fs=fs) 

3010 

3011 

3012def cheby2(N, rs, Wn, btype='low', analog=False, output='ba', fs=None): 

3013 """ 

3014 Chebyshev type II digital and analog filter design. 

3015 

3016 Design an Nth-order digital or analog Chebyshev type II filter and 

3017 return the filter coefficients. 

3018 

3019 Parameters 

3020 ---------- 

3021 N : int 

3022 The order of the filter. 

3023 rs : float 

3024 The minimum attenuation required in the stop band. 

3025 Specified in decibels, as a positive number. 

3026 Wn : array_like 

3027 A scalar or length-2 sequence giving the critical frequencies. 

3028 For Type II filters, this is the point in the transition band at which 

3029 the gain first reaches -`rs`. 

3030 

3031 For digital filters, `Wn` are in the same units as `fs`. By default, 

3032 `fs` is 2 half-cycles/sample, so these are normalized from 0 to 1, 

3033 where 1 is the Nyquist frequency. (`Wn` is thus in 

3034 half-cycles / sample.) 

3035 

3036 For analog filters, `Wn` is an angular frequency (e.g., rad/s). 

3037 btype : {'lowpass', 'highpass', 'bandpass', 'bandstop'}, optional 

3038 The type of filter. Default is 'lowpass'. 

3039 analog : bool, optional 

3040 When True, return an analog filter, otherwise a digital filter is 

3041 returned. 

3042 output : {'ba', 'zpk', 'sos'}, optional 

3043 Type of output: numerator/denominator ('ba'), pole-zero ('zpk'), or 

3044 second-order sections ('sos'). Default is 'ba' for backwards 

3045 compatibility, but 'sos' should be used for general-purpose filtering. 

3046 fs : float, optional 

3047 The sampling frequency of the digital system. 

3048 

3049 .. versionadded:: 1.2.0 

3050 

3051 Returns 

3052 ------- 

3053 b, a : ndarray, ndarray 

3054 Numerator (`b`) and denominator (`a`) polynomials of the IIR filter. 

3055 Only returned if ``output='ba'``. 

3056 z, p, k : ndarray, ndarray, float 

3057 Zeros, poles, and system gain of the IIR filter transfer 

3058 function. Only returned if ``output='zpk'``. 

3059 sos : ndarray 

3060 Second-order sections representation of the IIR filter. 

3061 Only returned if ``output=='sos'``. 

3062 

3063 See Also 

3064 -------- 

3065 cheb2ord, cheb2ap 

3066 

3067 Notes 

3068 ----- 

3069 The Chebyshev type II filter maximizes the rate of cutoff between the 

3070 frequency response's passband and stopband, at the expense of ripple in 

3071 the stopband and increased ringing in the step response. 

3072 

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

3074 

3075 The ``'sos'`` output parameter was added in 0.16.0. 

3076 

3077 Examples 

3078 -------- 

3079 Design an analog filter and plot its frequency response, showing the 

3080 critical points: 

3081 

3082 >>> from scipy import signal 

3083 >>> import matplotlib.pyplot as plt 

3084 

3085 >>> b, a = signal.cheby2(4, 40, 100, 'low', analog=True) 

3086 >>> w, h = signal.freqs(b, a) 

3087 >>> plt.semilogx(w, 20 * np.log10(abs(h))) 

3088 >>> plt.title('Chebyshev Type II frequency response (rs=40)') 

3089 >>> plt.xlabel('Frequency [radians / second]') 

3090 >>> plt.ylabel('Amplitude [dB]') 

3091 >>> plt.margins(0, 0.1) 

3092 >>> plt.grid(which='both', axis='both') 

3093 >>> plt.axvline(100, color='green') # cutoff frequency 

3094 >>> plt.axhline(-40, color='green') # rs 

3095 >>> plt.show() 

3096 

3097 Generate a signal made up of 10 Hz and 20 Hz, sampled at 1 kHz 

3098 

3099 >>> t = np.linspace(0, 1, 1000, False) # 1 second 

3100 >>> sig = np.sin(2*np.pi*10*t) + np.sin(2*np.pi*20*t) 

3101 >>> fig, (ax1, ax2) = plt.subplots(2, 1, sharex=True) 

3102 >>> ax1.plot(t, sig) 

3103 >>> ax1.set_title('10 Hz and 20 Hz sinusoids') 

3104 >>> ax1.axis([0, 1, -2, 2]) 

3105 

3106 Design a digital high-pass filter at 17 Hz to remove the 10 Hz tone, and 

3107 apply it to the signal. (It's recommended to use second-order sections 

3108 format when filtering, to avoid numerical error with transfer function 

3109 (``ba``) format): 

3110 

3111 >>> sos = signal.cheby2(12, 20, 17, 'hp', fs=1000, output='sos') 

3112 >>> filtered = signal.sosfilt(sos, sig) 

3113 >>> ax2.plot(t, filtered) 

3114 >>> ax2.set_title('After 17 Hz high-pass filter') 

3115 >>> ax2.axis([0, 1, -2, 2]) 

3116 >>> ax2.set_xlabel('Time [seconds]') 

3117 >>> plt.show() 

3118 """ 

3119 return iirfilter(N, Wn, rs=rs, btype=btype, analog=analog, 

3120 output=output, ftype='cheby2', fs=fs) 

3121 

3122 

3123def ellip(N, rp, rs, Wn, btype='low', analog=False, output='ba', fs=None): 

3124 """ 

3125 Elliptic (Cauer) digital and analog filter design. 

3126 

3127 Design an Nth-order digital or analog elliptic filter and return 

3128 the filter coefficients. 

3129 

3130 Parameters 

3131 ---------- 

3132 N : int 

3133 The order of the filter. 

3134 rp : float 

3135 The maximum ripple allowed below unity gain in the passband. 

3136 Specified in decibels, as a positive number. 

3137 rs : float 

3138 The minimum attenuation required in the stop band. 

3139 Specified in decibels, as a positive number. 

3140 Wn : array_like 

3141 A scalar or length-2 sequence giving the critical frequencies. 

3142 For elliptic filters, this is the point in the transition band at 

3143 which the gain first drops below -`rp`. 

3144 

3145 For digital filters, `Wn` are in the same units as `fs`. By default, 

3146 `fs` is 2 half-cycles/sample, so these are normalized from 0 to 1, 

3147 where 1 is the Nyquist frequency. (`Wn` is thus in 

3148 half-cycles / sample.) 

3149 

3150 For analog filters, `Wn` is an angular frequency (e.g., rad/s). 

3151 btype : {'lowpass', 'highpass', 'bandpass', 'bandstop'}, optional 

3152 The type of filter. Default is 'lowpass'. 

3153 analog : bool, optional 

3154 When True, return an analog filter, otherwise a digital filter is 

3155 returned. 

3156 output : {'ba', 'zpk', 'sos'}, optional 

3157 Type of output: numerator/denominator ('ba'), pole-zero ('zpk'), or 

3158 second-order sections ('sos'). Default is 'ba' for backwards 

3159 compatibility, but 'sos' should be used for general-purpose filtering. 

3160 fs : float, optional 

3161 The sampling frequency of the digital system. 

3162 

3163 .. versionadded:: 1.2.0 

3164 

3165 Returns 

3166 ------- 

3167 b, a : ndarray, ndarray 

3168 Numerator (`b`) and denominator (`a`) polynomials of the IIR filter. 

3169 Only returned if ``output='ba'``. 

3170 z, p, k : ndarray, ndarray, float 

3171 Zeros, poles, and system gain of the IIR filter transfer 

3172 function. Only returned if ``output='zpk'``. 

3173 sos : ndarray 

3174 Second-order sections representation of the IIR filter. 

3175 Only returned if ``output=='sos'``. 

3176 

3177 See Also 

3178 -------- 

3179 ellipord, ellipap 

3180 

3181 Notes 

3182 ----- 

3183 Also known as Cauer or Zolotarev filters, the elliptical filter maximizes 

3184 the rate of transition between the frequency response's passband and 

3185 stopband, at the expense of ripple in both, and increased ringing in the 

3186 step response. 

3187 

3188 As `rp` approaches 0, the elliptical filter becomes a Chebyshev 

3189 type II filter (`cheby2`). As `rs` approaches 0, it becomes a Chebyshev 

3190 type I filter (`cheby1`). As both approach 0, it becomes a Butterworth 

3191 filter (`butter`). 

3192 

3193 The equiripple passband has N maxima or minima (for example, a 

3194 5th-order filter has 3 maxima and 2 minima). Consequently, the DC gain is 

3195 unity for odd-order filters, or -rp dB for even-order filters. 

3196 

3197 The ``'sos'`` output parameter was added in 0.16.0. 

3198 

3199 Examples 

3200 -------- 

3201 Design an analog filter and plot its frequency response, showing the 

3202 critical points: 

3203 

3204 >>> from scipy import signal 

3205 >>> import matplotlib.pyplot as plt 

3206 

3207 >>> b, a = signal.ellip(4, 5, 40, 100, 'low', analog=True) 

3208 >>> w, h = signal.freqs(b, a) 

3209 >>> plt.semilogx(w, 20 * np.log10(abs(h))) 

3210 >>> plt.title('Elliptic filter frequency response (rp=5, rs=40)') 

3211 >>> plt.xlabel('Frequency [radians / second]') 

3212 >>> plt.ylabel('Amplitude [dB]') 

3213 >>> plt.margins(0, 0.1) 

3214 >>> plt.grid(which='both', axis='both') 

3215 >>> plt.axvline(100, color='green') # cutoff frequency 

3216 >>> plt.axhline(-40, color='green') # rs 

3217 >>> plt.axhline(-5, color='green') # rp 

3218 >>> plt.show() 

3219 

3220 Generate a signal made up of 10 Hz and 20 Hz, sampled at 1 kHz 

3221 

3222 >>> t = np.linspace(0, 1, 1000, False) # 1 second 

3223 >>> sig = np.sin(2*np.pi*10*t) + np.sin(2*np.pi*20*t) 

3224 >>> fig, (ax1, ax2) = plt.subplots(2, 1, sharex=True) 

3225 >>> ax1.plot(t, sig) 

3226 >>> ax1.set_title('10 Hz and 20 Hz sinusoids') 

3227 >>> ax1.axis([0, 1, -2, 2]) 

3228 

3229 Design a digital high-pass filter at 17 Hz to remove the 10 Hz tone, and 

3230 apply it to the signal. (It's recommended to use second-order sections 

3231 format when filtering, to avoid numerical error with transfer function 

3232 (``ba``) format): 

3233 

3234 >>> sos = signal.ellip(8, 1, 100, 17, 'hp', fs=1000, output='sos') 

3235 >>> filtered = signal.sosfilt(sos, sig) 

3236 >>> ax2.plot(t, filtered) 

3237 >>> ax2.set_title('After 17 Hz high-pass filter') 

3238 >>> ax2.axis([0, 1, -2, 2]) 

3239 >>> ax2.set_xlabel('Time [seconds]') 

3240 >>> plt.tight_layout() 

3241 >>> plt.show() 

3242 """ 

3243 return iirfilter(N, Wn, rs=rs, rp=rp, btype=btype, analog=analog, 

3244 output=output, ftype='elliptic', fs=fs) 

3245 

3246 

3247def bessel(N, Wn, btype='low', analog=False, output='ba', norm='phase', 

3248 fs=None): 

3249 """ 

3250 Bessel/Thomson digital and analog filter design. 

3251 

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

3253 filter coefficients. 

3254 

3255 Parameters 

3256 ---------- 

3257 N : int 

3258 The order of the filter. 

3259 Wn : array_like 

3260 A scalar or length-2 sequence giving the critical frequencies (defined 

3261 by the `norm` parameter). 

3262 For analog filters, `Wn` is an angular frequency (e.g., rad/s). 

3263 

3264 For digital filters, `Wn` are in the same units as `fs`. By default, 

3265 `fs` is 2 half-cycles/sample, so these are normalized from 0 to 1, 

3266 where 1 is the Nyquist frequency. (`Wn` is thus in 

3267 half-cycles / sample.) 

3268 btype : {'lowpass', 'highpass', 'bandpass', 'bandstop'}, optional 

3269 The type of filter. Default is 'lowpass'. 

3270 analog : bool, optional 

3271 When True, return an analog filter, otherwise a digital filter is 

3272 returned. (See Notes.) 

3273 output : {'ba', 'zpk', 'sos'}, optional 

3274 Type of output: numerator/denominator ('ba'), pole-zero ('zpk'), or 

3275 second-order sections ('sos'). Default is 'ba'. 

3276 norm : {'phase', 'delay', 'mag'}, optional 

3277 Critical frequency normalization: 

3278 

3279 ``phase`` 

3280 The filter is normalized such that the phase response reaches its 

3281 midpoint at angular (e.g. rad/s) frequency `Wn`. This happens for 

3282 both low-pass and high-pass filters, so this is the 

3283 "phase-matched" case. 

3284 

3285 The magnitude response asymptotes are the same as a Butterworth 

3286 filter of the same order with a cutoff of `Wn`. 

3287 

3288 This is the default, and matches MATLAB's implementation. 

3289 

3290 ``delay`` 

3291 The filter is normalized such that the group delay in the passband 

3292 is 1/`Wn` (e.g., seconds). This is the "natural" type obtained by 

3293 solving Bessel polynomials. 

3294 

3295 ``mag`` 

3296 The filter is normalized such that the gain magnitude is -3 dB at 

3297 angular frequency `Wn`. 

3298 

3299 .. versionadded:: 0.18.0 

3300 fs : float, optional 

3301 The sampling frequency of the digital system. 

3302 

3303 .. versionadded:: 1.2.0 

3304 

3305 Returns 

3306 ------- 

3307 b, a : ndarray, ndarray 

3308 Numerator (`b`) and denominator (`a`) polynomials of the IIR filter. 

3309 Only returned if ``output='ba'``. 

3310 z, p, k : ndarray, ndarray, float 

3311 Zeros, poles, and system gain of the IIR filter transfer 

3312 function. Only returned if ``output='zpk'``. 

3313 sos : ndarray 

3314 Second-order sections representation of the IIR filter. 

3315 Only returned if ``output=='sos'``. 

3316 

3317 Notes 

3318 ----- 

3319 Also known as a Thomson filter, the analog Bessel filter has maximally 

3320 flat group delay and maximally linear phase response, with very little 

3321 ringing in the step response. [1]_ 

3322 

3323 The Bessel is inherently an analog filter. This function generates digital 

3324 Bessel filters using the bilinear transform, which does not preserve the 

3325 phase response of the analog filter. As such, it is only approximately 

3326 correct at frequencies below about fs/4. To get maximally-flat group 

3327 delay at higher frequencies, the analog Bessel filter must be transformed 

3328 using phase-preserving techniques. 

3329 

3330 See `besselap` for implementation details and references. 

3331 

3332 The ``'sos'`` output parameter was added in 0.16.0. 

3333 

3334 Examples 

3335 -------- 

3336 Plot the phase-normalized frequency response, showing the relationship 

3337 to the Butterworth's cutoff frequency (green): 

3338 

3339 >>> from scipy import signal 

3340 >>> import matplotlib.pyplot as plt 

3341 

3342 >>> b, a = signal.butter(4, 100, 'low', analog=True) 

3343 >>> w, h = signal.freqs(b, a) 

3344 >>> plt.semilogx(w, 20 * np.log10(np.abs(h)), color='silver', ls='dashed') 

3345 >>> b, a = signal.bessel(4, 100, 'low', analog=True, norm='phase') 

3346 >>> w, h = signal.freqs(b, a) 

3347 >>> plt.semilogx(w, 20 * np.log10(np.abs(h))) 

3348 >>> plt.title('Bessel filter magnitude response (with Butterworth)') 

3349 >>> plt.xlabel('Frequency [radians / second]') 

3350 >>> plt.ylabel('Amplitude [dB]') 

3351 >>> plt.margins(0, 0.1) 

3352 >>> plt.grid(which='both', axis='both') 

3353 >>> plt.axvline(100, color='green') # cutoff frequency 

3354 >>> plt.show() 

3355 

3356 and the phase midpoint: 

3357 

3358 >>> plt.figure() 

3359 >>> plt.semilogx(w, np.unwrap(np.angle(h))) 

3360 >>> plt.axvline(100, color='green') # cutoff frequency 

3361 >>> plt.axhline(-np.pi, color='red') # phase midpoint 

3362 >>> plt.title('Bessel filter phase response') 

3363 >>> plt.xlabel('Frequency [radians / second]') 

3364 >>> plt.ylabel('Phase [radians]') 

3365 >>> plt.margins(0, 0.1) 

3366 >>> plt.grid(which='both', axis='both') 

3367 >>> plt.show() 

3368 

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

3370 

3371 >>> b, a = signal.bessel(3, 10, 'low', analog=True, norm='mag') 

3372 >>> w, h = signal.freqs(b, a) 

3373 >>> plt.semilogx(w, 20 * np.log10(np.abs(h))) 

3374 >>> plt.axhline(-3, color='red') # -3 dB magnitude 

3375 >>> plt.axvline(10, color='green') # cutoff frequency 

3376 >>> plt.title('Magnitude-normalized Bessel filter frequency response') 

3377 >>> plt.xlabel('Frequency [radians / second]') 

3378 >>> plt.ylabel('Amplitude [dB]') 

3379 >>> plt.margins(0, 0.1) 

3380 >>> plt.grid(which='both', axis='both') 

3381 >>> plt.show() 

3382 

3383 Plot the delay-normalized filter, showing the maximally-flat group delay 

3384 at 0.1 seconds: 

3385 

3386 >>> b, a = signal.bessel(5, 1/0.1, 'low', analog=True, norm='delay') 

3387 >>> w, h = signal.freqs(b, a) 

3388 >>> plt.figure() 

3389 >>> plt.semilogx(w[1:], -np.diff(np.unwrap(np.angle(h)))/np.diff(w)) 

3390 >>> plt.axhline(0.1, color='red') # 0.1 seconds group delay 

3391 >>> plt.title('Bessel filter group delay') 

3392 >>> plt.xlabel('Frequency [radians / second]') 

3393 >>> plt.ylabel('Group delay [seconds]') 

3394 >>> plt.margins(0, 0.1) 

3395 >>> plt.grid(which='both', axis='both') 

3396 >>> plt.show() 

3397 

3398 References 

3399 ---------- 

3400 .. [1] Thomson, W.E., "Delay Networks having Maximally Flat Frequency 

3401 Characteristics", Proceedings of the Institution of Electrical 

3402 Engineers, Part III, November 1949, Vol. 96, No. 44, pp. 487-490. 

3403 

3404 """ 

3405 return iirfilter(N, Wn, btype=btype, analog=analog, 

3406 output=output, ftype='bessel_'+norm, fs=fs) 

3407 

3408 

3409def maxflat(): 

3410 pass 

3411 

3412 

3413def yulewalk(): 

3414 pass 

3415 

3416 

3417def band_stop_obj(wp, ind, passb, stopb, gpass, gstop, type): 

3418 """ 

3419 Band Stop Objective Function for order minimization. 

3420 

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

3422 

3423 Parameters 

3424 ---------- 

3425 wp : scalar 

3426 Edge of passband `passb`. 

3427 ind : int, {0, 1} 

3428 Index specifying which `passb` edge to vary (0 or 1). 

3429 passb : ndarray 

3430 Two element sequence of fixed passband edges. 

3431 stopb : ndarray 

3432 Two element sequence of fixed stopband edges. 

3433 gstop : float 

3434 Amount of attenuation in stopband in dB. 

3435 gpass : float 

3436 Amount of ripple in the passband in dB. 

3437 type : {'butter', 'cheby', 'ellip'} 

3438 Type of filter. 

3439 

3440 Returns 

3441 ------- 

3442 n : scalar 

3443 Filter order (possibly non-integer). 

3444 

3445 """ 

3446 

3447 _validate_gpass_gstop(gpass, gstop) 

3448 

3449 passbC = passb.copy() 

3450 passbC[ind] = wp 

3451 nat = (stopb * (passbC[0] - passbC[1]) / 

3452 (stopb ** 2 - passbC[0] * passbC[1])) 

3453 nat = min(abs(nat)) 

3454 

3455 if type == 'butter': 

3456 GSTOP = 10 ** (0.1 * abs(gstop)) 

3457 GPASS = 10 ** (0.1 * abs(gpass)) 

3458 n = (log10((GSTOP - 1.0) / (GPASS - 1.0)) / (2 * log10(nat))) 

3459 elif type == 'cheby': 

3460 GSTOP = 10 ** (0.1 * abs(gstop)) 

3461 GPASS = 10 ** (0.1 * abs(gpass)) 

3462 n = arccosh(sqrt((GSTOP - 1.0) / (GPASS - 1.0))) / arccosh(nat) 

3463 elif type == 'ellip': 

3464 GSTOP = 10 ** (0.1 * gstop) 

3465 GPASS = 10 ** (0.1 * gpass) 

3466 arg1 = sqrt((GPASS - 1.0) / (GSTOP - 1.0)) 

3467 arg0 = 1.0 / nat 

3468 d0 = special.ellipk([arg0 ** 2, 1 - arg0 ** 2]) 

3469 d1 = special.ellipk([arg1 ** 2, 1 - arg1 ** 2]) 

3470 n = (d0[0] * d1[1] / (d0[1] * d1[0])) 

3471 else: 

3472 raise ValueError("Incorrect type: %s" % type) 

3473 return n 

3474 

3475 

3476def buttord(wp, ws, gpass, gstop, analog=False, fs=None): 

3477 """Butterworth filter order selection. 

3478 

3479 Return the order of the lowest order digital or analog Butterworth filter 

3480 that loses no more than `gpass` dB in the passband and has at least 

3481 `gstop` dB attenuation in the stopband. 

3482 

3483 Parameters 

3484 ---------- 

3485 wp, ws : float 

3486 Passband and stopband edge frequencies. 

3487 

3488 For digital filters, these are in the same units as `fs`. By default, 

3489 `fs` is 2 half-cycles/sample, so these are normalized from 0 to 1, 

3490 where 1 is the Nyquist frequency. (`wp` and `ws` are thus in 

3491 half-cycles / sample.) For example: 

3492 

3493 - Lowpass: wp = 0.2, ws = 0.3 

3494 - Highpass: wp = 0.3, ws = 0.2 

3495 - Bandpass: wp = [0.2, 0.5], ws = [0.1, 0.6] 

3496 - Bandstop: wp = [0.1, 0.6], ws = [0.2, 0.5] 

3497 

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

3499 gpass : float 

3500 The maximum loss in the passband (dB). 

3501 gstop : float 

3502 The minimum attenuation in the stopband (dB). 

3503 analog : bool, optional 

3504 When True, return an analog filter, otherwise a digital filter is 

3505 returned. 

3506 fs : float, optional 

3507 The sampling frequency of the digital system. 

3508 

3509 .. versionadded:: 1.2.0 

3510 

3511 Returns 

3512 ------- 

3513 ord : int 

3514 The lowest order for a Butterworth filter which meets specs. 

3515 wn : ndarray or float 

3516 The Butterworth natural frequency (i.e. the "3dB frequency"). Should 

3517 be used with `butter` to give filter results. If `fs` is specified, 

3518 this is in the same units, and `fs` must also be passed to `butter`. 

3519 

3520 See Also 

3521 -------- 

3522 butter : Filter design using order and critical points 

3523 cheb1ord : Find order and critical points from passband and stopband spec 

3524 cheb2ord, ellipord 

3525 iirfilter : General filter design using order and critical frequencies 

3526 iirdesign : General filter design using passband and stopband spec 

3527 

3528 Examples 

3529 -------- 

3530 Design an analog bandpass filter with passband within 3 dB from 20 to 

3531 50 rad/s, while rejecting at least -40 dB below 14 and above 60 rad/s. 

3532 Plot its frequency response, showing the passband and stopband 

3533 constraints in gray. 

3534 

3535 >>> from scipy import signal 

3536 >>> import matplotlib.pyplot as plt 

3537 

3538 >>> N, Wn = signal.buttord([20, 50], [14, 60], 3, 40, True) 

3539 >>> b, a = signal.butter(N, Wn, 'band', True) 

3540 >>> w, h = signal.freqs(b, a, np.logspace(1, 2, 500)) 

3541 >>> plt.semilogx(w, 20 * np.log10(abs(h))) 

3542 >>> plt.title('Butterworth bandpass filter fit to constraints') 

3543 >>> plt.xlabel('Frequency [radians / second]') 

3544 >>> plt.ylabel('Amplitude [dB]') 

3545 >>> plt.grid(which='both', axis='both') 

3546 >>> plt.fill([1, 14, 14, 1], [-40, -40, 99, 99], '0.9', lw=0) # stop 

3547 >>> plt.fill([20, 20, 50, 50], [-99, -3, -3, -99], '0.9', lw=0) # pass 

3548 >>> plt.fill([60, 60, 1e9, 1e9], [99, -40, -40, 99], '0.9', lw=0) # stop 

3549 >>> plt.axis([10, 100, -60, 3]) 

3550 >>> plt.show() 

3551 

3552 """ 

3553 

3554 _validate_gpass_gstop(gpass, gstop) 

3555 

3556 wp = atleast_1d(wp) 

3557 ws = atleast_1d(ws) 

3558 if fs is not None: 

3559 if analog: 

3560 raise ValueError("fs cannot be specified for an analog filter") 

3561 wp = 2*wp/fs 

3562 ws = 2*ws/fs 

3563 

3564 filter_type = 2 * (len(wp) - 1) 

3565 filter_type += 1 

3566 if wp[0] >= ws[0]: 

3567 filter_type += 1 

3568 

3569 # Pre-warp frequencies for digital filter design 

3570 if not analog: 

3571 passb = tan(pi * wp / 2.0) 

3572 stopb = tan(pi * ws / 2.0) 

3573 else: 

3574 passb = wp * 1.0 

3575 stopb = ws * 1.0 

3576 

3577 if filter_type == 1: # low 

3578 nat = stopb / passb 

3579 elif filter_type == 2: # high 

3580 nat = passb / stopb 

3581 elif filter_type == 3: # stop 

3582 wp0 = optimize.fminbound(band_stop_obj, passb[0], stopb[0] - 1e-12, 

3583 args=(0, passb, stopb, gpass, gstop, 

3584 'butter'), 

3585 disp=0) 

3586 passb[0] = wp0 

3587 wp1 = optimize.fminbound(band_stop_obj, stopb[1] + 1e-12, passb[1], 

3588 args=(1, passb, stopb, gpass, gstop, 

3589 'butter'), 

3590 disp=0) 

3591 passb[1] = wp1 

3592 nat = ((stopb * (passb[0] - passb[1])) / 

3593 (stopb ** 2 - passb[0] * passb[1])) 

3594 elif filter_type == 4: # pass 

3595 nat = ((stopb ** 2 - passb[0] * passb[1]) / 

3596 (stopb * (passb[0] - passb[1]))) 

3597 

3598 nat = min(abs(nat)) 

3599 

3600 GSTOP = 10 ** (0.1 * abs(gstop)) 

3601 GPASS = 10 ** (0.1 * abs(gpass)) 

3602 ord = int(ceil(log10((GSTOP - 1.0) / (GPASS - 1.0)) / (2 * log10(nat)))) 

3603 

3604 # Find the Butterworth natural frequency WN (or the "3dB" frequency") 

3605 # to give exactly gpass at passb. 

3606 try: 

3607 W0 = (GPASS - 1.0) ** (-1.0 / (2.0 * ord)) 

3608 except ZeroDivisionError: 

3609 W0 = 1.0 

3610 print("Warning, order is zero...check input parameters.") 

3611 

3612 # now convert this frequency back from lowpass prototype 

3613 # to the original analog filter 

3614 

3615 if filter_type == 1: # low 

3616 WN = W0 * passb 

3617 elif filter_type == 2: # high 

3618 WN = passb / W0 

3619 elif filter_type == 3: # stop 

3620 WN = numpy.zeros(2, float) 

3621 discr = sqrt((passb[1] - passb[0]) ** 2 + 

3622 4 * W0 ** 2 * passb[0] * passb[1]) 

3623 WN[0] = ((passb[1] - passb[0]) + discr) / (2 * W0) 

3624 WN[1] = ((passb[1] - passb[0]) - discr) / (2 * W0) 

3625 WN = numpy.sort(abs(WN)) 

3626 elif filter_type == 4: # pass 

3627 W0 = numpy.array([-W0, W0], float) 

3628 WN = (-W0 * (passb[1] - passb[0]) / 2.0 + 

3629 sqrt(W0 ** 2 / 4.0 * (passb[1] - passb[0]) ** 2 + 

3630 passb[0] * passb[1])) 

3631 WN = numpy.sort(abs(WN)) 

3632 else: 

3633 raise ValueError("Bad type: %s" % filter_type) 

3634 

3635 if not analog: 

3636 wn = (2.0 / pi) * arctan(WN) 

3637 else: 

3638 wn = WN 

3639 

3640 if len(wn) == 1: 

3641 wn = wn[0] 

3642 

3643 if fs is not None: 

3644 wn = wn*fs/2 

3645 

3646 return ord, wn 

3647 

3648 

3649def cheb1ord(wp, ws, gpass, gstop, analog=False, fs=None): 

3650 """Chebyshev type I filter order selection. 

3651 

3652 Return the order of the lowest order digital or analog Chebyshev Type I 

3653 filter that loses no more than `gpass` dB in the passband and has at 

3654 least `gstop` dB attenuation in the stopband. 

3655 

3656 Parameters 

3657 ---------- 

3658 wp, ws : float 

3659 Passband and stopband edge frequencies. 

3660 

3661 For digital filters, these are in the same units as `fs`. By default, 

3662 `fs` is 2 half-cycles/sample, so these are normalized from 0 to 1, 

3663 where 1 is the Nyquist frequency. (`wp` and `ws` are thus in 

3664 half-cycles / sample.) For example: 

3665 

3666 - Lowpass: wp = 0.2, ws = 0.3 

3667 - Highpass: wp = 0.3, ws = 0.2 

3668 - Bandpass: wp = [0.2, 0.5], ws = [0.1, 0.6] 

3669 - Bandstop: wp = [0.1, 0.6], ws = [0.2, 0.5] 

3670 

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

3672 gpass : float 

3673 The maximum loss in the passband (dB). 

3674 gstop : float 

3675 The minimum attenuation in the stopband (dB). 

3676 analog : bool, optional 

3677 When True, return an analog filter, otherwise a digital filter is 

3678 returned. 

3679 fs : float, optional 

3680 The sampling frequency of the digital system. 

3681 

3682 .. versionadded:: 1.2.0 

3683 

3684 Returns 

3685 ------- 

3686 ord : int 

3687 The lowest order for a Chebyshev type I filter that meets specs. 

3688 wn : ndarray or float 

3689 The Chebyshev natural frequency (the "3dB frequency") for use with 

3690 `cheby1` to give filter results. If `fs` is specified, 

3691 this is in the same units, and `fs` must also be passed to `cheby1`. 

3692 

3693 See Also 

3694 -------- 

3695 cheby1 : Filter design using order and critical points 

3696 buttord : Find order and critical points from passband and stopband spec 

3697 cheb2ord, ellipord 

3698 iirfilter : General filter design using order and critical frequencies 

3699 iirdesign : General filter design using passband and stopband spec 

3700 

3701 Examples 

3702 -------- 

3703 Design a digital lowpass filter such that the passband is within 3 dB up 

3704 to 0.2*(fs/2), while rejecting at least -40 dB above 0.3*(fs/2). Plot its 

3705 frequency response, showing the passband and stopband constraints in gray. 

3706 

3707 >>> from scipy import signal 

3708 >>> import matplotlib.pyplot as plt 

3709 

3710 >>> N, Wn = signal.cheb1ord(0.2, 0.3, 3, 40) 

3711 >>> b, a = signal.cheby1(N, 3, Wn, 'low') 

3712 >>> w, h = signal.freqz(b, a) 

3713 >>> plt.semilogx(w / np.pi, 20 * np.log10(abs(h))) 

3714 >>> plt.title('Chebyshev I lowpass filter fit to constraints') 

3715 >>> plt.xlabel('Normalized frequency') 

3716 >>> plt.ylabel('Amplitude [dB]') 

3717 >>> plt.grid(which='both', axis='both') 

3718 >>> plt.fill([.01, 0.2, 0.2, .01], [-3, -3, -99, -99], '0.9', lw=0) # stop 

3719 >>> plt.fill([0.3, 0.3, 2, 2], [ 9, -40, -40, 9], '0.9', lw=0) # pass 

3720 >>> plt.axis([0.08, 1, -60, 3]) 

3721 >>> plt.show() 

3722 

3723 """ 

3724 

3725 _validate_gpass_gstop(gpass, gstop) 

3726 

3727 wp = atleast_1d(wp) 

3728 ws = atleast_1d(ws) 

3729 if fs is not None: 

3730 if analog: 

3731 raise ValueError("fs cannot be specified for an analog filter") 

3732 wp = 2*wp/fs 

3733 ws = 2*ws/fs 

3734 

3735 filter_type = 2 * (len(wp) - 1) 

3736 if wp[0] < ws[0]: 

3737 filter_type += 1 

3738 else: 

3739 filter_type += 2 

3740 

3741 # Pre-warp frequencies for digital filter design 

3742 if not analog: 

3743 passb = tan(pi * wp / 2.0) 

3744 stopb = tan(pi * ws / 2.0) 

3745 else: 

3746 passb = wp * 1.0 

3747 stopb = ws * 1.0 

3748 

3749 if filter_type == 1: # low 

3750 nat = stopb / passb 

3751 elif filter_type == 2: # high 

3752 nat = passb / stopb 

3753 elif filter_type == 3: # stop 

3754 wp0 = optimize.fminbound(band_stop_obj, passb[0], stopb[0] - 1e-12, 

3755 args=(0, passb, stopb, gpass, gstop, 'cheby'), 

3756 disp=0) 

3757 passb[0] = wp0 

3758 wp1 = optimize.fminbound(band_stop_obj, stopb[1] + 1e-12, passb[1], 

3759 args=(1, passb, stopb, gpass, gstop, 'cheby'), 

3760 disp=0) 

3761 passb[1] = wp1 

3762 nat = ((stopb * (passb[0] - passb[1])) / 

3763 (stopb ** 2 - passb[0] * passb[1])) 

3764 elif filter_type == 4: # pass 

3765 nat = ((stopb ** 2 - passb[0] * passb[1]) / 

3766 (stopb * (passb[0] - passb[1]))) 

3767 

3768 nat = min(abs(nat)) 

3769 

3770 GSTOP = 10 ** (0.1 * abs(gstop)) 

3771 GPASS = 10 ** (0.1 * abs(gpass)) 

3772 ord = int(ceil(arccosh(sqrt((GSTOP - 1.0) / (GPASS - 1.0))) / 

3773 arccosh(nat))) 

3774 

3775 # Natural frequencies are just the passband edges 

3776 if not analog: 

3777 wn = (2.0 / pi) * arctan(passb) 

3778 else: 

3779 wn = passb 

3780 

3781 if len(wn) == 1: 

3782 wn = wn[0] 

3783 

3784 if fs is not None: 

3785 wn = wn*fs/2 

3786 

3787 return ord, wn 

3788 

3789 

3790def cheb2ord(wp, ws, gpass, gstop, analog=False, fs=None): 

3791 """Chebyshev type II filter order selection. 

3792 

3793 Return the order of the lowest order digital or analog Chebyshev Type II 

3794 filter that loses no more than `gpass` dB in the passband and has at least 

3795 `gstop` dB attenuation in the stopband. 

3796 

3797 Parameters 

3798 ---------- 

3799 wp, ws : float 

3800 Passband and stopband edge frequencies. 

3801 

3802 For digital filters, these are in the same units as `fs`. By default, 

3803 `fs` is 2 half-cycles/sample, so these are normalized from 0 to 1, 

3804 where 1 is the Nyquist frequency. (`wp` and `ws` are thus in 

3805 half-cycles / sample.) For example: 

3806 

3807 - Lowpass: wp = 0.2, ws = 0.3 

3808 - Highpass: wp = 0.3, ws = 0.2 

3809 - Bandpass: wp = [0.2, 0.5], ws = [0.1, 0.6] 

3810 - Bandstop: wp = [0.1, 0.6], ws = [0.2, 0.5] 

3811 

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

3813 gpass : float 

3814 The maximum loss in the passband (dB). 

3815 gstop : float 

3816 The minimum attenuation in the stopband (dB). 

3817 analog : bool, optional 

3818 When True, return an analog filter, otherwise a digital filter is 

3819 returned. 

3820 fs : float, optional 

3821 The sampling frequency of the digital system. 

3822 

3823 .. versionadded:: 1.2.0 

3824 

3825 Returns 

3826 ------- 

3827 ord : int 

3828 The lowest order for a Chebyshev type II filter that meets specs. 

3829 wn : ndarray or float 

3830 The Chebyshev natural frequency (the "3dB frequency") for use with 

3831 `cheby2` to give filter results. If `fs` is specified, 

3832 this is in the same units, and `fs` must also be passed to `cheby2`. 

3833 

3834 See Also 

3835 -------- 

3836 cheby2 : Filter design using order and critical points 

3837 buttord : Find order and critical points from passband and stopband spec 

3838 cheb1ord, ellipord 

3839 iirfilter : General filter design using order and critical frequencies 

3840 iirdesign : General filter design using passband and stopband spec 

3841 

3842 Examples 

3843 -------- 

3844 Design a digital bandstop filter which rejects -60 dB from 0.2*(fs/2) to 

3845 0.5*(fs/2), while staying within 3 dB below 0.1*(fs/2) or above 

3846 0.6*(fs/2). Plot its frequency response, showing the passband and 

3847 stopband constraints in gray. 

3848 

3849 >>> from scipy import signal 

3850 >>> import matplotlib.pyplot as plt 

3851 

3852 >>> N, Wn = signal.cheb2ord([0.1, 0.6], [0.2, 0.5], 3, 60) 

3853 >>> b, a = signal.cheby2(N, 60, Wn, 'stop') 

3854 >>> w, h = signal.freqz(b, a) 

3855 >>> plt.semilogx(w / np.pi, 20 * np.log10(abs(h))) 

3856 >>> plt.title('Chebyshev II bandstop filter fit to constraints') 

3857 >>> plt.xlabel('Normalized frequency') 

3858 >>> plt.ylabel('Amplitude [dB]') 

3859 >>> plt.grid(which='both', axis='both') 

3860 >>> plt.fill([.01, .1, .1, .01], [-3, -3, -99, -99], '0.9', lw=0) # stop 

3861 >>> plt.fill([.2, .2, .5, .5], [ 9, -60, -60, 9], '0.9', lw=0) # pass 

3862 >>> plt.fill([.6, .6, 2, 2], [-99, -3, -3, -99], '0.9', lw=0) # stop 

3863 >>> plt.axis([0.06, 1, -80, 3]) 

3864 >>> plt.show() 

3865 

3866 """ 

3867 

3868 _validate_gpass_gstop(gpass, gstop) 

3869 

3870 wp = atleast_1d(wp) 

3871 ws = atleast_1d(ws) 

3872 if fs is not None: 

3873 if analog: 

3874 raise ValueError("fs cannot be specified for an analog filter") 

3875 wp = 2*wp/fs 

3876 ws = 2*ws/fs 

3877 

3878 filter_type = 2 * (len(wp) - 1) 

3879 if wp[0] < ws[0]: 

3880 filter_type += 1 

3881 else: 

3882 filter_type += 2 

3883 

3884 # Pre-warp frequencies for digital filter design 

3885 if not analog: 

3886 passb = tan(pi * wp / 2.0) 

3887 stopb = tan(pi * ws / 2.0) 

3888 else: 

3889 passb = wp * 1.0 

3890 stopb = ws * 1.0 

3891 

3892 if filter_type == 1: # low 

3893 nat = stopb / passb 

3894 elif filter_type == 2: # high 

3895 nat = passb / stopb 

3896 elif filter_type == 3: # stop 

3897 wp0 = optimize.fminbound(band_stop_obj, passb[0], stopb[0] - 1e-12, 

3898 args=(0, passb, stopb, gpass, gstop, 'cheby'), 

3899 disp=0) 

3900 passb[0] = wp0 

3901 wp1 = optimize.fminbound(band_stop_obj, stopb[1] + 1e-12, passb[1], 

3902 args=(1, passb, stopb, gpass, gstop, 'cheby'), 

3903 disp=0) 

3904 passb[1] = wp1 

3905 nat = ((stopb * (passb[0] - passb[1])) / 

3906 (stopb ** 2 - passb[0] * passb[1])) 

3907 elif filter_type == 4: # pass 

3908 nat = ((stopb ** 2 - passb[0] * passb[1]) / 

3909 (stopb * (passb[0] - passb[1]))) 

3910 

3911 nat = min(abs(nat)) 

3912 

3913 GSTOP = 10 ** (0.1 * abs(gstop)) 

3914 GPASS = 10 ** (0.1 * abs(gpass)) 

3915 ord = int(ceil(arccosh(sqrt((GSTOP - 1.0) / (GPASS - 1.0))) / 

3916 arccosh(nat))) 

3917 

3918 # Find frequency where analog response is -gpass dB. 

3919 # Then convert back from low-pass prototype to the original filter. 

3920 

3921 new_freq = cosh(1.0 / ord * arccosh(sqrt((GSTOP - 1.0) / (GPASS - 1.0)))) 

3922 new_freq = 1.0 / new_freq 

3923 

3924 if filter_type == 1: 

3925 nat = passb / new_freq 

3926 elif filter_type == 2: 

3927 nat = passb * new_freq 

3928 elif filter_type == 3: 

3929 nat = numpy.zeros(2, float) 

3930 nat[0] = (new_freq / 2.0 * (passb[0] - passb[1]) + 

3931 sqrt(new_freq ** 2 * (passb[1] - passb[0]) ** 2 / 4.0 + 

3932 passb[1] * passb[0])) 

3933 nat[1] = passb[1] * passb[0] / nat[0] 

3934 elif filter_type == 4: 

3935 nat = numpy.zeros(2, float) 

3936 nat[0] = (1.0 / (2.0 * new_freq) * (passb[0] - passb[1]) + 

3937 sqrt((passb[1] - passb[0]) ** 2 / (4.0 * new_freq ** 2) + 

3938 passb[1] * passb[0])) 

3939 nat[1] = passb[0] * passb[1] / nat[0] 

3940 

3941 if not analog: 

3942 wn = (2.0 / pi) * arctan(nat) 

3943 else: 

3944 wn = nat 

3945 

3946 if len(wn) == 1: 

3947 wn = wn[0] 

3948 

3949 if fs is not None: 

3950 wn = wn*fs/2 

3951 

3952 return ord, wn 

3953 

3954 

3955def ellipord(wp, ws, gpass, gstop, analog=False, fs=None): 

3956 """Elliptic (Cauer) filter order selection. 

3957 

3958 Return the order of the lowest order digital or analog elliptic filter 

3959 that loses no more than `gpass` dB in the passband and has at least 

3960 `gstop` dB attenuation in the stopband. 

3961 

3962 Parameters 

3963 ---------- 

3964 wp, ws : float 

3965 Passband and stopband edge frequencies. 

3966 

3967 For digital filters, these are in the same units as `fs`. By default, 

3968 `fs` is 2 half-cycles/sample, so these are normalized from 0 to 1, 

3969 where 1 is the Nyquist frequency. (`wp` and `ws` are thus in 

3970 half-cycles / sample.) For example: 

3971 

3972 - Lowpass: wp = 0.2, ws = 0.3 

3973 - Highpass: wp = 0.3, ws = 0.2 

3974 - Bandpass: wp = [0.2, 0.5], ws = [0.1, 0.6] 

3975 - Bandstop: wp = [0.1, 0.6], ws = [0.2, 0.5] 

3976 

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

3978 gpass : float 

3979 The maximum loss in the passband (dB). 

3980 gstop : float 

3981 The minimum attenuation in the stopband (dB). 

3982 analog : bool, optional 

3983 When True, return an analog filter, otherwise a digital filter is 

3984 returned. 

3985 fs : float, optional 

3986 The sampling frequency of the digital system. 

3987 

3988 .. versionadded:: 1.2.0 

3989 

3990 Returns 

3991 ------- 

3992 ord : int 

3993 The lowest order for an Elliptic (Cauer) filter that meets specs. 

3994 wn : ndarray or float 

3995 The Chebyshev natural frequency (the "3dB frequency") for use with 

3996 `ellip` to give filter results. If `fs` is specified, 

3997 this is in the same units, and `fs` must also be passed to `ellip`. 

3998 

3999 See Also 

4000 -------- 

4001 ellip : Filter design using order and critical points 

4002 buttord : Find order and critical points from passband and stopband spec 

4003 cheb1ord, cheb2ord 

4004 iirfilter : General filter design using order and critical frequencies 

4005 iirdesign : General filter design using passband and stopband spec 

4006 

4007 Examples 

4008 -------- 

4009 Design an analog highpass filter such that the passband is within 3 dB 

4010 above 30 rad/s, while rejecting -60 dB at 10 rad/s. Plot its 

4011 frequency response, showing the passband and stopband constraints in gray. 

4012 

4013 >>> from scipy import signal 

4014 >>> import matplotlib.pyplot as plt 

4015 

4016 >>> N, Wn = signal.ellipord(30, 10, 3, 60, True) 

4017 >>> b, a = signal.ellip(N, 3, 60, Wn, 'high', True) 

4018 >>> w, h = signal.freqs(b, a, np.logspace(0, 3, 500)) 

4019 >>> plt.semilogx(w, 20 * np.log10(abs(h))) 

4020 >>> plt.title('Elliptical highpass filter fit to constraints') 

4021 >>> plt.xlabel('Frequency [radians / second]') 

4022 >>> plt.ylabel('Amplitude [dB]') 

4023 >>> plt.grid(which='both', axis='both') 

4024 >>> plt.fill([.1, 10, 10, .1], [1e4, 1e4, -60, -60], '0.9', lw=0) # stop 

4025 >>> plt.fill([30, 30, 1e9, 1e9], [-99, -3, -3, -99], '0.9', lw=0) # pass 

4026 >>> plt.axis([1, 300, -80, 3]) 

4027 >>> plt.show() 

4028 

4029 """ 

4030 

4031 _validate_gpass_gstop(gpass, gstop) 

4032 

4033 wp = atleast_1d(wp) 

4034 ws = atleast_1d(ws) 

4035 if fs is not None: 

4036 if analog: 

4037 raise ValueError("fs cannot be specified for an analog filter") 

4038 wp = 2*wp/fs 

4039 ws = 2*ws/fs 

4040 

4041 filter_type = 2 * (len(wp) - 1) 

4042 filter_type += 1 

4043 if wp[0] >= ws[0]: 

4044 filter_type += 1 

4045 

4046 # Pre-warp frequencies for digital filter design 

4047 if not analog: 

4048 passb = tan(pi * wp / 2.0) 

4049 stopb = tan(pi * ws / 2.0) 

4050 else: 

4051 passb = wp * 1.0 

4052 stopb = ws * 1.0 

4053 

4054 if filter_type == 1: # low 

4055 nat = stopb / passb 

4056 elif filter_type == 2: # high 

4057 nat = passb / stopb 

4058 elif filter_type == 3: # stop 

4059 wp0 = optimize.fminbound(band_stop_obj, passb[0], stopb[0] - 1e-12, 

4060 args=(0, passb, stopb, gpass, gstop, 'ellip'), 

4061 disp=0) 

4062 passb[0] = wp0 

4063 wp1 = optimize.fminbound(band_stop_obj, stopb[1] + 1e-12, passb[1], 

4064 args=(1, passb, stopb, gpass, gstop, 'ellip'), 

4065 disp=0) 

4066 passb[1] = wp1 

4067 nat = ((stopb * (passb[0] - passb[1])) / 

4068 (stopb ** 2 - passb[0] * passb[1])) 

4069 elif filter_type == 4: # pass 

4070 nat = ((stopb ** 2 - passb[0] * passb[1]) / 

4071 (stopb * (passb[0] - passb[1]))) 

4072 

4073 nat = min(abs(nat)) 

4074 

4075 GSTOP = 10 ** (0.1 * gstop) 

4076 GPASS = 10 ** (0.1 * gpass) 

4077 arg1 = sqrt((GPASS - 1.0) / (GSTOP - 1.0)) 

4078 arg0 = 1.0 / nat 

4079 d0 = special.ellipk([arg0 ** 2, 1 - arg0 ** 2]) 

4080 d1 = special.ellipk([arg1 ** 2, 1 - arg1 ** 2]) 

4081 ord = int(ceil(d0[0] * d1[1] / (d0[1] * d1[0]))) 

4082 

4083 if not analog: 

4084 wn = arctan(passb) * 2.0 / pi 

4085 else: 

4086 wn = passb 

4087 

4088 if len(wn) == 1: 

4089 wn = wn[0] 

4090 

4091 if fs is not None: 

4092 wn = wn*fs/2 

4093 

4094 return ord, wn 

4095 

4096 

4097def buttap(N): 

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

4099 

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

4101 

4102 See Also 

4103 -------- 

4104 butter : Filter design function using this prototype 

4105 

4106 """ 

4107 if abs(int(N)) != N: 

4108 raise ValueError("Filter order must be a nonnegative integer") 

4109 z = numpy.array([]) 

4110 m = numpy.arange(-N+1, N, 2) 

4111 # Middle value is 0 to ensure an exactly real pole 

4112 p = -numpy.exp(1j * pi * m / (2 * N)) 

4113 k = 1 

4114 return z, p, k 

4115 

4116 

4117def cheb1ap(N, rp): 

4118 """ 

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

4120 

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

4122 

4123 The filter's angular (e.g. rad/s) cutoff frequency is normalized to 1, 

4124 defined as the point at which the gain first drops below ``-rp``. 

4125 

4126 See Also 

4127 -------- 

4128 cheby1 : Filter design function using this prototype 

4129 

4130 """ 

4131 if abs(int(N)) != N: 

4132 raise ValueError("Filter order must be a nonnegative integer") 

4133 elif N == 0: 

4134 # Avoid divide-by-zero error 

4135 # Even order filters have DC gain of -rp dB 

4136 return numpy.array([]), numpy.array([]), 10**(-rp/20) 

4137 z = numpy.array([]) 

4138 

4139 # Ripple factor (epsilon) 

4140 eps = numpy.sqrt(10 ** (0.1 * rp) - 1.0) 

4141 mu = 1.0 / N * arcsinh(1 / eps) 

4142 

4143 # Arrange poles in an ellipse on the left half of the S-plane 

4144 m = numpy.arange(-N+1, N, 2) 

4145 theta = pi * m / (2*N) 

4146 p = -sinh(mu + 1j*theta) 

4147 

4148 k = numpy.prod(-p, axis=0).real 

4149 if N % 2 == 0: 

4150 k = k / sqrt((1 + eps * eps)) 

4151 

4152 return z, p, k 

4153 

4154 

4155def cheb2ap(N, rs): 

4156 """ 

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

4158 

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

4160 

4161 The filter's angular (e.g. rad/s) cutoff frequency is normalized to 1, 

4162 defined as the point at which the gain first reaches ``-rs``. 

4163 

4164 See Also 

4165 -------- 

4166 cheby2 : Filter design function using this prototype 

4167 

4168 """ 

4169 if abs(int(N)) != N: 

4170 raise ValueError("Filter order must be a nonnegative integer") 

4171 elif N == 0: 

4172 # Avoid divide-by-zero warning 

4173 return numpy.array([]), numpy.array([]), 1 

4174 

4175 # Ripple factor (epsilon) 

4176 de = 1.0 / sqrt(10 ** (0.1 * rs) - 1) 

4177 mu = arcsinh(1.0 / de) / N 

4178 

4179 if N % 2: 

4180 m = numpy.concatenate((numpy.arange(-N+1, 0, 2), 

4181 numpy.arange(2, N, 2))) 

4182 else: 

4183 m = numpy.arange(-N+1, N, 2) 

4184 

4185 z = -conjugate(1j / sin(m * pi / (2.0 * N))) 

4186 

4187 # Poles around the unit circle like Butterworth 

4188 p = -exp(1j * pi * numpy.arange(-N+1, N, 2) / (2 * N)) 

4189 # Warp into Chebyshev II 

4190 p = sinh(mu) * p.real + 1j * cosh(mu) * p.imag 

4191 p = 1.0 / p 

4192 

4193 k = (numpy.prod(-p, axis=0) / numpy.prod(-z, axis=0)).real 

4194 return z, p, k 

4195 

4196 

4197EPSILON = 2e-16 

4198 

4199 

4200def _vratio(u, ineps, mp): 

4201 [s, c, d, phi] = special.ellipj(u, mp) 

4202 ret = abs(ineps - s / c) 

4203 return ret 

4204 

4205 

4206def _kratio(m, k_ratio): 

4207 m = float(m) 

4208 if m < 0: 

4209 m = 0.0 

4210 if m > 1: 

4211 m = 1.0 

4212 if abs(m) > EPSILON and (abs(m) + EPSILON) < 1: 

4213 k = special.ellipk([m, 1 - m]) 

4214 r = k[0] / k[1] - k_ratio 

4215 elif abs(m) > EPSILON: 

4216 r = -k_ratio 

4217 else: 

4218 r = 1e20 

4219 return abs(r) 

4220 

4221 

4222def ellipap(N, rp, rs): 

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

4224 

4225 The filter is a normalized prototype that has `rp` decibels of ripple 

4226 in the passband and a stopband `rs` decibels down. 

4227 

4228 The filter's angular (e.g., rad/s) cutoff frequency is normalized to 1, 

4229 defined as the point at which the gain first drops below ``-rp``. 

4230 

4231 See Also 

4232 -------- 

4233 ellip : Filter design function using this prototype 

4234 

4235 References 

4236 ---------- 

4237 .. [1] Lutova, Tosic, and Evans, "Filter Design for Signal Processing", 

4238 Chapters 5 and 12. 

4239 

4240 """ 

4241 if abs(int(N)) != N: 

4242 raise ValueError("Filter order must be a nonnegative integer") 

4243 elif N == 0: 

4244 # Avoid divide-by-zero warning 

4245 # Even order filters have DC gain of -rp dB 

4246 return numpy.array([]), numpy.array([]), 10**(-rp/20) 

4247 elif N == 1: 

4248 p = -sqrt(1.0 / (10 ** (0.1 * rp) - 1.0)) 

4249 k = -p 

4250 z = [] 

4251 return asarray(z), asarray(p), k 

4252 

4253 eps = numpy.sqrt(10 ** (0.1 * rp) - 1) 

4254 ck1 = eps / numpy.sqrt(10 ** (0.1 * rs) - 1) 

4255 ck1p = numpy.sqrt(1 - ck1 * ck1) 

4256 if ck1p == 1: 

4257 raise ValueError("Cannot design a filter with given rp and rs" 

4258 " specifications.") 

4259 

4260 val = special.ellipk([ck1 * ck1, ck1p * ck1p]) 

4261 if abs(1 - ck1p * ck1p) < EPSILON: 

4262 krat = 0 

4263 else: 

4264 krat = N * val[0] / val[1] 

4265 

4266 m = optimize.fmin(_kratio, [0.5], args=(krat,), maxfun=250, maxiter=250, 

4267 disp=0) 

4268 if m < 0 or m > 1: 

4269 m = optimize.fminbound(_kratio, 0, 1, args=(krat,), maxfun=250, 

4270 disp=0) 

4271 

4272 capk = special.ellipk(m) 

4273 

4274 j = numpy.arange(1 - N % 2, N, 2) 

4275 jj = len(j) 

4276 

4277 [s, c, d, phi] = special.ellipj(j * capk / N, m * numpy.ones(jj)) 

4278 snew = numpy.compress(abs(s) > EPSILON, s, axis=-1) 

4279 z = 1.0 / (sqrt(m) * snew) 

4280 z = 1j * z 

4281 z = numpy.concatenate((z, conjugate(z))) 

4282 

4283 r = optimize.fmin(_vratio, special.ellipk(m), args=(1. / eps, ck1p * ck1p), 

4284 maxfun=250, maxiter=250, disp=0) 

4285 v0 = capk * r / (N * val[0]) 

4286 

4287 [sv, cv, dv, phi] = special.ellipj(v0, 1 - m) 

4288 p = -(c * d * sv * cv + 1j * s * dv) / (1 - (d * sv) ** 2.0) 

4289 

4290 if N % 2: 

4291 newp = numpy.compress(abs(p.imag) > EPSILON * 

4292 numpy.sqrt(numpy.sum(p * numpy.conjugate(p), 

4293 axis=0).real), 

4294 p, axis=-1) 

4295 p = numpy.concatenate((p, conjugate(newp))) 

4296 else: 

4297 p = numpy.concatenate((p, conjugate(p))) 

4298 

4299 k = (numpy.prod(-p, axis=0) / numpy.prod(-z, axis=0)).real 

4300 if N % 2 == 0: 

4301 k = k / numpy.sqrt((1 + eps * eps)) 

4302 

4303 return z, p, k 

4304 

4305 

4306# TODO: Make this a real public function scipy.misc.ff 

4307def _falling_factorial(x, n): 

4308 r""" 

4309 Return the factorial of `x` to the `n` falling. 

4310 

4311 This is defined as: 

4312 

4313 .. math:: x^\underline n = (x)_n = x (x-1) \cdots (x-n+1) 

4314 

4315 This can more efficiently calculate ratios of factorials, since: 

4316 

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

4318 

4319 where n >= m 

4320 

4321 skipping the factors that cancel out 

4322 

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

4324 """ 

4325 val = 1 

4326 for k in range(x - n + 1, x + 1): 

4327 val *= k 

4328 return val 

4329 

4330 

4331def _bessel_poly(n, reverse=False): 

4332 """ 

4333 Return the coefficients of Bessel polynomial of degree `n` 

4334 

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

4336 

4337 Output is a list of coefficients: 

4338 [1] = 1 

4339 [1, 1] = 1*s + 1 

4340 [1, 3, 3] = 1*s^2 + 3*s + 3 

4341 [1, 6, 15, 15] = 1*s^3 + 6*s^2 + 15*s + 15 

4342 [1, 10, 45, 105, 105] = 1*s^4 + 10*s^3 + 45*s^2 + 105*s + 105 

4343 etc. 

4344 

4345 Output is a Python list of arbitrary precision long ints, so n is only 

4346 limited by your hardware's memory. 

4347 

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

4349 match http://oeis.org/A001498/b001498.txt : 

4350 

4351 >>> i = 0 

4352 >>> for n in range(51): 

4353 ... for x in _bessel_poly(n, reverse=True): 

4354 ... print(i, x) 

4355 ... i += 1 

4356 

4357 """ 

4358 if abs(int(n)) != n: 

4359 raise ValueError("Polynomial order must be a nonnegative integer") 

4360 else: 

4361 n = int(n) # np.int32 doesn't work, for instance 

4362 

4363 out = [] 

4364 for k in range(n + 1): 

4365 num = _falling_factorial(2*n - k, n) 

4366 den = 2**(n - k) * factorial(k, exact=True) 

4367 out.append(num // den) 

4368 

4369 if reverse: 

4370 return out[::-1] 

4371 else: 

4372 return out 

4373 

4374 

4375def _campos_zeros(n): 

4376 """ 

4377 Return approximate zero locations of Bessel polynomials y_n(x) for order 

4378 `n` using polynomial fit (Campos-Calderon 2011) 

4379 """ 

4380 if n == 1: 

4381 return asarray([-1+0j]) 

4382 

4383 s = npp_polyval(n, [0, 0, 2, 0, -3, 1]) 

4384 b3 = npp_polyval(n, [16, -8]) / s 

4385 b2 = npp_polyval(n, [-24, -12, 12]) / s 

4386 b1 = npp_polyval(n, [8, 24, -12, -2]) / s 

4387 b0 = npp_polyval(n, [0, -6, 0, 5, -1]) / s 

4388 

4389 r = npp_polyval(n, [0, 0, 2, 1]) 

4390 a1 = npp_polyval(n, [-6, -6]) / r 

4391 a2 = 6 / r 

4392 

4393 k = np.arange(1, n+1) 

4394 x = npp_polyval(k, [0, a1, a2]) 

4395 y = npp_polyval(k, [b0, b1, b2, b3]) 

4396 

4397 return x + 1j*y 

4398 

4399 

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

4401 """ 

4402 Given a function `f`, its first derivative `fp`, and a set of initial 

4403 guesses `x0`, simultaneously find the roots of the polynomial using the 

4404 Aberth-Ehrlich method. 

4405 

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

4407 

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

4409 """ 

4410 

4411 N = len(x0) 

4412 

4413 x = array(x0, complex) 

4414 beta = np.empty_like(x0) 

4415 

4416 for iteration in range(maxiter): 

4417 alpha = -f(x) / fp(x) # Newton's method 

4418 

4419 # Model "repulsion" between zeros 

4420 for k in range(N): 

4421 beta[k] = np.sum(1/(x[k] - x[k+1:])) 

4422 beta[k] += np.sum(1/(x[k] - x[:k])) 

4423 

4424 x += alpha / (1 + alpha * beta) 

4425 

4426 if not all(np.isfinite(x)): 

4427 raise RuntimeError('Root-finding calculation failed') 

4428 

4429 # Mekwi: The iterative process can be stopped when |hn| has become 

4430 # less than the largest error one is willing to permit in the root. 

4431 if all(abs(alpha) <= tol): 

4432 break 

4433 else: 

4434 raise Exception('Zeros failed to converge') 

4435 

4436 return x 

4437 

4438 

4439def _bessel_zeros(N): 

4440 """ 

4441 Find zeros of ordinary Bessel polynomial of order `N`, by root-finding of 

4442 modified Bessel function of the second kind 

4443 """ 

4444 if N == 0: 

4445 return asarray([]) 

4446 

4447 # Generate starting points 

4448 x0 = _campos_zeros(N) 

4449 

4450 # Zeros are the same for exp(1/x)*K_{N+0.5}(1/x) and Nth-order ordinary 

4451 # Bessel polynomial y_N(x) 

4452 def f(x): 

4453 return special.kve(N+0.5, 1/x) 

4454 

4455 # First derivative of above 

4456 def fp(x): 

4457 return (special.kve(N-0.5, 1/x)/(2*x**2) - 

4458 special.kve(N+0.5, 1/x)/(x**2) + 

4459 special.kve(N+1.5, 1/x)/(2*x**2)) 

4460 

4461 # Starting points converge to true zeros 

4462 x = _aberth(f, fp, x0) 

4463 

4464 # Improve precision using Newton's method on each 

4465 for i in range(len(x)): 

4466 x[i] = optimize.newton(f, x[i], fp, tol=1e-15) 

4467 

4468 # Average complex conjugates to make them exactly symmetrical 

4469 x = np.mean((x, x[::-1].conj()), 0) 

4470 

4471 # Zeros should sum to -1 

4472 if abs(np.sum(x) + 1) > 1e-15: 

4473 raise RuntimeError('Generated zeros are inaccurate') 

4474 

4475 return x 

4476 

4477 

4478def _norm_factor(p, k): 

4479 """ 

4480 Numerically find frequency shift to apply to delay-normalized filter such 

4481 that -3 dB point is at 1 rad/sec. 

4482 

4483 `p` is an array_like of polynomial poles 

4484 `k` is a float gain 

4485 

4486 First 10 values are listed in "Bessel Scale Factors" table, 

4487 "Bessel Filters Polynomials, Poles and Circuit Elements 2003, C. Bond." 

4488 """ 

4489 p = asarray(p, dtype=complex) 

4490 

4491 def G(w): 

4492 """ 

4493 Gain of filter 

4494 """ 

4495 return abs(k / prod(1j*w - p)) 

4496 

4497 def cutoff(w): 

4498 """ 

4499 When gain = -3 dB, return 0 

4500 """ 

4501 return G(w) - 1/np.sqrt(2) 

4502 

4503 return optimize.newton(cutoff, 1.5) 

4504 

4505 

4506def besselap(N, norm='phase'): 

4507 """ 

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

4509 

4510 Parameters 

4511 ---------- 

4512 N : int 

4513 The order of the filter. 

4514 norm : {'phase', 'delay', 'mag'}, optional 

4515 Frequency normalization: 

4516 

4517 ``phase`` 

4518 The filter is normalized such that the phase response reaches its 

4519 midpoint at an angular (e.g., rad/s) cutoff frequency of 1. This 

4520 happens for both low-pass and high-pass filters, so this is the 

4521 "phase-matched" case. [6]_ 

4522 

4523 The magnitude response asymptotes are the same as a Butterworth 

4524 filter of the same order with a cutoff of `Wn`. 

4525 

4526 This is the default, and matches MATLAB's implementation. 

4527 

4528 ``delay`` 

4529 The filter is normalized such that the group delay in the passband 

4530 is 1 (e.g., 1 second). This is the "natural" type obtained by 

4531 solving Bessel polynomials 

4532 

4533 ``mag`` 

4534 The filter is normalized such that the gain magnitude is -3 dB at 

4535 angular frequency 1. This is called "frequency normalization" by 

4536 Bond. [1]_ 

4537 

4538 .. versionadded:: 0.18.0 

4539 

4540 Returns 

4541 ------- 

4542 z : ndarray 

4543 Zeros of the transfer function. Is always an empty array. 

4544 p : ndarray 

4545 Poles of the transfer function. 

4546 k : scalar 

4547 Gain of the transfer function. For phase-normalized, this is always 1. 

4548 

4549 See Also 

4550 -------- 

4551 bessel : Filter design function using this prototype 

4552 

4553 Notes 

4554 ----- 

4555 To find the pole locations, approximate starting points are generated [2]_ 

4556 for the zeros of the ordinary Bessel polynomial [3]_, then the 

4557 Aberth-Ehrlich method [4]_ [5]_ is used on the Kv(x) Bessel function to 

4558 calculate more accurate zeros, and these locations are then inverted about 

4559 the unit circle. 

4560 

4561 References 

4562 ---------- 

4563 .. [1] C.R. Bond, "Bessel Filter Constants", 

4564 http://www.crbond.com/papers/bsf.pdf 

4565 .. [2] Campos and Calderon, "Approximate closed-form formulas for the 

4566 zeros of the Bessel Polynomials", :arXiv:`1105.0957`. 

4567 .. [3] Thomson, W.E., "Delay Networks having Maximally Flat Frequency 

4568 Characteristics", Proceedings of the Institution of Electrical 

4569 Engineers, Part III, November 1949, Vol. 96, No. 44, pp. 487-490. 

4570 .. [4] Aberth, "Iteration Methods for Finding all Zeros of a Polynomial 

4571 Simultaneously", Mathematics of Computation, Vol. 27, No. 122, 

4572 April 1973 

4573 .. [5] Ehrlich, "A modified Newton method for polynomials", Communications 

4574 of the ACM, Vol. 10, Issue 2, pp. 107-108, Feb. 1967, 

4575 :DOI:`10.1145/363067.363115` 

4576 .. [6] Miller and Bohn, "A Bessel Filter Crossover, and Its Relation to 

4577 Others", RaneNote 147, 1998, https://www.ranecommercial.com/legacy/note147.html 

4578 

4579 """ 

4580 if abs(int(N)) != N: 

4581 raise ValueError("Filter order must be a nonnegative integer") 

4582 if N == 0: 

4583 p = [] 

4584 k = 1 

4585 else: 

4586 # Find roots of reverse Bessel polynomial 

4587 p = 1/_bessel_zeros(N) 

4588 

4589 a_last = _falling_factorial(2*N, N) // 2**N 

4590 

4591 # Shift them to a different normalization if required 

4592 if norm in ('delay', 'mag'): 

4593 # Normalized for group delay of 1 

4594 k = a_last 

4595 if norm == 'mag': 

4596 # -3 dB magnitude point is at 1 rad/sec 

4597 norm_factor = _norm_factor(p, k) 

4598 p /= norm_factor 

4599 k = norm_factor**-N * a_last 

4600 elif norm == 'phase': 

4601 # Phase-matched (1/2 max phase shift at 1 rad/sec) 

4602 # Asymptotes are same as Butterworth filter 

4603 p *= 10**(-math.log10(a_last)/N) 

4604 k = 1 

4605 else: 

4606 raise ValueError('normalization not understood') 

4607 

4608 return asarray([]), asarray(p, dtype=complex), float(k) 

4609 

4610 

4611def iirnotch(w0, Q, fs=2.0): 

4612 """ 

4613 Design second-order IIR notch digital filter. 

4614 

4615 A notch filter is a band-stop filter with a narrow bandwidth 

4616 (high quality factor). It rejects a narrow frequency band and 

4617 leaves the rest of the spectrum little changed. 

4618 

4619 Parameters 

4620 ---------- 

4621 w0 : float 

4622 Frequency to remove from a signal. If `fs` is specified, this is in 

4623 the same units as `fs`. By default, it is a normalized scalar that must 

4624 satisfy ``0 < w0 < 1``, with ``w0 = 1`` corresponding to half of the 

4625 sampling frequency. 

4626 Q : float 

4627 Quality factor. Dimensionless parameter that characterizes 

4628 notch filter -3 dB bandwidth ``bw`` relative to its center 

4629 frequency, ``Q = w0/bw``. 

4630 fs : float, optional 

4631 The sampling frequency of the digital system. 

4632 

4633 .. versionadded:: 1.2.0 

4634 

4635 Returns 

4636 ------- 

4637 b, a : ndarray, ndarray 

4638 Numerator (``b``) and denominator (``a``) polynomials 

4639 of the IIR filter. 

4640 

4641 See Also 

4642 -------- 

4643 iirpeak 

4644 

4645 Notes 

4646 ----- 

4647 .. versionadded:: 0.19.0 

4648 

4649 References 

4650 ---------- 

4651 .. [1] Sophocles J. Orfanidis, "Introduction To Signal Processing", 

4652 Prentice-Hall, 1996 

4653 

4654 Examples 

4655 -------- 

4656 Design and plot filter to remove the 60 Hz component from a 

4657 signal sampled at 200 Hz, using a quality factor Q = 30 

4658 

4659 >>> from scipy import signal 

4660 >>> import matplotlib.pyplot as plt 

4661 

4662 >>> fs = 200.0 # Sample frequency (Hz) 

4663 >>> f0 = 60.0 # Frequency to be removed from signal (Hz) 

4664 >>> Q = 30.0 # Quality factor 

4665 >>> # Design notch filter 

4666 >>> b, a = signal.iirnotch(f0, Q, fs) 

4667 

4668 >>> # Frequency response 

4669 >>> freq, h = signal.freqz(b, a, fs=fs) 

4670 >>> # Plot 

4671 >>> fig, ax = plt.subplots(2, 1, figsize=(8, 6)) 

4672 >>> ax[0].plot(freq, 20*np.log10(abs(h)), color='blue') 

4673 >>> ax[0].set_title("Frequency Response") 

4674 >>> ax[0].set_ylabel("Amplitude (dB)", color='blue') 

4675 >>> ax[0].set_xlim([0, 100]) 

4676 >>> ax[0].set_ylim([-25, 10]) 

4677 >>> ax[0].grid() 

4678 >>> ax[1].plot(freq, np.unwrap(np.angle(h))*180/np.pi, color='green') 

4679 >>> ax[1].set_ylabel("Angle (degrees)", color='green') 

4680 >>> ax[1].set_xlabel("Frequency (Hz)") 

4681 >>> ax[1].set_xlim([0, 100]) 

4682 >>> ax[1].set_yticks([-90, -60, -30, 0, 30, 60, 90]) 

4683 >>> ax[1].set_ylim([-90, 90]) 

4684 >>> ax[1].grid() 

4685 >>> plt.show() 

4686 """ 

4687 

4688 return _design_notch_peak_filter(w0, Q, "notch", fs) 

4689 

4690 

4691def iirpeak(w0, Q, fs=2.0): 

4692 """ 

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

4694 

4695 A peak filter is a band-pass filter with a narrow bandwidth 

4696 (high quality factor). It rejects components outside a narrow 

4697 frequency band. 

4698 

4699 Parameters 

4700 ---------- 

4701 w0 : float 

4702 Frequency to be retained in a signal. If `fs` is specified, this is in 

4703 the same units as `fs`. By default, it is a normalized scalar that must 

4704 satisfy ``0 < w0 < 1``, with ``w0 = 1`` corresponding to half of the 

4705 sampling frequency. 

4706 Q : float 

4707 Quality factor. Dimensionless parameter that characterizes 

4708 peak filter -3 dB bandwidth ``bw`` relative to its center 

4709 frequency, ``Q = w0/bw``. 

4710 fs : float, optional 

4711 The sampling frequency of the digital system. 

4712 

4713 .. versionadded:: 1.2.0 

4714 

4715 Returns 

4716 ------- 

4717 b, a : ndarray, ndarray 

4718 Numerator (``b``) and denominator (``a``) polynomials 

4719 of the IIR filter. 

4720 

4721 See Also 

4722 -------- 

4723 iirnotch 

4724 

4725 Notes 

4726 ----- 

4727 .. versionadded:: 0.19.0 

4728 

4729 References 

4730 ---------- 

4731 .. [1] Sophocles J. Orfanidis, "Introduction To Signal Processing", 

4732 Prentice-Hall, 1996 

4733 

4734 Examples 

4735 -------- 

4736 Design and plot filter to remove the frequencies other than the 300 Hz 

4737 component from a signal sampled at 1000 Hz, using a quality factor Q = 30 

4738 

4739 >>> from scipy import signal 

4740 >>> import matplotlib.pyplot as plt 

4741 

4742 >>> fs = 1000.0 # Sample frequency (Hz) 

4743 >>> f0 = 300.0 # Frequency to be retained (Hz) 

4744 >>> Q = 30.0 # Quality factor 

4745 >>> # Design peak filter 

4746 >>> b, a = signal.iirpeak(f0, Q, fs) 

4747 

4748 >>> # Frequency response 

4749 >>> freq, h = signal.freqz(b, a, fs=fs) 

4750 >>> # Plot 

4751 >>> fig, ax = plt.subplots(2, 1, figsize=(8, 6)) 

4752 >>> ax[0].plot(freq, 20*np.log10(np.maximum(abs(h), 1e-5)), color='blue') 

4753 >>> ax[0].set_title("Frequency Response") 

4754 >>> ax[0].set_ylabel("Amplitude (dB)", color='blue') 

4755 >>> ax[0].set_xlim([0, 500]) 

4756 >>> ax[0].set_ylim([-50, 10]) 

4757 >>> ax[0].grid() 

4758 >>> ax[1].plot(freq, np.unwrap(np.angle(h))*180/np.pi, color='green') 

4759 >>> ax[1].set_ylabel("Angle (degrees)", color='green') 

4760 >>> ax[1].set_xlabel("Frequency (Hz)") 

4761 >>> ax[1].set_xlim([0, 500]) 

4762 >>> ax[1].set_yticks([-90, -60, -30, 0, 30, 60, 90]) 

4763 >>> ax[1].set_ylim([-90, 90]) 

4764 >>> ax[1].grid() 

4765 >>> plt.show() 

4766 """ 

4767 

4768 return _design_notch_peak_filter(w0, Q, "peak", fs) 

4769 

4770 

4771def _design_notch_peak_filter(w0, Q, ftype, fs=2.0): 

4772 """ 

4773 Design notch or peak digital filter. 

4774 

4775 Parameters 

4776 ---------- 

4777 w0 : float 

4778 Normalized frequency to remove from a signal. If `fs` is specified, 

4779 this is in the same units as `fs`. By default, it is a normalized 

4780 scalar that must satisfy ``0 < w0 < 1``, with ``w0 = 1`` 

4781 corresponding to half of the sampling frequency. 

4782 Q : float 

4783 Quality factor. Dimensionless parameter that characterizes 

4784 notch filter -3 dB bandwidth ``bw`` relative to its center 

4785 frequency, ``Q = w0/bw``. 

4786 ftype : str 

4787 The type of IIR filter to design: 

4788 

4789 - notch filter : ``notch`` 

4790 - peak filter : ``peak`` 

4791 fs : float, optional 

4792 The sampling frequency of the digital system. 

4793 

4794 .. versionadded:: 1.2.0: 

4795 

4796 Returns 

4797 ------- 

4798 b, a : ndarray, ndarray 

4799 Numerator (``b``) and denominator (``a``) polynomials 

4800 of the IIR filter. 

4801 """ 

4802 

4803 # Guarantee that the inputs are floats 

4804 w0 = float(w0) 

4805 Q = float(Q) 

4806 w0 = 2*w0/fs 

4807 

4808 # Checks if w0 is within the range 

4809 if w0 > 1.0 or w0 < 0.0: 

4810 raise ValueError("w0 should be such that 0 < w0 < 1") 

4811 

4812 # Get bandwidth 

4813 bw = w0/Q 

4814 

4815 # Normalize inputs 

4816 bw = bw*np.pi 

4817 w0 = w0*np.pi 

4818 

4819 # Compute -3dB atenuation 

4820 gb = 1/np.sqrt(2) 

4821 

4822 if ftype == "notch": 

4823 # Compute beta: formula 11.3.4 (p.575) from reference [1] 

4824 beta = (np.sqrt(1.0-gb**2.0)/gb)*np.tan(bw/2.0) 

4825 elif ftype == "peak": 

4826 # Compute beta: formula 11.3.19 (p.579) from reference [1] 

4827 beta = (gb/np.sqrt(1.0-gb**2.0))*np.tan(bw/2.0) 

4828 else: 

4829 raise ValueError("Unknown ftype.") 

4830 

4831 # Compute gain: formula 11.3.6 (p.575) from reference [1] 

4832 gain = 1.0/(1.0+beta) 

4833 

4834 # Compute numerator b and denominator a 

4835 # formulas 11.3.7 (p.575) and 11.3.21 (p.579) 

4836 # from reference [1] 

4837 if ftype == "notch": 

4838 b = gain*np.array([1.0, -2.0*np.cos(w0), 1.0]) 

4839 else: 

4840 b = (1.0-gain)*np.array([1.0, 0.0, -1.0]) 

4841 a = np.array([1.0, -2.0*gain*np.cos(w0), (2.0*gain-1.0)]) 

4842 

4843 return b, a 

4844 

4845 

4846filter_dict = {'butter': [buttap, buttord], 

4847 'butterworth': [buttap, buttord], 

4848 

4849 'cauer': [ellipap, ellipord], 

4850 'elliptic': [ellipap, ellipord], 

4851 'ellip': [ellipap, ellipord], 

4852 

4853 'bessel': [besselap], 

4854 'bessel_phase': [besselap], 

4855 'bessel_delay': [besselap], 

4856 'bessel_mag': [besselap], 

4857 

4858 'cheby1': [cheb1ap, cheb1ord], 

4859 'chebyshev1': [cheb1ap, cheb1ord], 

4860 'chebyshevi': [cheb1ap, cheb1ord], 

4861 

4862 'cheby2': [cheb2ap, cheb2ord], 

4863 'chebyshev2': [cheb2ap, cheb2ord], 

4864 'chebyshevii': [cheb2ap, cheb2ord], 

4865 } 

4866 

4867band_dict = {'band': 'bandpass', 

4868 'bandpass': 'bandpass', 

4869 'pass': 'bandpass', 

4870 'bp': 'bandpass', 

4871 

4872 'bs': 'bandstop', 

4873 'bandstop': 'bandstop', 

4874 'bands': 'bandstop', 

4875 'stop': 'bandstop', 

4876 

4877 'l': 'lowpass', 

4878 'low': 'lowpass', 

4879 'lowpass': 'lowpass', 

4880 'lp': 'lowpass', 

4881 

4882 'high': 'highpass', 

4883 'highpass': 'highpass', 

4884 'h': 'highpass', 

4885 'hp': 'highpass', 

4886 } 

4887 

4888bessel_norms = {'bessel': 'phase', 

4889 'bessel_phase': 'phase', 

4890 'bessel_delay': 'delay', 

4891 'bessel_mag': 'mag'}