view scripts/linear-algebra/qzhess.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 4d917a6a858b
children 53eaa83e4181
line wrap: on
line source

## Copyright (C) 1993-2012 John W. Eaton
##
## 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
## <http://www.gnu.org/licenses/>.

## -*- texinfo -*-
## @deftypefn {Function File} {[@var{aa}, @var{bb}, @var{q}, @var{z}] =} qzhess (@var{A}, @var{B})
## Compute the Hessenberg-triangular decomposition of the matrix pencil
## @code{(@var{A}, @var{B})}, returning
## @code{@var{aa} = @var{q} * @var{A} * @var{z}},
## @code{@var{bb} = @var{q} * @var{B} * @var{z}}, with @var{q} and @var{z}
## orthogonal.  For example:
##
## @example
## @group
## [aa, bb, q, z] = qzhess ([1, 2; 3, 4], [5, 6; 7, 8])
##      @result{} aa = [ -3.02244, -4.41741;  0.92998,  0.69749 ]
##      @result{} bb = [ -8.60233, -9.99730;  0.00000, -0.23250 ]
##      @result{}  q = [ -0.58124, -0.81373; -0.81373,  0.58124 ]
##      @result{}  z = [ 1, 0; 0, 1 ]
## @end group
## @end example
##
## The Hessenberg-triangular decomposition is the first step in
## Moler and Stewart's QZ@tie{}decomposition algorithm.
##
## Algorithm taken from Golub and Van Loan,
## @cite{Matrix Computations, 2nd edition}.
## @end deftypefn

## Author: A. S. Hodel <scotte@eng.auburn.edu>
## Created: August 1993
## Adapted-By: jwe

function [aa, bb, q, z] = qzhess (A, B)

  if (nargin != 2)
    print_usage ();
  endif

  [na, ma] = size (A);
  [nb, mb] = size (B);
  if (na != ma || na != nb || nb != mb)
    error ("qzhess: incompatible dimensions");
  endif

  ## Reduce to hessenberg-triangular form.

  [q, bb] = qr (B);
  aa = q' * A;
  q = q';
  z = eye (na);
  for j = 1:(na-2)
    for i = na:-1:(j+2)

      ## disp (["zero out aa(", num2str(i), ",", num2str(j), ")"])

      rot = givens (aa (i-1, j), aa (i, j));
      aa ((i-1):i, :) = rot *aa ((i-1):i, :);
      bb ((i-1):i, :) = rot *bb ((i-1):i, :);
      q  ((i-1):i, :) = rot *q  ((i-1):i, :);

      ## disp (["now zero out bb(", num2str(i), ",", num2str(i-1), ")"])

      rot = givens (bb (i, i), bb (i, i-1))';
      bb (:, (i-1):i) = bb (:, (i-1):i) * rot';
      aa (:, (i-1):i) = aa (:, (i-1):i) * rot';
      z  (:, (i-1):i) = z  (:, (i-1):i) * rot';

    endfor
  endfor

  bb (2, 1) = 0.0;
  for i = 3:na
    bb (i, 1:(i-1)) = zeros (1, i-1);
    aa (i, 1:(i-2)) = zeros (1, i-2);
  endfor

endfunction


%!test
%! a = [1 2 1 3;
%!      2 5 3 2;
%!      5 5 1 0;
%!      4 0 3 2];
%! b = [0 4 2 1;
%!      2 3 1 1;
%!      1 0 2 1;
%!      2 5 3 2];
%! mask = [0 0 0 0;
%!         0 0 0 0;
%!         1 0 0 0;
%!         1 1 0 0];
%! [aa, bb, q, z] = qzhess (a, b);
%! assert (inv (q) - q', zeros (4), 2e-8);
%! assert (inv (z) - z', zeros (4), 2e-8);
%! assert (q * a * z, aa, 2e-8);
%! assert (aa .* mask, zeros (4), 2e-8);
%! assert (q * b * z, bb, 2e-8);
%! assert (bb .* mask, zeros (4), 2e-8);

%!test
%! a = [1 2 3 4 5;
%!      3 2 3 1 0;
%!      4 3 2 1 1;
%!      0 1 0 1 0;
%!      3 2 1 0 5];
%! b = [5 0 4 0 1;
%!      1 1 1 2 5;
%!      0 3 2 1 0;
%!      4 3 0 3 5;
%!      2 1 2 1 3];
%! mask = [0 0 0 0 0;
%!         0 0 0 0 0;
%!         1 0 0 0 0;
%!         1 1 0 0 0;
%!         1 1 1 0 0];
%! [aa, bb, q, z] = qzhess (a, b);
%! assert (inv (q) - q', zeros (5), 2e-8);
%! assert (inv (z) - z', zeros (5), 2e-8);
%! assert (q * a * z, aa, 2e-8);
%! assert (aa .* mask, zeros (5), 2e-8);
%! assert (q * b * z, bb, 2e-8);
%! assert (bb .* mask, zeros (5), 2e-8);

%!error qzhess ([0])
%!error qzhess ()