Mercurial > hg > octave-lyh
view scripts/sparse/gmres.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 | b76f0740940e |
line wrap: on
line source
## Copyright (C) 2009-2012 Carlo de Falco ## ## 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{x} =} gmres (@var{A}, @var{b}, @var{m}, @var{rtol}, @var{maxit}, @var{M1}, @var{M2}, @var{x0}) ## @deftypefnx {Function File} {@var{x} =} gmres (@var{A}, @var{b}, @var{m}, @var{rtol}, @var{maxit}, @var{P}) ## @deftypefnx {Function File} {[@var{x}, @var{flag}, @var{relres}, @var{iter}, @var{resvec}] =} gmres (@dots{}) ## Solve @code{A x = b} using the Preconditioned GMRES iterative method ## with restart, a.k.a. PGMRES(m). ## ## @itemize @minus ## @item @var{rtol} is the relative tolerance, ## if not given or set to [] the default value 1e-6 is used. ## ## @item @var{maxit} is the maximum number of outer iterations, ## if not given or set to [] the default value ## @code{min (10, numel (b) / restart)} is used. ## ## @item @var{x0} is the initial guess, ## if not given or set to [] the default value @code{zeros(size (b))} is used. ## ## @item @var{m} is the restart parameter, ## if not given or set to [] the default value @code{numel (b)} is used. ## @end itemize ## ## Argument @var{A} can be passed as a matrix, function handle, or ## inline function @code{f} such that @code{f(x) = A*x}. ## ## The preconditioner @var{P} is given as @code{P = M1 * M2}. ## Both @var{M1} and @var{M2} can be passed as a matrix, function handle, or ## inline function @code{g} such that @code{g(x) = M1\x} or @code{g(x) = M2\x}. ## ## Besides the vector @var{x}, additional outputs are: ## ## @itemize @minus ## @item @var{flag} indicates the exit status: ## @table @asis ## @item 0 : iteration converged to within the specified tolerance ## ## @item 1 : maximum number of iterations exceeded ## ## @item 2 : unused, but skipped for compatibility ## ## @item 3 : algorithm reached stagnation ## @end table ## ## @item @var{relres} is the final value of the relative residual. ## ## @item @var{iter} is a vector containing the number of outer iterations and ## total iterations performed. ## ## @item @var{resvec} is a vector containing the relative residual at each ## iteration. ## @end itemize ## ## @seealso{bicg, bicgstab, cgs, pcg} ## @end deftypefn function [x, flag, presn, it, resids] = gmres (A, b, restart, rtol, maxit, M1, M2, x0) if (nargin < 2 || nargin > 8) print_usage (); endif if (ischar (A)) Ax = str2func (A); elseif (ismatrix (A)) Ax = @(x) A*x; elseif (isa (A, "function_handle")) Ax = A; else error ("gmres: A must be a function or matrix"); endif if (nargin < 3 || isempty (restart)) restart = rows (b); endif if (nargin < 4 || isempty (rtol)) rtol = 1e-6; endif if (nargin < 5 || isempty (maxit)) maxit = min (rows (b)/restart, 10); endif if (nargin < 6 || isempty (M1)) M1m1x = @(x) x; elseif (ischar (M1)) M1m1x = str2func (M1); elseif (ismatrix (M1)) M1m1x = @(x) M1 \ x; elseif (isa (M1, "function_handle")) M1m1x = M1; else error ("gmres: preconditioner M1 must be a function or matrix"); endif if (nargin < 7 || isempty (M2)) M2m1x = @(x) x; elseif (ischar (M2)) M2m1x = str2func (M2); elseif (ismatrix (M2)) M2m1x = @(x) M2 \ x; elseif (isa (M2, "function_handle")) M2m1x = M2; else error ("gmres: preconditioner M2 must be a function or matrix"); endif Pm1x = @(x) M2m1x (M1m1x (x)); if (nargin < 8 || isempty (x0)) x0 = zeros (size (b)); endif x_old = x0; x = x_old; prec_res = Pm1x (b - Ax (x_old)); presn = norm (prec_res, 2); B = zeros (restart + 1, 1); V = zeros (rows (x), restart); H = zeros (restart + 1, restart); ## begin loop iter = 1; restart_it = restart + 1; resids = zeros (maxit, 1); resids(1) = presn; prec_b_norm = norm (Pm1x (b), 2); flag = 1; while (iter <= maxit * restart && presn > rtol * prec_b_norm) ## restart if (restart_it > restart) restart_it = 1; x_old = x; prec_res = Pm1x (b - Ax (x_old)); presn = norm (prec_res, 2); B(1) = presn; H(:) = 0; V(:, 1) = prec_res / presn; endif ## basic iteration tmp = Pm1x (Ax (V(:, restart_it))); [V(:,restart_it+1), H(1:restart_it+1, restart_it)] = ... mgorth (tmp, V(:,1:restart_it)); Y = (H(1:restart_it+1, 1:restart_it) \ B (1:restart_it+1)); little_res = B(1:restart_it+1) - ... H(1:restart_it+1, 1:restart_it) * Y(1:restart_it); presn = norm (little_res, 2); x = x_old + V(:, 1:restart_it) * Y(1:restart_it); resids(iter) = presn; if (norm (x - x_old, inf) <= eps) flag = 3; break endif restart_it++ ; iter++; endwhile if (presn > rtol * prec_b_norm) flag = 0; endif resids = resids(1:iter-1); it = [ceil(iter / restart), rem(iter, restart)]; endfunction %!shared A, b, dim %! dim = 100; %!test %! A = spdiags ([-ones(dim,1) 2*ones(dim,1) ones(dim,1)], [-1:1], dim, dim); %! b = ones (dim, 1); %! x = gmres (A, b, 10, 1e-10, dim, @(x) x ./ diag (A), [], b); %! assert (x, A\b, 1e-9*norm (x, Inf)); %! %!test %! x = gmres (A, b, dim, 1e-10, 1e4, @(x) diag (diag (A)) \ x, [], b); %! assert(x, A\b, 1e-7*norm (x, Inf)); %! %!test %! A = spdiags ([[1./(2:2:2*(dim-1)) 0]; 1./(1:2:2*dim-1); [0 1./(2:2:2*(dim-1))]]', -1:1, dim, dim); %! A = A'*A; %! b = rand (dim, 1); %! [x, resids] = gmres (@(x) A*x, b, dim, 1e-10, dim, @(x) x./diag (A), [], []); %! assert (x, A\b, 1e-9*norm (x, Inf)); %! x = gmres (@(x) A*x, b, dim, 1e-10, 1e6, @(x) diag (diag (A)) \ x, [], []); %! assert (x, A\b, 1e-9*norm (x, Inf)); %!test %! x = gmres (@(x) A*x, b, dim, 1e-10, 1e6, @(x) x./diag(A), [], []); %! assert (x, A\b, 1e-7*norm (x, Inf));