Mercurial > hg > octave-nkf
annotate scripts/linear-algebra/rref.m @ 14868:5d3a684236b0
maint: Use Octave coding conventions for cuddling parentheses in scripts directory
* lin2mu.m, loadaudio.m, wavread.m, accumarray.m, bicubic.m, celldisp.m,
colon.m, cplxpair.m, dblquad.m, divergence.m, genvarname.m, gradient.m,
int2str.m, interp1.m, interp1q.m, interp2.m, interpn.m, loadobj.m, nthargout.m,
__isequal__.m, __splinen__.m, quadgk.m, quadl.m, quadv.m, rat.m, rot90.m,
rotdim.m, saveobj.m, subsindex.m, triplequad.m, delaunay3.m, griddata.m,
inpolygon.m, tsearchn.m, voronoi.m, get_first_help_sentence.m, which.m,
gray2ind.m, pink.m, dlmwrite.m, strread.m, textread.m, textscan.m, housh.m,
ishermitian.m, issymmetric.m, krylov.m, logm.m, null.m, rref.m,
compare_versions.m, copyfile.m, dump_prefs.m, edit.m, fileparts.m,
getappdata.m, isappdata.m, movefile.m, orderfields.m, parseparams.m,
__xzip__.m, rmappdata.m, setappdata.m, swapbytes.m, unpack.m, ver.m, fminbnd.m,
fminunc.m, fsolve.m, glpk.m, lsqnonneg.m, qp.m, sqp.m, configure_make.m,
copy_files.m, describe.m, get_description.m, get_forge_pkg.m, install.m,
installed_packages.m, is_architecture_dependent.m, load_package_dirs.m,
print_package_description.m, rebuild.m, repackage.m, save_order.m, shell.m,
allchild.m, ancestor.m, area.m, axes.m, axis.m, clabel.m, close.m, colorbar.m,
comet.m, comet3.m, contour.m, cylinder.m, ezmesh.m, ezsurf.m, findobj.m,
fplot.m, hist.m, isocolors.m, isonormals.m, isosurface.m, isprop.m, legend.m,
mesh.m, meshz.m, pareto.m, pcolor.m, peaks.m, plot3.m, plotmatrix.m, plotyy.m,
polar.m, print.m, __add_datasource__.m, __add_default_menu__.m,
__axes_limits__.m, __bar__.m, __clabel__.m, __contour__.m, __errcomm__.m,
__errplot__.m, __ezplot__.m, __file_filter__.m, __fltk_print__.m,
__ghostscript__.m, __gnuplot_print__.m, __go_draw_axes__.m,
__go_draw_figure__.m, __interp_cube__.m, __marching_cube__.m, __patch__.m,
__pie__.m, __plt__.m, __print_parse_opts__.m, __quiver__.m, __scatter__.m,
__stem__.m, __tight_eps_bbox__.m, __uigetdir_fltk__.m, __uigetfile_fltk__.m,
__uiputfile_fltk__.m, quiver.m, quiver3.m, rectangle.m, refreshdata.m,
ribbon.m, scatter.m, semilogy.m, shading.m, slice.m, subplot.m, surface.m,
surfl.m, surfnorm.m, text.m, uigetfile.m, uiputfile.m, whitebg.m, deconv.m,
mkpp.m, pchip.m, polyaffine.m, polyder.m, polygcd.m, polyout.m, polyval.m,
ppint.m, ppjumps.m, ppval.m, residue.m, roots.m, spline.m, splinefit.m,
addpref.m, getpref.m, setpref.m, ismember.m, setxor.m, arch_fit.m, arch_rnd.m,
arch_test.m, autoreg_matrix.m, diffpara.m, fftconv.m, filter2.m, hanning.m,
hurst.m, periodogram.m, triangle_sw.m, sinc.m, spectral_xdf.m, spencer.m,
stft.m, synthesis.m, unwrap.m, yulewalker.m, bicgstab.m, gmres.m, pcg.m, pcr.m,
__sprand_impl__.m, speye.m, spfun.m, sprandn.m, spstats.m, svds.m,
treelayout.m, treeplot.m, bessel.m, factor.m, legendre.m, perms.m, primes.m,
magic.m, toeplitz.m, corr.m, cov.m, mean.m, median.m, mode.m, qqplot.m,
quantile.m, ranks.m, zscore.m, logistic_regression_likelihood.m,
bartlett_test.m, chisquare_test_homogeneity.m, chisquare_test_independence.m,
kolmogorov_smirnov_test.m, run_test.m, u_test.m, wilcoxon_test.m, z_test.m,
z_test_2.m, bin2dec.m, dec2base.m, mat2str.m, strcat.m, strchr.m, strjust.m,
strtok.m, substr.m, untabify.m, assert.m, demo.m, example.m, fail.m, speed.m,
test.m, now.m: Use Octave coding conventions for cuddling parentheses in
scripts directory.
author | Rik <octave@nomad.inbox5.com> |
---|---|
date | Tue, 17 Jul 2012 07:08:39 -0700 |
parents | f3d52523cde1 |
children | d63878346099 |
rev | line source |
---|---|
14138
72c96de7a403
maint: update copyright notices for 2012
John W. Eaton <jwe@octave.org>
parents:
13048
diff
changeset
|
1 ## Copyright (C) 2000-2012 Paul Kienzle |
5827 | 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 | |
7016 | 7 ## the Free Software Foundation; either version 3 of the License, or (at |
8 ## your option) any later version. | |
5827 | 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 | |
7016 | 16 ## along with Octave; see the file COPYING. If not, see |
17 ## <http://www.gnu.org/licenses/>. | |
5827 | 18 |
19 ## -*- texinfo -*- | |
12584
7ef7e20057fa
Improve documentation strings in Linear Algebra chapter.
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
20 ## @deftypefn {Function File} {} rref (@var{A}) |
7ef7e20057fa
Improve documentation strings in Linear Algebra chapter.
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
21 ## @deftypefnx {Function File} {} rref (@var{A}, @var{tol}) |
7ef7e20057fa
Improve documentation strings in Linear Algebra chapter.
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
22 ## @deftypefnx {Function File} {[@var{r}, @var{k}] =} rref (@dots{}) |
7ef7e20057fa
Improve documentation strings in Linear Algebra chapter.
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
23 ## Return the reduced row echelon form of @var{A}. @var{tol} defaults |
11469
c776f063fefe
Overhaul m-script files to use common variable name between code and documentation.
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
24 ## to @code{eps * max (size (@var{A})) * norm (@var{A}, inf)}. |
5827 | 25 ## |
26 ## Called with two return arguments, @var{k} returns the vector of | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
27 ## "bound variables", which are those columns on which elimination |
5827 | 28 ## has been performed. |
29 ## | |
30 ## @end deftypefn | |
31 | |
32 ## Author: Paul Kienzle <pkienzle@users.sf.net> | |
11469
c776f063fefe
Overhaul m-script files to use common variable name between code and documentation.
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
33 ## (based on an anonymous source from the public domain) |
5827 | 34 |
11469
c776f063fefe
Overhaul m-script files to use common variable name between code and documentation.
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
35 function [A, k] = rref (A, tol) |
5827 | 36 |
37 if (nargin < 1 || nargin > 2) | |
38 print_usage (); | |
39 endif | |
40 | |
41 if (ndims (A) > 2) | |
42 error ("rref: expecting matrix argument"); | |
43 endif | |
44 | |
45 [rows, cols] = size (A); | |
46 | |
47 if (nargin < 2) | |
7795
df9519e9990c
Handle single precision eps values
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
48 if (isa (A, "single")) |
11469
c776f063fefe
Overhaul m-script files to use common variable name between code and documentation.
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
49 tol = eps ("single") * max (rows, cols) * norm (A, inf ("single")); |
7795
df9519e9990c
Handle single precision eps values
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
50 else |
11469
c776f063fefe
Overhaul m-script files to use common variable name between code and documentation.
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
51 tol = eps * max (rows, cols) * norm (A, inf); |
7795
df9519e9990c
Handle single precision eps values
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
52 endif |
5827 | 53 endif |
54 | |
55 used = zeros (1, cols); | |
56 r = 1; | |
57 for c = 1:cols | |
58 ## Find the pivot row | |
59 [m, pivot] = max (abs (A(r:rows,c))); | |
60 pivot = r + pivot - 1; | |
61 | |
11469
c776f063fefe
Overhaul m-script files to use common variable name between code and documentation.
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
62 if (m <= tol) |
5827 | 63 ## Skip column c, making sure the approximately zero terms are |
64 ## actually zero. | |
14868
5d3a684236b0
maint: Use Octave coding conventions for cuddling parentheses in scripts directory
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
65 A(r:rows, c) = zeros (rows-r+1, 1); |
5827 | 66 else |
67 ## keep track of bound variables | |
68 used (1, c) = 1; | |
69 | |
70 ## Swap current row and pivot row | |
14868
5d3a684236b0
maint: Use Octave coding conventions for cuddling parentheses in scripts directory
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
71 A([pivot, r], c:cols) = A([r, pivot], c:cols); |
5827 | 72 |
73 ## Normalize pivot row | |
14868
5d3a684236b0
maint: Use Octave coding conventions for cuddling parentheses in scripts directory
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
74 A(r, c:cols) = A(r, c:cols) / A(r, c); |
5827 | 75 |
76 ## Eliminate the current column | |
77 ridx = [1:r-1, r+1:rows]; | |
14868
5d3a684236b0
maint: Use Octave coding conventions for cuddling parentheses in scripts directory
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
78 A(ridx, c:cols) = A(ridx, c:cols) - A(ridx, c) * A(r, c:cols); |
5827 | 79 |
80 ## Check if done | |
81 if (r++ == rows) | |
10549 | 82 break; |
5827 | 83 endif |
84 endif | |
85 endfor | |
86 k = find (used); | |
87 | |
88 endfunction | |
13048
c5c94b63931f
codesprint: linear algebra tests: cross, housh, planerot, qzhess, rref
Roman Belov <romblv@gmail.com>
parents:
12584
diff
changeset
|
89 |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
90 |
13048
c5c94b63931f
codesprint: linear algebra tests: cross, housh, planerot, qzhess, rref
Roman Belov <romblv@gmail.com>
parents:
12584
diff
changeset
|
91 %!test |
c5c94b63931f
codesprint: linear algebra tests: cross, housh, planerot, qzhess, rref
Roman Belov <romblv@gmail.com>
parents:
12584
diff
changeset
|
92 %! a = [1]; |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
93 %! [r k] = rref (a); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
94 %! assert (r, [1], 2e-8); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
95 %! assert (k, [1], 2e-8); |
13048
c5c94b63931f
codesprint: linear algebra tests: cross, housh, planerot, qzhess, rref
Roman Belov <romblv@gmail.com>
parents:
12584
diff
changeset
|
96 |
c5c94b63931f
codesprint: linear algebra tests: cross, housh, planerot, qzhess, rref
Roman Belov <romblv@gmail.com>
parents:
12584
diff
changeset
|
97 %!test |
c5c94b63931f
codesprint: linear algebra tests: cross, housh, planerot, qzhess, rref
Roman Belov <romblv@gmail.com>
parents:
12584
diff
changeset
|
98 %! a = [1 3; 4 5]; |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
99 %! [r k] = rref (a); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
100 %! assert (rank (a), rank (r), 2e-8); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
101 %! assert (r, eye (2), 2e-8); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
102 %! assert (k == [1, 2] || k == [2, 1]); |
13048
c5c94b63931f
codesprint: linear algebra tests: cross, housh, planerot, qzhess, rref
Roman Belov <romblv@gmail.com>
parents:
12584
diff
changeset
|
103 |
c5c94b63931f
codesprint: linear algebra tests: cross, housh, planerot, qzhess, rref
Roman Belov <romblv@gmail.com>
parents:
12584
diff
changeset
|
104 %!test |
c5c94b63931f
codesprint: linear algebra tests: cross, housh, planerot, qzhess, rref
Roman Belov <romblv@gmail.com>
parents:
12584
diff
changeset
|
105 %! a = [1 3; 4 5; 7 9]; |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
106 %! [r k] = rref (a); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
107 %! assert (rank (a), rank (r), 2e-8); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
108 %! assert (r, eye(3)(:,1:2), 2e-8); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
109 %! assert (k, [1 2], 2e-8); |
13048
c5c94b63931f
codesprint: linear algebra tests: cross, housh, planerot, qzhess, rref
Roman Belov <romblv@gmail.com>
parents:
12584
diff
changeset
|
110 |
c5c94b63931f
codesprint: linear algebra tests: cross, housh, planerot, qzhess, rref
Roman Belov <romblv@gmail.com>
parents:
12584
diff
changeset
|
111 %!test |
c5c94b63931f
codesprint: linear algebra tests: cross, housh, planerot, qzhess, rref
Roman Belov <romblv@gmail.com>
parents:
12584
diff
changeset
|
112 %! a = [1 2 3; 2 4 6; 7 2 0]; |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
113 %! [r k] = rref (a); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
114 %! assert (rank (a), rank (r), 2e-8); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
115 %! assert (r, [1 0 (3-7/2); 0 1 (7/4); 0 0 0], 2e-8); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
116 %! assert (k, [1 2], 2e-8); |
13048
c5c94b63931f
codesprint: linear algebra tests: cross, housh, planerot, qzhess, rref
Roman Belov <romblv@gmail.com>
parents:
12584
diff
changeset
|
117 |
c5c94b63931f
codesprint: linear algebra tests: cross, housh, planerot, qzhess, rref
Roman Belov <romblv@gmail.com>
parents:
12584
diff
changeset
|
118 %!test |
c5c94b63931f
codesprint: linear algebra tests: cross, housh, planerot, qzhess, rref
Roman Belov <romblv@gmail.com>
parents:
12584
diff
changeset
|
119 %! a = [1 2 1; 2 4 2.01; 2 4 2.1]; |
c5c94b63931f
codesprint: linear algebra tests: cross, housh, planerot, qzhess, rref
Roman Belov <romblv@gmail.com>
parents:
12584
diff
changeset
|
120 %! tol = 0.02; |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
121 %! [r k] = rref (a, tol); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
122 %! assert (rank (a, tol), rank (r, tol), 2e-8); |
13048
c5c94b63931f
codesprint: linear algebra tests: cross, housh, planerot, qzhess, rref
Roman Belov <romblv@gmail.com>
parents:
12584
diff
changeset
|
123 %! tol = 0.2; |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
124 %! [r k] = rref (a, tol); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
125 %! assert (rank (a, tol), rank (r, tol), 2e-8); |
13048
c5c94b63931f
codesprint: linear algebra tests: cross, housh, planerot, qzhess, rref
Roman Belov <romblv@gmail.com>
parents:
12584
diff
changeset
|
126 |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
127 %!error rref (); |
13048
c5c94b63931f
codesprint: linear algebra tests: cross, housh, planerot, qzhess, rref
Roman Belov <romblv@gmail.com>
parents:
12584
diff
changeset
|
128 |