1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237
|
########################################################################
##
## Copyright (C) 1994-2024 The Octave Project Developers
##
## See the file COPYRIGHT.md in the top-level directory of this
## distribution or <https://octave.org/copyright/>.
##
## This file is part of Octave.
##
## Octave is free software: you can redistribute it and/or modify it
## under the terms of the GNU General Public License as published by
## the Free Software Foundation, either version 3 of the License, or
## (at your option) any later version.
##
## Octave is distributed in the hope that it will be useful, but
## WITHOUT ANY WARRANTY; without even the implied warranty of
## MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
## GNU General Public License for more details.
##
## You should have received a copy of the GNU General Public License
## along with Octave; see the file COPYING. If not, see
## <https://www.gnu.org/licenses/>.
##
########################################################################
## -*- texinfo -*-
## @deftypefn {} {[@var{h}, @var{w}] =} freqz (@var{b}, @var{a}, @var{n}, "whole")
## @deftypefnx {} {[@var{h}, @var{w}] =} freqz (@var{b})
## @deftypefnx {} {[@var{h}, @var{w}] =} freqz (@var{b}, @var{a})
## @deftypefnx {} {[@var{h}, @var{w}] =} freqz (@var{b}, @var{a}, @var{n})
## @deftypefnx {} {@var{h} =} freqz (@var{b}, @var{a}, @var{w})
## @deftypefnx {} {[@var{h}, @var{w}] =} freqz (@dots{}, @var{Fs})
## @deftypefnx {} {} freqz (@dots{})
##
## Return the complex frequency response @var{h} of the rational IIR filter
## whose numerator and denominator coefficients are @var{b} and @var{a},
## respectively.
##
## The response is evaluated at @var{n} angular frequencies between 0 and
## @ifnottex
## 2*pi.
## @end ifnottex
## @tex
## $2\pi$.
## @end tex
##
## @noindent
## The output value @var{w} is a vector of the frequencies.
##
## If @var{a} is omitted, the denominator is assumed to be 1 (this
## corresponds to a simple FIR filter).
##
## If @var{n} is omitted, a value of 512 is assumed. For fastest computation,
## @var{n} should factor into a small number of small primes.
##
## If the fourth argument, @qcode{"whole"}, is omitted the response is
## evaluated at frequencies between 0 and
## @ifnottex
## pi.
## @end ifnottex
## @tex
## $\pi$.
## @end tex
##
## @code{freqz (@var{b}, @var{a}, @var{w})}
##
## Evaluate the response at the specific frequencies in the vector @var{w}.
## The values for @var{w} are measured in radians.
##
## @code{[@dots{}] = freqz (@dots{}, @var{Fs})}
##
## Return frequencies in Hz instead of radians assuming a sampling rate
## @var{Fs}. If you are evaluating the response at specific frequencies
## @var{w}, those frequencies should be requested in Hz rather than radians.
##
## @code{freqz (@dots{})}
##
## Plot the magnitude and phase response of @var{h} rather than returning them.
##
## @seealso{freqz_plot}
## @end deftypefn
function [h_r, f_r] = freqz (b, a, n, region, Fs)
if (nargin < 1)
print_usage ();
elseif (nargin == 1)
## Response of an FIR filter.
a = n = region = Fs = [];
elseif (nargin == 2)
## Response of an IIR filter
n = region = Fs = [];
elseif (nargin == 3)
region = Fs = [];
elseif (nargin == 4)
Fs = [];
if (! ischar (region) && ! isempty (region))
Fs = region;
region = [];
endif
endif
if (isempty (b))
b = 1;
elseif (! isvector (b))
error ("freqz: B must be a vector");
endif
if (isempty (a))
a = 1;
elseif (! isvector (a))
error ("freqz: A must be a vector");
endif
if (isempty (n))
n = 512;
elseif (isscalar (n) && n < 1)
error ("freqz: N must be a positive integer");
endif
if (isempty (region))
if (isreal (b) && isreal (a))
region = "half";
else
region = "whole";
endif
endif
if (isempty (Fs))
freq_norm = true;
if (nargout == 0)
Fs = 2;
else
Fs = 2*pi;
endif
else
freq_norm = false;
endif
plot_output = (nargout == 0);
whole_region = strcmp (region, "whole");
a = a(:);
b = b(:);
if (! isscalar (n))
## Explicit frequency vector given
w = f = n;
if (nargin == 4)
## Sampling rate Fs was specified
w = 2*pi*f/Fs;
endif
k = max (length (b), length (a));
hb = polyval (postpad (b, k), exp (j*w));
ha = polyval (postpad (a, k), exp (j*w));
else
## polyval(fliplr(P),exp(jw)) is O(p n) and fft(x) is O(n log(n)),
## where p is the order of the polynomial P. For small p it
## would be faster to use polyval but in practice the overhead for
## polyval is much higher and the little bit of time saved isn't
## worth the extra code.
k = max (length (b), length (a));
if (k > n/2 && nargout == 0)
## Ensure a causal phase response.
n *= 2 .^ ceil (log2 (2*k/n));
endif
if (whole_region)
N = n;
if (plot_output)
f = Fs * (0:n).' / N; # do 1 more for the plot
else
f = Fs * (0:n-1).' / N;
endif
else
N = 2*n;
if (plot_output)
n += 1;
endif
f = Fs * (0:n-1).' / N;
endif
pad_sz = N*ceil (k/N);
b = postpad (b, pad_sz);
a = postpad (a, pad_sz);
hb = zeros (n, 1);
ha = zeros (n, 1);
for i = 1:N:pad_sz
hb += fft (postpad (b(i:i+N-1), N))(1:n);
ha += fft (postpad (a(i:i+N-1), N))(1:n);
endfor
endif
h = hb ./ ha;
if (plot_output)
## Plot and don't return values.
if (whole_region && isscalar (n))
h(end+1) = h(1); # Solution is periodic. Copy first value to end.
endif
freqz_plot (f, h, freq_norm);
else
## Return values and don't plot.
h_r = h;
f_r = f;
endif
endfunction
%!testif HAVE_FFTW # correct values and fft-polyval consistency
%! ## butterworth filter, order 2, cutoff pi/2 radians
%! b = [0.292893218813452 0.585786437626905 0.292893218813452];
%! a = [1 0 0.171572875253810];
%! [h,w] = freqz (b,a,32);
%! assert (h(1),1,10*eps);
%! assert (abs (h(17)).^2,0.5,10*eps);
%! assert (h,freqz (b,a,w),10*eps); # fft should be consistent with polyval
%!testif HAVE_FFTW # whole-half consistency
%! b = [1 1 1]/3; # 3-sample average
%! [h,w] = freqz (b,1,32,"whole");
%! assert (h(2:16),conj (h(32:-1:18)),20*eps);
%! [h2,w2] = freqz (b,1,16,"half");
%! assert (h(1:16),h2,20*eps);
%! assert (w(1:16),w2,20*eps);
%!testif HAVE_FFTW # Sampling frequency properly interpreted
%! b = [1 1 1]/3; a = [1 0.2];
%! [h,f] = freqz (b,a,16,320);
%! assert (f,[0:15]'*10,10*eps);
%! [h2,f2] = freqz (b,a,[0:15]*10,320);
%! assert (f2,[0:15]*10,10*eps);
%! assert (h,h2.',20*eps);
%! [h3,f3] = freqz (b,a,32,"whole",320);
%! assert (f3,[0:31]'*10,10*eps);
## Test input validation
## FIXME: Need to put tests here and simplify input validation in the main code.
|