annotate scripts/signal/freqz.m @ 14363:f3d52523cde1

Use Octave coding conventions in all m-file %!test blocks * wavread.m, acosd.m, acot.m, acotd.m, acoth.m, acsc.m, acscd.m, acsch.m, asec.m, asecd.m, asech.m, asind.m, atand.m, cosd.m, cot.m, cotd.m, coth.m, csc.m, cscd.m, csch.m, sec.m, secd.m, sech.m, sind.m, tand.m, accumarray.m, accumdim.m, bitcmp.m, bitget.m, bitset.m, blkdiag.m, cart2pol.m, cart2sph.m, celldisp.m, chop.m, circshift.m, colon.m, common_size.m, cplxpair.m, cumtrapz.m, curl.m, dblquad.m, deal.m, divergence.m, flipdim.m, fliplr.m, flipud.m, genvarname.m, gradient.m, idivide.m, int2str.m, interp1.m, interp1q.m, interp2.m, interp3.m, interpft.m, interpn.m, isa.m, isdir.m, isequal.m, isequalwithequalnans.m, issquare.m, logspace.m, nargchk.m, narginchk.m, nargoutchk.m, nextpow2.m, nthargout.m, num2str.m, pol2cart.m, polyarea.m, postpad.m, prepad.m, profile.m, profshow.m, quadgk.m, quadv.m, randi.m, rat.m, repmat.m, rot90.m, rotdim.m, shift.m, shiftdim.m, sph2cart.m, structfun.m, trapz.m, triplequad.m, convhull.m, dsearch.m, dsearchn.m, griddata3.m, griddatan.m, rectint.m, tsearchn.m, __makeinfo__.m, doc.m, get_first_help_sentence.m, help.m, type.m, unimplemented.m, which.m, imread.m, imwrite.m, dlmwrite.m, fileread.m, is_valid_file_id.m, strread.m, textread.m, textscan.m, commutation_matrix.m, cond.m, condest.m, cross.m, duplication_matrix.m, expm.m, housh.m, isdefinite.m, ishermitian.m, issymmetric.m, logm.m, normest.m, null.m, onenormest.m, orth.m, planerot.m, qzhess.m, rank.m, rref.m, trace.m, vech.m, ans.m, bincoeff.m, bug_report.m, bzip2.m, comma.m, compare_versions.m, computer.m, edit.m, fileparts.m, fullfile.m, getfield.m, gzip.m, info.m, inputname.m, isappdata.m, isdeployed.m, ismac.m, ispc.m, isunix.m, list_primes.m, ls.m, mexext.m, namelengthmax.m, news.m, orderfields.m, paren.m, recycle.m, rmappdata.m, semicolon.m, setappdata.m, setfield.m, substruct.m, symvar.m, ver.m, version.m, warning_ids.m, xor.m, fminbnd.m, fsolve.m, fzero.m, lsqnonneg.m, optimset.m, pqpnonneg.m, sqp.m, matlabroot.m, __gnuplot_drawnow__.m, __plt_get_axis_arg__.m, ancestor.m, cla.m, clf.m, close.m, colorbar.m, colstyle.m, comet3.m, contourc.m, figure.m, gca.m, gcbf.m, gcbo.m, gcf.m, ginput.m, graphics_toolkit.m, gtext.m, hggroup.m, hist.m, hold.m, isfigure.m, ishghandle.m, ishold.m, isocolors.m, isonormals.m, isosurface.m, isprop.m, legend.m, line.m, loglog.m, loglogerr.m, meshgrid.m, ndgrid.m, newplot.m, orient.m, patch.m, plot3.m, plotyy.m, __print_parse_opts__.m, quiver3.m, refreshdata.m, ribbon.m, semilogx.m, semilogxerr.m, semilogy.m, stem.m, stem3.m, subplot.m, title.m, uigetfile.m, view.m, whitebg.m, compan.m, conv.m, deconv.m, mkpp.m, mpoles.m, pchip.m, poly.m, polyaffine.m, polyder.m, polyfit.m, polygcd.m, polyint.m, polyout.m, polyval.m, polyvalm.m, ppder.m, ppint.m, ppjumps.m, ppval.m, residue.m, roots.m, spline.m, intersect.m, ismember.m, powerset.m, setdiff.m, setxor.m, union.m, unique.m, autoreg_matrix.m, bartlett.m, blackman.m, detrend.m, fftconv.m, fftfilt.m, fftshift.m, freqz.m, hamming.m, hanning.m, ifftshift.m, sinc.m, sinetone.m, sinewave.m, unwrap.m, bicg.m, bicgstab.m, gmres.m, gplot.m, nonzeros.m, pcg.m, pcr.m, spaugment.m, spconvert.m, spdiags.m, speye.m, spfun.m, spones.m, sprand.m, sprandsym.m, spstats.m, spy.m, svds.m, treelayout.m, bessel.m, beta.m, betaln.m, factor.m, factorial.m, isprime.m, lcm.m, legendre.m, nchoosek.m, nthroot.m, perms.m, pow2.m, primes.m, reallog.m, realpow.m, realsqrt.m, hadamard.m, hankel.m, hilb.m, invhilb.m, magic.m, rosser.m, vander.m, __finish__.m, center.m, cloglog.m, corr.m, cov.m, gls.m, histc.m, iqr.m, kendall.m, kurtosis.m, logit.m, mahalanobis.m, mean.m, meansq.m, median.m, mode.m, moment.m, ols.m, ppplot.m, prctile.m, probit.m, quantile.m, range.m, ranks.m, run_count.m, runlength.m, skewness.m, spearman.m, statistics.m, std.m, table.m, var.m, zscore.m, betacdf.m, betainv.m, betapdf.m, betarnd.m, binocdf.m, binoinv.m, binopdf.m, binornd.m, cauchy_cdf.m, cauchy_inv.m, cauchy_pdf.m, cauchy_rnd.m, chi2cdf.m, chi2inv.m, chi2pdf.m, chi2rnd.m, discrete_cdf.m, discrete_inv.m, discrete_pdf.m, discrete_rnd.m, empirical_cdf.m, empirical_inv.m, empirical_pdf.m, empirical_rnd.m, expcdf.m, expinv.m, exppdf.m, exprnd.m, fcdf.m, finv.m, fpdf.m, frnd.m, gamcdf.m, gaminv.m, gampdf.m, gamrnd.m, geocdf.m, geoinv.m, geopdf.m, geornd.m, hygecdf.m, hygeinv.m, hygepdf.m, hygernd.m, kolmogorov_smirnov_cdf.m, laplace_cdf.m, laplace_inv.m, laplace_pdf.m, laplace_rnd.m, logistic_cdf.m, logistic_inv.m, logistic_pdf.m, logistic_rnd.m, logncdf.m, logninv.m, lognpdf.m, lognrnd.m, nbincdf.m, nbininv.m, nbinpdf.m, nbinrnd.m, normcdf.m, norminv.m, normpdf.m, normrnd.m, poisscdf.m, poissinv.m, poisspdf.m, poissrnd.m, stdnormal_cdf.m, stdnormal_inv.m, stdnormal_pdf.m, stdnormal_rnd.m, tcdf.m, tinv.m, tpdf.m, trnd.m, unidcdf.m, unidinv.m, unidpdf.m, unidrnd.m, unifcdf.m, unifinv.m, unifpdf.m, unifrnd.m, wblcdf.m, wblinv.m, wblpdf.m, wblrnd.m, kolmogorov_smirnov_test.m, kruskal_wallis_test.m, base2dec.m, bin2dec.m, blanks.m, cstrcat.m, deblank.m, dec2base.m, dec2bin.m, dec2hex.m, findstr.m, hex2dec.m, index.m, isletter.m, mat2str.m, rindex.m, str2num.m, strcat.m, strjust.m, strmatch.m, strsplit.m, strtok.m, strtrim.m, strtrunc.m, substr.m, validatestring.m, demo.m, example.m, fail.m, speed.m, addtodate.m, asctime.m, clock.m, ctime.m, date.m, datenum.m, datetick.m, datevec.m, eomday.m, etime.m, is_leap_year.m, now.m: Use Octave coding conventions in all m-file %!test blocks
author Rik <octave@nomad.inbox5.com>
date Mon, 13 Feb 2012 07:29:44 -0800
parents 72c96de7a403
children 36dba9be680b
Ignore whitespace changes - Everywhere: Within whitespace: At end of lines:
rev   line source
14138
72c96de7a403 maint: update copyright notices for 2012
John W. Eaton <jwe@octave.org>
parents: 11587
diff changeset
1 ## Copyright (C) 1994-2012 John W. Eaton
2313
5ca126254d15 [project @ 1996-07-11 21:25:22 by jwe]
jwe
parents: 2311
diff changeset
2 ##
5ca126254d15 [project @ 1996-07-11 21:25:22 by jwe]
jwe
parents: 2311
diff changeset
3 ## This file is part of Octave.
5ca126254d15 [project @ 1996-07-11 21:25:22 by jwe]
jwe
parents: 2311
diff changeset
4 ##
5ca126254d15 [project @ 1996-07-11 21:25:22 by jwe]
jwe
parents: 2311
diff changeset
5 ## Octave is free software; you can redistribute it and/or modify it
5ca126254d15 [project @ 1996-07-11 21:25:22 by jwe]
jwe
parents: 2311
diff changeset
6 ## under the terms of the GNU General Public License as published by
7016
93c65f2a5668 [project @ 2007-10-12 06:40:56 by jwe]
jwe
parents: 6653
diff changeset
7 ## the Free Software Foundation; either version 3 of the License, or (at
93c65f2a5668 [project @ 2007-10-12 06:40:56 by jwe]
jwe
parents: 6653
diff changeset
8 ## your option) any later version.
2313
5ca126254d15 [project @ 1996-07-11 21:25:22 by jwe]
jwe
parents: 2311
diff changeset
9 ##
5ca126254d15 [project @ 1996-07-11 21:25:22 by jwe]
jwe
parents: 2311
diff changeset
10 ## Octave is distributed in the hope that it will be useful, but
5ca126254d15 [project @ 1996-07-11 21:25:22 by jwe]
jwe
parents: 2311
diff changeset
11 ## WITHOUT ANY WARRANTY; without even the implied warranty of
5ca126254d15 [project @ 1996-07-11 21:25:22 by jwe]
jwe
parents: 2311
diff changeset
12 ## MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
5ca126254d15 [project @ 1996-07-11 21:25:22 by jwe]
jwe
parents: 2311
diff changeset
13 ## General Public License for more details.
5ca126254d15 [project @ 1996-07-11 21:25:22 by jwe]
jwe
parents: 2311
diff changeset
14 ##
5ca126254d15 [project @ 1996-07-11 21:25:22 by jwe]
jwe
parents: 2311
diff changeset
15 ## You should have received a copy of the GNU General Public License
7016
93c65f2a5668 [project @ 2007-10-12 06:40:56 by jwe]
jwe
parents: 6653
diff changeset
16 ## along with Octave; see the file COPYING. If not, see
93c65f2a5668 [project @ 2007-10-12 06:40:56 by jwe]
jwe
parents: 6653
diff changeset
17 ## <http://www.gnu.org/licenses/>.
2303
5cffc4b8de57 [project @ 1996-06-24 09:15:24 by jwe]
jwe
parents: 1337
diff changeset
18
3367
0748b03c3510 [project @ 1999-11-20 14:52:38 by jwe]
jwe
parents: 2847
diff changeset
19 ## -*- texinfo -*-
0748b03c3510 [project @ 1999-11-20 14:52:38 by jwe]
jwe
parents: 2847
diff changeset
20 ## @deftypefn {Function File} {[@var{h}, @var{w}] =} freqz (@var{b}, @var{a}, @var{n}, "whole")
0748b03c3510 [project @ 1999-11-20 14:52:38 by jwe]
jwe
parents: 2847
diff changeset
21 ## Return the complex frequency response @var{h} of the rational IIR filter
0748b03c3510 [project @ 1999-11-20 14:52:38 by jwe]
jwe
parents: 2847
diff changeset
22 ## whose numerator and denominator coefficients are @var{b} and @var{a},
0748b03c3510 [project @ 1999-11-20 14:52:38 by jwe]
jwe
parents: 2847
diff changeset
23 ## respectively. The response is evaluated at @var{n} angular frequencies
0748b03c3510 [project @ 1999-11-20 14:52:38 by jwe]
jwe
parents: 2847
diff changeset
24 ## between 0 and
8517
81d6ab3ac93c Allow documentation tobe built for other formats than tex and info
sh@sh-laptop
parents: 8506
diff changeset
25 ## @ifnottex
3367
0748b03c3510 [project @ 1999-11-20 14:52:38 by jwe]
jwe
parents: 2847
diff changeset
26 ## 2*pi.
8517
81d6ab3ac93c Allow documentation tobe built for other formats than tex and info
sh@sh-laptop
parents: 8506
diff changeset
27 ## @end ifnottex
3367
0748b03c3510 [project @ 1999-11-20 14:52:38 by jwe]
jwe
parents: 2847
diff changeset
28 ## @tex
0748b03c3510 [project @ 1999-11-20 14:52:38 by jwe]
jwe
parents: 2847
diff changeset
29 ## $2\pi$.
0748b03c3510 [project @ 1999-11-20 14:52:38 by jwe]
jwe
parents: 2847
diff changeset
30 ## @end tex
3426
f8dde1807dee [project @ 2000-01-13 08:40:00 by jwe]
jwe
parents: 3367
diff changeset
31 ##
3367
0748b03c3510 [project @ 1999-11-20 14:52:38 by jwe]
jwe
parents: 2847
diff changeset
32 ## @noindent
0748b03c3510 [project @ 1999-11-20 14:52:38 by jwe]
jwe
parents: 2847
diff changeset
33 ## The output value @var{w} is a vector of the frequencies.
3426
f8dde1807dee [project @ 2000-01-13 08:40:00 by jwe]
jwe
parents: 3367
diff changeset
34 ##
3367
0748b03c3510 [project @ 1999-11-20 14:52:38 by jwe]
jwe
parents: 2847
diff changeset
35 ## If the fourth argument is omitted, the response is evaluated at
0748b03c3510 [project @ 1999-11-20 14:52:38 by jwe]
jwe
parents: 2847
diff changeset
36 ## frequencies between 0 and
8517
81d6ab3ac93c Allow documentation tobe built for other formats than tex and info
sh@sh-laptop
parents: 8506
diff changeset
37 ## @ifnottex
3367
0748b03c3510 [project @ 1999-11-20 14:52:38 by jwe]
jwe
parents: 2847
diff changeset
38 ## pi.
8517
81d6ab3ac93c Allow documentation tobe built for other formats than tex and info
sh@sh-laptop
parents: 8506
diff changeset
39 ## @end ifnottex
3367
0748b03c3510 [project @ 1999-11-20 14:52:38 by jwe]
jwe
parents: 2847
diff changeset
40 ## @tex
0748b03c3510 [project @ 1999-11-20 14:52:38 by jwe]
jwe
parents: 2847
diff changeset
41 ## $\pi$.
0748b03c3510 [project @ 1999-11-20 14:52:38 by jwe]
jwe
parents: 2847
diff changeset
42 ## @end tex
3426
f8dde1807dee [project @ 2000-01-13 08:40:00 by jwe]
jwe
parents: 3367
diff changeset
43 ##
3367
0748b03c3510 [project @ 1999-11-20 14:52:38 by jwe]
jwe
parents: 2847
diff changeset
44 ## If @var{n} is omitted, a value of 512 is assumed.
3426
f8dde1807dee [project @ 2000-01-13 08:40:00 by jwe]
jwe
parents: 3367
diff changeset
45 ##
3367
0748b03c3510 [project @ 1999-11-20 14:52:38 by jwe]
jwe
parents: 2847
diff changeset
46 ## If @var{a} is omitted, the denominator is assumed to be 1 (this
0748b03c3510 [project @ 1999-11-20 14:52:38 by jwe]
jwe
parents: 2847
diff changeset
47 ## corresponds to a simple FIR filter).
3426
f8dde1807dee [project @ 2000-01-13 08:40:00 by jwe]
jwe
parents: 3367
diff changeset
48 ##
3367
0748b03c3510 [project @ 1999-11-20 14:52:38 by jwe]
jwe
parents: 2847
diff changeset
49 ## For fastest computation, @var{n} should factor into a small number of
0748b03c3510 [project @ 1999-11-20 14:52:38 by jwe]
jwe
parents: 2847
diff changeset
50 ## small primes.
3893
abd8659eea11 [project @ 2002-04-09 21:36:31 by jwe]
jwe
parents: 3457
diff changeset
51 ##
abd8659eea11 [project @ 2002-04-09 21:36:31 by jwe]
jwe
parents: 3457
diff changeset
52 ## @deftypefnx {Function File} {@var{h} =} freqz (@var{b}, @var{a}, @var{w})
abd8659eea11 [project @ 2002-04-09 21:36:31 by jwe]
jwe
parents: 3457
diff changeset
53 ## Evaluate the response at the specific frequencies in the vector @var{w}.
abd8659eea11 [project @ 2002-04-09 21:36:31 by jwe]
jwe
parents: 3457
diff changeset
54 ## The values for @var{w} are measured in radians.
abd8659eea11 [project @ 2002-04-09 21:36:31 by jwe]
jwe
parents: 3457
diff changeset
55 ##
abd8659eea11 [project @ 2002-04-09 21:36:31 by jwe]
jwe
parents: 3457
diff changeset
56 ## @deftypefnx {Function File} {[@dots{}] =} freqz (@dots{}, @var{Fs})
abd8659eea11 [project @ 2002-04-09 21:36:31 by jwe]
jwe
parents: 3457
diff changeset
57 ## Return frequencies in Hz instead of radians assuming a sampling rate
11587
c792872f8942 all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents: 11523
diff changeset
58 ## @var{Fs}. If you are evaluating the response at specific frequencies
3893
abd8659eea11 [project @ 2002-04-09 21:36:31 by jwe]
jwe
parents: 3457
diff changeset
59 ## @var{w}, those frequencies should be requested in Hz rather than radians.
abd8659eea11 [project @ 2002-04-09 21:36:31 by jwe]
jwe
parents: 3457
diff changeset
60 ##
3920
87db95b22f8f [project @ 2002-05-01 04:07:31 by jwe]
jwe
parents: 3909
diff changeset
61 ## @deftypefnx {Function File} {} freqz (@dots{})
3907
437884fae441 [project @ 2002-04-24 19:33:08 by jwe]
jwe
parents: 3893
diff changeset
62 ## Plot the pass band, stop band and phase response of @var{h} rather
437884fae441 [project @ 2002-04-24 19:33:08 by jwe]
jwe
parents: 3893
diff changeset
63 ## than returning them.
3367
0748b03c3510 [project @ 1999-11-20 14:52:38 by jwe]
jwe
parents: 2847
diff changeset
64 ## @end deftypefn
904
3470f1e25a79 [project @ 1994-11-09 21:22:15 by jwe]
jwe
parents: 566
diff changeset
65
3367
0748b03c3510 [project @ 1999-11-20 14:52:38 by jwe]
jwe
parents: 2847
diff changeset
66 ## Author: jwe ???
2314
949ab8eba8bc [project @ 1996-07-12 03:58:02 by jwe]
jwe
parents: 2313
diff changeset
67
5377
bd4ee620c38a [project @ 2005-06-02 15:42:39 by jwe]
jwe
parents: 5307
diff changeset
68 function [h_r, f_r] = freqz (b, a, n, region, Fs)
566
8529a21443fa [project @ 1994-07-26 02:07:38 by jwe]
jwe
parents:
diff changeset
69
3893
abd8659eea11 [project @ 2002-04-09 21:36:31 by jwe]
jwe
parents: 3457
diff changeset
70 if (nargin < 1 || nargin > 5)
6046
34f96dd5441b [project @ 2006-10-10 16:10:25 by jwe]
jwe
parents: 5986
diff changeset
71 print_usage ();
3893
abd8659eea11 [project @ 2002-04-09 21:36:31 by jwe]
jwe
parents: 3457
diff changeset
72 elseif (nargin == 1)
2303
5cffc4b8de57 [project @ 1996-06-24 09:15:24 by jwe]
jwe
parents: 1337
diff changeset
73 ## Response of an FIR filter.
3893
abd8659eea11 [project @ 2002-04-09 21:36:31 by jwe]
jwe
parents: 3457
diff changeset
74 a = n = region = Fs = [];
566
8529a21443fa [project @ 1994-07-26 02:07:38 by jwe]
jwe
parents:
diff changeset
75 elseif (nargin == 2)
2303
5cffc4b8de57 [project @ 1996-06-24 09:15:24 by jwe]
jwe
parents: 1337
diff changeset
76 ## Response of an IIR filter
3893
abd8659eea11 [project @ 2002-04-09 21:36:31 by jwe]
jwe
parents: 3457
diff changeset
77 n = region = Fs = [];
566
8529a21443fa [project @ 1994-07-26 02:07:38 by jwe]
jwe
parents:
diff changeset
78 elseif (nargin == 3)
3893
abd8659eea11 [project @ 2002-04-09 21:36:31 by jwe]
jwe
parents: 3457
diff changeset
79 region = Fs = [];
566
8529a21443fa [project @ 1994-07-26 02:07:38 by jwe]
jwe
parents:
diff changeset
80 elseif (nargin == 4)
3893
abd8659eea11 [project @ 2002-04-09 21:36:31 by jwe]
jwe
parents: 3457
diff changeset
81 Fs = [];
5443
ec8c33dcd1bf [project @ 2005-09-08 01:40:57 by jwe]
jwe
parents: 5378
diff changeset
82 if (! ischar (region) && ! isempty (region))
11587
c792872f8942 all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents: 11523
diff changeset
83 Fs = region;
3893
abd8659eea11 [project @ 2002-04-09 21:36:31 by jwe]
jwe
parents: 3457
diff changeset
84 region = [];
abd8659eea11 [project @ 2002-04-09 21:36:31 by jwe]
jwe
parents: 3457
diff changeset
85 endif
abd8659eea11 [project @ 2002-04-09 21:36:31 by jwe]
jwe
parents: 3457
diff changeset
86 endif
abd8659eea11 [project @ 2002-04-09 21:36:31 by jwe]
jwe
parents: 3457
diff changeset
87
5377
bd4ee620c38a [project @ 2005-06-02 15:42:39 by jwe]
jwe
parents: 5307
diff changeset
88 if (isempty (b))
bd4ee620c38a [project @ 2005-06-02 15:42:39 by jwe]
jwe
parents: 5307
diff changeset
89 b = 1;
bd4ee620c38a [project @ 2005-06-02 15:42:39 by jwe]
jwe
parents: 5307
diff changeset
90 endif
11587
c792872f8942 all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents: 11523
diff changeset
91 if (isempty (a))
c792872f8942 all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents: 11523
diff changeset
92 a = 1;
3893
abd8659eea11 [project @ 2002-04-09 21:36:31 by jwe]
jwe
parents: 3457
diff changeset
93 endif
abd8659eea11 [project @ 2002-04-09 21:36:31 by jwe]
jwe
parents: 3457
diff changeset
94 if (isempty (n))
11587
c792872f8942 all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents: 11523
diff changeset
95 n = 512;
3893
abd8659eea11 [project @ 2002-04-09 21:36:31 by jwe]
jwe
parents: 3457
diff changeset
96 endif
abd8659eea11 [project @ 2002-04-09 21:36:31 by jwe]
jwe
parents: 3457
diff changeset
97 if (isempty (region))
abd8659eea11 [project @ 2002-04-09 21:36:31 by jwe]
jwe
parents: 3457
diff changeset
98 if (isreal (b) && isreal (a))
abd8659eea11 [project @ 2002-04-09 21:36:31 by jwe]
jwe
parents: 3457
diff changeset
99 region = "half";
abd8659eea11 [project @ 2002-04-09 21:36:31 by jwe]
jwe
parents: 3457
diff changeset
100 else
abd8659eea11 [project @ 2002-04-09 21:36:31 by jwe]
jwe
parents: 3457
diff changeset
101 region = "whole";
abd8659eea11 [project @ 2002-04-09 21:36:31 by jwe]
jwe
parents: 3457
diff changeset
102 endif
abd8659eea11 [project @ 2002-04-09 21:36:31 by jwe]
jwe
parents: 3457
diff changeset
103 endif
11587
c792872f8942 all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents: 11523
diff changeset
104 if (isempty (Fs))
c792872f8942 all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents: 11523
diff changeset
105 if (nargout == 0)
c792872f8942 all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents: 11523
diff changeset
106 Fs = 2;
c792872f8942 all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents: 11523
diff changeset
107 else
c792872f8942 all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents: 11523
diff changeset
108 Fs = 2*pi;
3893
abd8659eea11 [project @ 2002-04-09 21:36:31 by jwe]
jwe
parents: 3457
diff changeset
109 endif
566
8529a21443fa [project @ 1994-07-26 02:07:38 by jwe]
jwe
parents:
diff changeset
110 endif
8529a21443fa [project @ 1994-07-26 02:07:38 by jwe]
jwe
parents:
diff changeset
111
5583
8eebdcfde94e [project @ 2005-12-15 01:16:26 by jwe]
jwe
parents: 5443
diff changeset
112 a = a(:);
8eebdcfde94e [project @ 2005-12-15 01:16:26 by jwe]
jwe
parents: 5443
diff changeset
113 b = b(:);
566
8529a21443fa [project @ 1994-07-26 02:07:38 by jwe]
jwe
parents:
diff changeset
114
8506
bc982528de11 comment style fixes
John W. Eaton <jwe@octave.org>
parents: 7151
diff changeset
115 if (! isscalar (n))
bc982528de11 comment style fixes
John W. Eaton <jwe@octave.org>
parents: 7151
diff changeset
116 ## Explicit frequency vector given
5377
bd4ee620c38a [project @ 2005-06-02 15:42:39 by jwe]
jwe
parents: 5307
diff changeset
117 w = f = n;
8506
bc982528de11 comment style fixes
John W. Eaton <jwe@octave.org>
parents: 7151
diff changeset
118 if (nargin == 4)
bc982528de11 comment style fixes
John W. Eaton <jwe@octave.org>
parents: 7151
diff changeset
119 ## Sampling rate Fs was specified
5377
bd4ee620c38a [project @ 2005-06-02 15:42:39 by jwe]
jwe
parents: 5307
diff changeset
120 w = 2*pi*f/Fs;
566
8529a21443fa [project @ 1994-07-26 02:07:38 by jwe]
jwe
parents:
diff changeset
121 endif
5986
14078139f941 [project @ 2006-09-12 02:46:33 by jwe]
jwe
parents: 5583
diff changeset
122 k = max (length (b), length (a));
14078139f941 [project @ 2006-09-12 02:46:33 by jwe]
jwe
parents: 5583
diff changeset
123 hb = polyval (postpad (b, k), exp (j*w));
14078139f941 [project @ 2006-09-12 02:46:33 by jwe]
jwe
parents: 5583
diff changeset
124 ha = polyval (postpad (a, k), exp (j*w));
8667
a89198168175 freqz.m: freqz.m: fix for long input
Ben Abbott <bpabbott@mac.com>
parents: 8517
diff changeset
125 else
5583
8eebdcfde94e [project @ 2005-12-15 01:16:26 by jwe]
jwe
parents: 5443
diff changeset
126 ## polyval(fliplr(P),exp(jw)) is O(p n) and fft(x) is O(n log(n)),
6653
673686daec87 [project @ 2007-05-22 15:36:09 by jwe]
jwe
parents: 6046
diff changeset
127 ## where p is the order of the polynomial P. For small p it
5377
bd4ee620c38a [project @ 2005-06-02 15:42:39 by jwe]
jwe
parents: 5307
diff changeset
128 ## would be faster to use polyval but in practice the overhead for
bd4ee620c38a [project @ 2005-06-02 15:42:39 by jwe]
jwe
parents: 5307
diff changeset
129 ## polyval is much higher and the little bit of time saved isn't
bd4ee620c38a [project @ 2005-06-02 15:42:39 by jwe]
jwe
parents: 5307
diff changeset
130 ## worth the extra code.
8667
a89198168175 freqz.m: freqz.m: fix for long input
Ben Abbott <bpabbott@mac.com>
parents: 8517
diff changeset
131 k = max (length (b), length (a));
a89198168175 freqz.m: freqz.m: fix for long input
Ben Abbott <bpabbott@mac.com>
parents: 8517
diff changeset
132 if (k > n/2 && nargout == 0)
a89198168175 freqz.m: freqz.m: fix for long input
Ben Abbott <bpabbott@mac.com>
parents: 8517
diff changeset
133 ## Ensure a causal phase response.
a89198168175 freqz.m: freqz.m: fix for long input
Ben Abbott <bpabbott@mac.com>
parents: 8517
diff changeset
134 n = n * 2 .^ ceil (log2 (2*k/n));
a89198168175 freqz.m: freqz.m: fix for long input
Ben Abbott <bpabbott@mac.com>
parents: 8517
diff changeset
135 endif
a89198168175 freqz.m: freqz.m: fix for long input
Ben Abbott <bpabbott@mac.com>
parents: 8517
diff changeset
136
a89198168175 freqz.m: freqz.m: fix for long input
Ben Abbott <bpabbott@mac.com>
parents: 8517
diff changeset
137 if (strcmp (region, "whole"))
a89198168175 freqz.m: freqz.m: fix for long input
Ben Abbott <bpabbott@mac.com>
parents: 8517
diff changeset
138 N = n;
a89198168175 freqz.m: freqz.m: fix for long input
Ben Abbott <bpabbott@mac.com>
parents: 8517
diff changeset
139 else
a89198168175 freqz.m: freqz.m: fix for long input
Ben Abbott <bpabbott@mac.com>
parents: 8517
diff changeset
140 N = 2*n;
a89198168175 freqz.m: freqz.m: fix for long input
Ben Abbott <bpabbott@mac.com>
parents: 8517
diff changeset
141 endif
a89198168175 freqz.m: freqz.m: fix for long input
Ben Abbott <bpabbott@mac.com>
parents: 8517
diff changeset
142
a89198168175 freqz.m: freqz.m: fix for long input
Ben Abbott <bpabbott@mac.com>
parents: 8517
diff changeset
143 f = Fs * (0:n-1).' / N;
a89198168175 freqz.m: freqz.m: fix for long input
Ben Abbott <bpabbott@mac.com>
parents: 8517
diff changeset
144
a89198168175 freqz.m: freqz.m: fix for long input
Ben Abbott <bpabbott@mac.com>
parents: 8517
diff changeset
145 pad_sz = N*ceil (k/N);
a89198168175 freqz.m: freqz.m: fix for long input
Ben Abbott <bpabbott@mac.com>
parents: 8517
diff changeset
146 b = postpad (b, pad_sz);
a89198168175 freqz.m: freqz.m: fix for long input
Ben Abbott <bpabbott@mac.com>
parents: 8517
diff changeset
147 a = postpad (a, pad_sz);
a89198168175 freqz.m: freqz.m: fix for long input
Ben Abbott <bpabbott@mac.com>
parents: 8517
diff changeset
148
a89198168175 freqz.m: freqz.m: fix for long input
Ben Abbott <bpabbott@mac.com>
parents: 8517
diff changeset
149 hb = zeros (n, 1);
a89198168175 freqz.m: freqz.m: fix for long input
Ben Abbott <bpabbott@mac.com>
parents: 8517
diff changeset
150 ha = zeros (n, 1);
a89198168175 freqz.m: freqz.m: fix for long input
Ben Abbott <bpabbott@mac.com>
parents: 8517
diff changeset
151
a89198168175 freqz.m: freqz.m: fix for long input
Ben Abbott <bpabbott@mac.com>
parents: 8517
diff changeset
152 for i = 1:N:pad_sz
a89198168175 freqz.m: freqz.m: fix for long input
Ben Abbott <bpabbott@mac.com>
parents: 8517
diff changeset
153 hb = hb + fft (postpad (b(i:i+N-1), N))(1:n);
a89198168175 freqz.m: freqz.m: fix for long input
Ben Abbott <bpabbott@mac.com>
parents: 8517
diff changeset
154 ha = ha + fft (postpad (a(i:i+N-1), N))(1:n);
a89198168175 freqz.m: freqz.m: fix for long input
Ben Abbott <bpabbott@mac.com>
parents: 8517
diff changeset
155 endfor
a89198168175 freqz.m: freqz.m: fix for long input
Ben Abbott <bpabbott@mac.com>
parents: 8517
diff changeset
156
566
8529a21443fa [project @ 1994-07-26 02:07:38 by jwe]
jwe
parents:
diff changeset
157 endif
8529a21443fa [project @ 1994-07-26 02:07:38 by jwe]
jwe
parents:
diff changeset
158
3893
abd8659eea11 [project @ 2002-04-09 21:36:31 by jwe]
jwe
parents: 3457
diff changeset
159 h = hb ./ ha;
abd8659eea11 [project @ 2002-04-09 21:36:31 by jwe]
jwe
parents: 3457
diff changeset
160
8506
bc982528de11 comment style fixes
John W. Eaton <jwe@octave.org>
parents: 7151
diff changeset
161 if (nargout != 0)
bc982528de11 comment style fixes
John W. Eaton <jwe@octave.org>
parents: 7151
diff changeset
162 ## Return values and don't plot.
3907
437884fae441 [project @ 2002-04-24 19:33:08 by jwe]
jwe
parents: 3893
diff changeset
163 h_r = h;
5377
bd4ee620c38a [project @ 2005-06-02 15:42:39 by jwe]
jwe
parents: 5307
diff changeset
164 f_r = f;
8506
bc982528de11 comment style fixes
John W. Eaton <jwe@octave.org>
parents: 7151
diff changeset
165 else
bc982528de11 comment style fixes
John W. Eaton <jwe@octave.org>
parents: 7151
diff changeset
166 ## Plot and don't return values.
5377
bd4ee620c38a [project @ 2005-06-02 15:42:39 by jwe]
jwe
parents: 5307
diff changeset
167 freqz_plot (f, h);
7151
aeeb646f6538 [project @ 2007-11-09 19:34:17 by jwe]
jwe
parents: 7017
diff changeset
168 endif
3907
437884fae441 [project @ 2002-04-24 19:33:08 by jwe]
jwe
parents: 3893
diff changeset
169
566
8529a21443fa [project @ 1994-07-26 02:07:38 by jwe]
jwe
parents:
diff changeset
170 endfunction
5377
bd4ee620c38a [project @ 2005-06-02 15:42:39 by jwe]
jwe
parents: 5307
diff changeset
171
14363
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
172
5377
bd4ee620c38a [project @ 2005-06-02 15:42:39 by jwe]
jwe
parents: 5307
diff changeset
173 %!test # correct values and fft-polyval consistency
bd4ee620c38a [project @ 2005-06-02 15:42:39 by jwe]
jwe
parents: 5307
diff changeset
174 %! # butterworth filter, order 2, cutoff pi/2 radians
bd4ee620c38a [project @ 2005-06-02 15:42:39 by jwe]
jwe
parents: 5307
diff changeset
175 %! b = [0.292893218813452 0.585786437626905 0.292893218813452];
bd4ee620c38a [project @ 2005-06-02 15:42:39 by jwe]
jwe
parents: 5307
diff changeset
176 %! a = [1 0 0.171572875253810];
14363
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
177 %! [h,w] = freqz (b,a,32);
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
178 %! assert (h(1),1,10*eps);
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
179 %! assert (abs (h(17)).^2,0.5,10*eps);
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
180 %! assert (h,freqz (b,a,w),10*eps); # fft should be consistent with polyval
5377
bd4ee620c38a [project @ 2005-06-02 15:42:39 by jwe]
jwe
parents: 5307
diff changeset
181
bd4ee620c38a [project @ 2005-06-02 15:42:39 by jwe]
jwe
parents: 5307
diff changeset
182 %!test # whole-half consistency
bd4ee620c38a [project @ 2005-06-02 15:42:39 by jwe]
jwe
parents: 5307
diff changeset
183 %! b = [1 1 1]/3; # 3-sample average
14363
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
184 %! [h,w] = freqz (b,1,32,"whole");
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
185 %! assert (h(2:16),conj (h(32:-1:18)),20*eps);
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
186 %! [h2,w2] = freqz (b,1,16,"half");
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
187 %! assert (h(1:16),h2,20*eps);
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
188 %! assert (w(1:16),w2,20*eps);
5377
bd4ee620c38a [project @ 2005-06-02 15:42:39 by jwe]
jwe
parents: 5307
diff changeset
189
bd4ee620c38a [project @ 2005-06-02 15:42:39 by jwe]
jwe
parents: 5307
diff changeset
190 %!test # Sampling frequency properly interpreted
5986
14078139f941 [project @ 2006-09-12 02:46:33 by jwe]
jwe
parents: 5583
diff changeset
191 %! b = [1 1 1]/3; a = [1 0.2];
14363
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
192 %! [h,f] = freqz (b,a,16,320);
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
193 %! assert (f,[0:15]'*10,10*eps);
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
194 %! [h2,f2] = freqz (b,a,[0:15]*10,320);
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
195 %! assert (f2,[0:15]*10,10*eps);
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
196 %! assert (h,h2.',20*eps);
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
197 %! [h3,f3] = freqz (b,a,32,"whole",320);
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
198 %! assert (f3,[0:31]'*10,10*eps);