1 ## Copyright (C) 2001 Paulo Neis <p_neis@yahoo.com.br>
2 ## Copyright (C) 2003 Doug Stewart <dastew@sympatico.ca>
4 ## This program is free software; you can redistribute it and/or modify it under
5 ## the terms of the GNU General Public License as published by the Free Software
6 ## Foundation; either version 3 of the License, or (at your option) any later
9 ## This program is distributed in the hope that it will be useful, but WITHOUT
10 ## ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
11 ## FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more
14 ## You should have received a copy of the GNU General Public License along with
15 ## this program; if not, see <http://www.gnu.org/licenses/>.
18 ##usage: [Zz, Zp, Zg] = ellip(n, Rp, Rs, Wp, stype,'s')
20 ## Generate an Elliptic or Cauer filter (discrete and contnuious).
22 ## [b,a] = ellip(n, Rp, Rs, Wp)
23 ## low pass filter with order n, cutoff pi*Wp radians, Rp decibels
24 ## of ripple in the passband and a stopband Rs decibels down.
26 ## [b,a] = ellip(n, Rp, Rs, Wp, 'high')
27 ## high pass filter with cutoff pi*Wp...
29 ## [b,a] = ellip(n, Rp, Rs, [Wl, Wh])
30 ## band pass filter with band pass edges pi*Wl and pi*Wh ...
32 ## [b,a] = ellip(n, Rp, Rs, [Wl, Wh], 'stop')
33 ## band reject filter with edges pi*Wl and pi*Wh, ...
35 ## [z,p,g] = ellip(...)
36 ## return filter as zero-pole-gain.
38 ## [...] = ellip(...,'s')
39 ## return a Laplace space filter, W can be larger than 1.
41 ## [a,b,c,d] = ellip(...)
42 ## return state-space matrices
46 ## - Oppenheim, Alan V., Discrete Time Signal Processing, Hardcover, 1999.
47 ## - Parente Ribeiro, E., Notas de aula da disciplina TE498 - Processamento
48 ## Digital de Sinais, UFPR, 2001/2002.
49 ## - Kienzle, Paul, functions from Octave-Forge, 1999 (http://octave.sf.net).
52 function [a,b,c,d] = ellip(n, Rp, Rs, W, varargin)
54 if (nargin>6 || nargin<4) || (nargout>4 || nargout<2)
58 ## interpret the input parameters
59 if (!(length(n)==1 && n == round(n) && n > 0))
60 error ("ellip: filter order n must be a positive integer");
66 for i=1:length(varargin)
68 case 's', digital = 0;
69 case 'z', digital = 1;
70 case { 'high', 'stop' }, stop = 1;
71 case { 'low', 'pass' }, stop = 0;
72 otherwise, error ("ellip: expected [high|stop] or [s|z]");
77 if (!(length(W)<=2 && (r==1 || c==1)))
78 error ("ellip: frequency must be given as w0 or [w0, w1]");
79 elseif (!(length(W)==1 || length(W) == 2))
80 error ("ellip: only one filter band allowed");
81 elseif (length(W)==2 && !(W(1) < W(2)))
82 error ("ellip: first band edge must be smaller than second");
85 if ( digital && !all(W >= 0 & W <= 1))
86 error ("ellip: critical frequencies must be in (0 1)");
87 elseif ( !digital && !all(W >= 0 ))
88 error ("ellip: critical frequencies must be in (0 inf)");
92 error("ellip: passband ripple must be positive decibels");
96 error("ellip: stopband ripple must be positive decibels");
100 ##Prewarp the digital frequencies
102 T = 2; # sampling frequency of 2 Hz
106 ##Generate s-plane poles, zeros and gain
107 [zero, pole, gain] = ncauer(Rp, Rs, n);
109 ## splane frequency transform
110 [zero, pole, gain] = sftrans(zero, pole, gain, W, stop);
112 ## Use bilinear transform to convert poles to the z plane
114 [zero, pole, gain] = bilinear(zero, pole, gain, T);
118 ## convert to the correct output form
120 a = real(gain*poly(zero));
121 b = real(poly(pole));
128 [a, b, c, d] = zp2ss (zero, pole, gain);
135 %! disp('---------------------------> NELLIP 0.2 EXAMPLE <-------------------------')
136 %! x=input("Let's calculate the filter order: [ENTER]");
138 %! x=input("[n, Ws] = ellipord([.1 .2],.4,1,90); [ENTER]");
139 %! [n, Ws] = ellipord([.1 .2],.4,1,90)
141 %! x=input("Let's calculate the filter: [ENTER]");
143 %! x=input("[b,a]=ellip(5,1,90,[.1,.2]); [ENTER]");
144 %! [b,a]=ellip(5,1,90,[.1,.2])
146 %! x=input("Let's calculate the frequency response: [ENTER]");
148 %! x=input("[h,w]=freqz(b,a); [ENTER]");
151 %! xlabel("Frequency");
152 %! ylabel("abs(H[w])[dB]");
153 %! axis([0,1,-100,0]);
154 %! plot(w./pi, 20*log10(abs(h)), ';;')
157 %! x=ones(1,length(h));
158 %! plot(w./pi, x.*-1, ';-1 dB;')
159 %! plot(w./pi, x.*-90, ';-90 dB;')