Skip to content

Commit

Permalink
ENH: Fix remez to recognize 'hilbert' option. Thanks to kumanna.
Browse files Browse the repository at this point in the history
  • Loading branch information
charris committed Aug 18, 2010
1 parent 1c3f75f commit 23135f4
Show file tree
Hide file tree
Showing 2 changed files with 39 additions and 5 deletions.
7 changes: 5 additions & 2 deletions scipy/signal/signaltools.py
Expand Up @@ -514,6 +514,9 @@ def remez(numtaps, bands, desired, weight=None, Hz=1, type='bandpass',
type --- The type of filter:
'bandpass' : flat response in bands.
'differentiator' : frequency proportional response in bands.
'hilbert' : filter with odd symmetry, that is, type III
(for even order) or type IV (for odd order)
linear phase filters
Outputs: (out,)
Expand All @@ -523,9 +526,9 @@ def remez(numtaps, bands, desired, weight=None, Hz=1, type='bandpass',
"""
# Convert type
try:
tnum = {'bandpass':1, 'differentiator':2}[type]
tnum = {'bandpass':1, 'differentiator':2, 'hilbert':3}[type]
except KeyError:
raise ValueError, "Type must be 'bandpass', or 'differentiator'"
raise ValueError("Type must be 'bandpass', 'differentiator', or 'hilbert'")

# Convert weight
if weight is None:
Expand Down
37 changes: 34 additions & 3 deletions scipy/signal/tests/test_filter_design.py
@@ -1,9 +1,9 @@
import warnings

import numpy as np
from numpy.testing import TestCase, assert_array_almost_equal
from numpy.testing import *

from scipy.signal import tf2zpk, bessel, BadCoefficients, kaiserord, firwin, freqz
from scipy.signal import tf2zpk, bessel, BadCoefficients, kaiserord, firwin, freqz, remez


class TestTf2zpk(TestCase):
Expand Down Expand Up @@ -40,7 +40,7 @@ def test_bad_filter(self):


class TestFirWin(TestCase):

def test_lowpass(self):
width = 0.04
ntaps, beta = kaiserord(120, width)
Expand All @@ -49,3 +49,34 @@ def test_lowpass(self):
freqs, response = freqz(taps, worN=np.pi*freq_samples)
assert_array_almost_equal(np.abs(response),
[1.0, 1.0, 1.0, 0.0, 0.0, 0.0], decimal=5)

class TestRemez(TestCase):

def test_hilbert(self):
N = 11 # number of taps in the filter
a = 0.1 # width of the transition band

# design an unity gain hilbert bandpass filter from w to 0.5-w
h = remez(11, [ a, 0.5-a ], [ 1 ], type='hilbert')

# make sure the filter has correct # of taps
assert_(len(h) == N, "Number of Taps")

# make sure it is type III (anti-symmtric tap coefficients)
assert_array_almost_equal(h[:(N-1)/2], -h[:-(N-1)/2-1:-1])

# Since the requested response is symmetric, all even coeffcients
# should be zero (or in this case really small)
assert_((abs(h[1::2]) < 1e-15).all(), "Even Coefficients Equal Zero")

# now check the frequency response
w, H = freqz(h, 1)
f = w/2/np.pi
Hmag = abs(H)

# should have a zero at 0 and pi (in this case close to zero)
assert_((Hmag[ [0,-1] ] < 0.02).all(), "Zero at zero and pi")

# check that the pass band is close to unity
idx = (f > a) * (f < 0.5-a)
assert_((abs(Hmag[idx] - 1) < 0.015).all(), "Pass Band Close To Unity")

0 comments on commit 23135f4

Please sign in to comment.