2847
|
1 ## Copyright (C) 1996, 1997 John W. Eaton |
2313
|
2 ## |
|
3 ## This file is part of Octave. |
|
4 ## |
|
5 ## Octave is free software; you can redistribute it and/or modify it |
|
6 ## under the terms of the GNU General Public License as published by |
|
7 ## the Free Software Foundation; either version 2, or (at your option) |
|
8 ## any later version. |
|
9 ## |
|
10 ## Octave is distributed in the hope that it will be useful, but |
|
11 ## WITHOUT ANY WARRANTY; without even the implied warranty of |
|
12 ## MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU |
|
13 ## General Public License for more details. |
|
14 ## |
|
15 ## You should have received a copy of the GNU General Public License |
|
16 ## along with Octave; see the file COPYING. If not, write to the Free |
|
17 ## Software Foundation, 59 Temple Place - Suite 330, Boston, MA |
|
18 ## 02111-1307, USA. |
2303
|
19 |
3367
|
20 ## -*- texinfo -*- |
|
21 ## @deftypefn {Function File} {[@var{h}, @var{w}] =} freqz (@var{b}, @var{a}, @var{n}, "whole") |
|
22 ## Return the complex frequency response @var{h} of the rational IIR filter |
|
23 ## whose numerator and denominator coefficients are @var{b} and @var{a}, |
|
24 ## respectively. The response is evaluated at @var{n} angular frequencies |
|
25 ## between 0 and |
|
26 ## @ifinfo |
|
27 ## 2*pi. |
|
28 ## @end ifinfo |
|
29 ## @iftex |
|
30 ## @tex |
|
31 ## $2\pi$. |
|
32 ## @end tex |
|
33 ## @end iftex |
3426
|
34 ## |
3367
|
35 ## @noindent |
|
36 ## The output value @var{w} is a vector of the frequencies. |
3426
|
37 ## |
3367
|
38 ## If the fourth argument is omitted, the response is evaluated at |
|
39 ## frequencies between 0 and |
|
40 ## @ifinfo |
|
41 ## pi. |
|
42 ## @end ifinfo |
|
43 ## @iftex |
|
44 ## @tex |
|
45 ## $\pi$. |
|
46 ## @end tex |
|
47 ## @end iftex |
3426
|
48 ## |
3367
|
49 ## If @var{n} is omitted, a value of 512 is assumed. |
3426
|
50 ## |
3367
|
51 ## If @var{a} is omitted, the denominator is assumed to be 1 (this |
|
52 ## corresponds to a simple FIR filter). |
3426
|
53 ## |
3367
|
54 ## For fastest computation, @var{n} should factor into a small number of |
|
55 ## small primes. |
3893
|
56 ## |
|
57 ## @deftypefnx {Function File} {@var{h} =} freqz (@var{b}, @var{a}, @var{w}) |
|
58 ## Evaluate the response at the specific frequencies in the vector @var{w}. |
|
59 ## The values for @var{w} are measured in radians. |
|
60 ## |
|
61 ## @deftypefnx {Function File} {[@dots{}] =} freqz (@dots{}, @var{Fs}) |
|
62 ## Return frequencies in Hz instead of radians assuming a sampling rate |
|
63 ## @var{Fs}. If you are evaluating the response at specific frequencies |
|
64 ## @var{w}, those frequencies should be requested in Hz rather than radians. |
|
65 ## |
3920
|
66 ## @deftypefnx {Function File} {} freqz (@dots{}) |
3907
|
67 ## Plot the pass band, stop band and phase response of @var{h} rather |
|
68 ## than returning them. |
3367
|
69 ## @end deftypefn |
904
|
70 |
3367
|
71 ## Author: jwe ??? |
2314
|
72 |
3907
|
73 function [h_r, w_r] = freqz (b, a, n, region, Fs) |
566
|
74 |
3893
|
75 if (nargin < 1 || nargin > 5) |
|
76 usage ("[h, w] = freqz (b, a, n [, \"whole\"] [, Fs])"); |
|
77 elseif (nargin == 1) |
2303
|
78 ## Response of an FIR filter. |
3893
|
79 a = n = region = Fs = []; |
566
|
80 elseif (nargin == 2) |
2303
|
81 ## Response of an IIR filter |
3893
|
82 n = region = Fs = []; |
566
|
83 elseif (nargin == 3) |
3893
|
84 region = Fs = []; |
566
|
85 elseif (nargin == 4) |
3893
|
86 Fs = []; |
|
87 if (! isstr (region) && ! isempty (region)) |
|
88 Fs = region; |
|
89 region = []; |
|
90 endif |
|
91 endif |
|
92 |
|
93 if (isempty (a)) |
|
94 a = 1; |
|
95 endif |
|
96 if (isempty (n)) |
|
97 n = 512; |
|
98 endif |
|
99 if (isempty (region)) |
|
100 if (isreal (b) && isreal (a)) |
|
101 region = "half"; |
|
102 else |
|
103 region = "whole"; |
|
104 endif |
|
105 endif |
|
106 if (isempty (Fs)) |
|
107 if (nargout == 0) |
|
108 Fs = 2; |
|
109 else |
|
110 Fs = 2*pi; |
|
111 endif |
566
|
112 endif |
|
113 |
3893
|
114 la = length (a); |
|
115 a = reshape (a, 1, la); |
|
116 lb = length (b); |
|
117 b = reshape (b, 1, lb); |
|
118 k = max ([la, lb]); |
566
|
119 |
3893
|
120 if (! is_scalar (n)) |
|
121 if (nargin == 4) ## Fs was specified |
|
122 w = 2*pi*n/Fs; |
566
|
123 else |
3893
|
124 w = n; |
566
|
125 endif |
3893
|
126 n = length (n); |
|
127 extent = 0; |
|
128 elseif (strcmp (region, "whole")) |
|
129 w = 2 * pi * (0:n-1) / n; |
|
130 extent = n; |
566
|
131 else |
3893
|
132 w = pi * (0:n-1) / n; |
|
133 extent = 2 * n; |
566
|
134 endif |
|
135 |
3893
|
136 if (length (b) == 1) |
|
137 if (length (a) == 1) |
|
138 hb = b * ones (1, n); |
|
139 else |
|
140 hb = b; |
|
141 endif |
|
142 elseif (extent >= k) |
|
143 hb = fft (postpad (b, extent)); |
|
144 hb = hb(1:n); |
|
145 else |
|
146 hb = polyval (postpad (b, k), exp (j*w)); |
|
147 endif |
|
148 if (length (a) == 1) |
|
149 ha = a; |
|
150 elseif (extent >= k) |
|
151 ha = fft (postpad (a, extent)); |
|
152 ha = ha(1:n); |
|
153 else |
|
154 ha = polyval (postpad (a, k), exp (j*w)); |
|
155 endif |
|
156 h = hb ./ ha; |
|
157 w = Fs * w / (2*pi); |
|
158 |
3907
|
159 if (nargout != 0), # return values and don't plot |
|
160 h_r = h; |
|
161 w_r = w; |
|
162 else # plot and don't return values |
|
163 freqz_plot (w, h); |
|
164 end |
|
165 |
566
|
166 endfunction |