Home > tt2 > solve > tt_gmres.m

tt_gmres

PURPOSE ^

TT-GMRES method

SYNOPSIS ^

function [x,RESVEC,rw,rx] = tt_gmres(A, b, tol, maxout, maxin, eps_x, eps_z, M1, M2, M3, x0, verbose, varargin)

DESCRIPTION ^

TT-GMRES method
   [x] = TT_GMRES(A, B, TOL, MAXOUT, MAXIN, EPS_X, 
   EPS_Z, M1, M2, M3, X0, VERBOSE, VARARGIN)
   GMRES in TT format. Solve a system A*x=b, up to the accuracy tol,
   with maximum number of outer iterations maxout, the size of Krylov basis
   maxin, compression of solution eps_x, compression of Krylov vectors
   eps_z, optional LEFT preconditioner: [M1*[M2]*[M3]], initial guess [x0]
   (default is 1e-308*b), [verbose] - if 1 or unspecified - print messages,
   if 0 - silent mode


 TT-Toolbox 2.2, 2009-2012

This is TT Toolbox, written by Ivan Oseledets et al.
Institute of Numerical Mathematics, Moscow, Russia
webpage: http://spring.inm.ras.ru/osel

For all questions, bugs and suggestions please mail
ivan.oseledets@gmail.com
---------------------------

CROSS-REFERENCE INFORMATION ^

This function calls: This function is called by:

SOURCE CODE ^

0001 
0002 
0003 function [x,RESVEC,rw,rx] = tt_gmres(A, b, tol, maxout, maxin, eps_x, eps_z, M1, M2, M3, x0, verbose, varargin)
0004 %TT-GMRES method
0005 %   [x] = TT_GMRES(A, B, TOL, MAXOUT, MAXIN, EPS_X,
0006 %   EPS_Z, M1, M2, M3, X0, VERBOSE, VARARGIN)
0007 %   GMRES in TT format. Solve a system A*x=b, up to the accuracy tol,
0008 %   with maximum number of outer iterations maxout, the size of Krylov basis
0009 %   maxin, compression of solution eps_x, compression of Krylov vectors
0010 %   eps_z, optional LEFT preconditioner: [M1*[M2]*[M3]], initial guess [x0]
0011 %   (default is 1e-308*b), [verbose] - if 1 or unspecified - print messages,
0012 %   if 0 - silent mode
0013 %
0014 %
0015 % TT-Toolbox 2.2, 2009-2012
0016 %
0017 %This is TT Toolbox, written by Ivan Oseledets et al.
0018 %Institute of Numerical Mathematics, Moscow, Russia
0019 %webpage: http://spring.inm.ras.ru/osel
0020 %
0021 %For all questions, bugs and suggestions please mail
0022 %ivan.oseledets@gmail.com
0023 %---------------------------
0024 
0025 [atype,afun,afcnstr] = tt_iterchk(A);
0026 
0027 %%%%%%%%%%%%% parameters %%%%%%%%%%%%%%%%%%%
0028 max_swp=12; % maximum sweeps in mvk2 allowed
0029 max_sp_factor = 1.5; % maximum stagpoints (as a part of Krylov basis size) allowed
0030 max_zrank_factor = 4; % maximum rank of Krylov vectors with respect to the rank of X.
0031 derr_tol_for_sp = 1.0; % minimal jump of residual at one iteration, which is considered as
0032                        % a stagnation.
0033 use_err_trunc = 1; % use compression accuracy eps_z/err for Krylov vectors
0034 err_trunc_power = 0.8; % theoretical estimate - 1 - sometimes sucks
0035 compute_real_res = 0; % Compute the real residual on each step (expensive!)
0036 % extra param in tt_iterapp.m: matvec_alg
0037 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0038 
0039 t0 = tic;
0040 
0041 existsM1=1;
0042 existsM2=1;
0043 existsM3=1;
0044 existsx0=1;
0045 if ((nargin<8)||(isempty(M1)))
0046     existsM1=0;
0047 end;
0048 if ((nargin<9)||(isempty(M2)))
0049     existsM2=0;
0050 end;
0051 if ((nargin<10)||(isempty(M3)))
0052     existsM3=0;
0053 end;
0054 if ((nargin<11)||(isempty(x0)))
0055     existsx0=0;
0056 end;
0057 
0058 if ((nargin<12)||(isempty(verbose)))
0059     verbose=1;
0060 end;
0061 
0062 pre_b = b;
0063 max_rank=[];
0064 
0065 if (compute_real_res==1)
0066     pre_b2 = b;
0067 end;
0068 
0069 if (existsM1)
0070     [m1type,m1fun,m1fcnstr] = tt_iterchk(M1);
0071     if (compute_real_res==1)
0072     pre_b2 = tt_iterapp('mldivide',m1fun,m1type,m1fcnstr,pre_b2,min(eps_z, 0.5), max_rank, max_swp,varargin{:});
0073     end;
0074 end;
0075 if (existsM2)
0076     [m2type,m2fun,m2fcnstr] = tt_iterchk(M2);
0077     if (compute_real_res==1)
0078     pre_b2 = tt_iterapp('mldivide',m2fun,m2type,m2fcnstr,pre_b2,min(eps_z, 0.5), max_rank, max_swp,varargin{:});    
0079     end;
0080 end;
0081 if (existsM3)
0082     [m3type,m3fun,m3fcnstr] = tt_iterchk(M3);
0083     if (compute_real_res==1)
0084     pre_b2 = tt_iterapp('mldivide',m3fun,m3type,m3fcnstr,pre_b2,min(eps_z, 0.5), max_rank, max_swp,varargin{:});    
0085     end;
0086 end;
0087 % pre_b = tt_stabilize(pre_b,0);
0088 
0089 norm_f = tt_dot2(pre_b,pre_b)*0.5;
0090 mod_norm_f = sign(norm_f)*mod(abs(norm_f), 10);
0091 order_norm_f = norm_f - mod_norm_f;
0092 cur_norm_f = exp(mod_norm_f);
0093 
0094 if (existsx0)
0095     x = x0;
0096 else
0097 %     x = pre_b;
0098 %     x = tt_scal2(b, log(1e-308), 1);
0099     x = tt_zeros(max(size(b)), tt_size(b));
0100 end;
0101 
0102 % x_old = x;
0103 
0104 H = zeros(maxin+1, maxin);
0105 v = cell(maxin, 1);
0106 
0107 if (nargout>1)
0108     RESVEC = ones(maxout, maxin);
0109 end;
0110 if (nargout>2)
0111     rw = zeros(maxout,maxin);
0112 end;
0113 if (nargout>2)
0114     rx = zeros(maxout,maxin);
0115 end;
0116 % z = cell(maxin, 1);
0117 
0118 err=2;
0119 max_err=Inf;
0120 old_err = 1;
0121 stagpoints=0;
0122 
0123 for nitout=1:maxout
0124     max_rank=[];
0125     err_for_trunc=1;    
0126     Ax = tt_iterapp('mtimes',afun,atype,afcnstr,x,min(eps_z/err_for_trunc, 0.5), max_rank, max_swp,varargin{:});
0127     r = tt_axpy2(0,1, pre_b, 0,-1, Ax, min(eps_z/err_for_trunc, 0.5), max_rank);
0128     if (existsM1)
0129         r = tt_iterapp('mldivide',m1fun,m1type,m1fcnstr,r,min(eps_z/err_for_trunc, 0.5), max_rank, max_swp,varargin{:});
0130 %         Ax = tt_stabilize(Ax,0);
0131         if (existsM2)
0132             r = tt_iterapp('mldivide',m2fun,m2type,m2fcnstr,r,min(eps_z/err_for_trunc, 0.5), max_rank, max_swp,varargin{:});
0133 %             Ax = tt_stabilize(Ax,0);
0134         end;
0135         if (existsM3)
0136             r = tt_iterapp('mldivide',m3fun,m3type,m3fcnstr,r,min(eps_z/err_for_trunc, 0.5), max_rank, max_swp,varargin{:});
0137 %             Ax = tt_stabilize(Ax,0);
0138         end;
0139     end;        
0140 %     Ax = tt_stabilize(Ax,0);
0141 %     r = tt_add(pre_b, tt_scal(Ax, -1));
0142 %     r = tt_compr2(r, min(eps_z/err, 0.5), max_rank);
0143 %     beta = sqrt(tt_dot(r, r))
0144     beta = tt_dot2(r,r)*0.5;
0145     mod_beta = sign(beta)*mod(abs(beta), 10);
0146     order_beta = beta - mod_beta;    
0147     cur_beta = exp(mod_beta);    
0148     if (verbose==1)
0149         real_beta = exp(beta);
0150         fprintf('   cur_beta = %g, real_beta = %g\n', cur_beta, real_beta);
0151     end;
0152 %     if (beta<tol) break; end;
0153     if (nitout==1) 
0154         cur_normb = cur_beta; 
0155         order_normb = order_beta;
0156 %         normb = beta;
0157     end;    
0158     v{1} = tt_scal2(r, -beta, 1);
0159 %     v{1}=tt_stabilize(v{1},0);
0160     for j=1:maxin
0161         max_w_rank = 0;
0162         w = tt_iterapp('mtimes',afun,atype,afcnstr,v{j},min(eps_z/err_for_trunc, 0.5), max_rank, max_swp,varargin{:});                
0163         max_w_rank = max([max_w_rank; tt_ranks(w)]);
0164 %         w=tt_stabilize(w,0);
0165 %         max_Mx_rank = max(tt_ranks(w))
0166         if (existsM1)
0167             w = tt_iterapp('mldivide',m1fun,m1type,m1fcnstr,w,min(eps_z/err_for_trunc, 0.5), max_rank, max_swp,varargin{:});
0168             max_w_rank = max([max_w_rank; tt_ranks(w)]);
0169 %             w=tt_stabilize(w,0);
0170             if (existsM2)
0171                 w = tt_iterapp('mldivide',m2fun,m2type,m2fcnstr,w,min(eps_z/err_for_trunc, 0.5), max_rank, max_swp,varargin{:});
0172                 max_w_rank = max([max_w_rank; tt_ranks(w)]);
0173 %                 w=tt_stabilize(w,0);
0174 %                 max_Px_rank = max(tt_ranks(w))
0175             end;
0176             if (existsM3)
0177                 w = tt_iterapp('mldivide',m3fun,m3type,m3fcnstr,w,min(eps_z/err_for_trunc, 0.5), max_rank, max_swp,varargin{:});
0178                 max_w_rank = max([max_w_rank; tt_ranks(w)]);
0179 %                 w=tt_stabilize(w,0);
0180             end;            
0181         end;
0182         
0183         max_wrank = max(tt_ranks(w));
0184                 
0185 %         wnew = w;
0186         for i=1:j
0187             H(i,j)=tt_dot(w, v{i});
0188 %             w = tt_axpy2(0,1, w, log(abs(H(i,j))+1e-308), -1*sign(H(i,j)), v{i}, min(eps_z/err_for_trunc, 0.5), max_rank);
0189             max_w_rank = max([max_w_rank; tt_ranks(w)]);
0190             w = tt_axpy2(0,1, w, log(abs(H(i,j))+1e-308), -1*sign(H(i,j)), v{i}, eps_z, max_rank);
0191             % Orthogonality muss sein!
0192 %             w = tt_axpy3(0,1, w, log(abs(H(i,j))+1e-308), -1*sign(H(i,j)), v{i});
0193         end;
0194 %         w = tt_wround([], w, min(eps_z/err_for_trunc,0.5), 'rmax', max_rank, 'nswp', max_swp, 'verb', 1);
0195 %         for i=1:j
0196 %             wv = abs(tt_dot(w, v{i}));
0197 %             if (wv>eps_z)
0198 %                 fprintf('Warning: dot(w,v{%d}) = %3.3e\n', i, wv);
0199 %             end;
0200 %         end;
0201 %         w=wnew;
0202         if (nargout>2)
0203             rw(nitout,j)=max_w_rank;
0204         end;
0205         
0206         H(j+1,j) = sqrt(tt_dot(w, w));
0207         if (j<maxin)
0208             v{j+1}=tt_scal2(w, -log(H(j+1,j)), 1);
0209 %             v{j+1}=tt_stabilize(v{j+1},0);
0210         end;        
0211         
0212         [UH,SH,VH]=svd(H(1:j+1, 1:j), 0);
0213         SH = diag(SH);
0214         sigma_min_H = min(SH);
0215         sigma_max_H = max(SH);
0216         if (verbose==1)       
0217             fprintf('   min(sigma(H)) = %g\n', sigma_min_H);
0218         end;
0219         SH(1:numel(find(SH>1e-100)))=1./SH(1:numel(find(SH>1e-100))); %pseudoinv(SH)
0220         SH = diag(SH);
0221         y = cur_beta*VH*SH*(UH(1,:)'); % y = pinv(H)*(beta*e1)
0222 %         y = beta*VH*SH*(UH(1,:)'); % y = pinv(H)*(beta*e1)
0223 
0224 %         err = norm(H(1:j+1, 1:j)*y-[beta zeros(1,j)]', 'fro')/normb;
0225 %         err = log(norm(H(1:j+1, 1:j)*y-[cur_beta zeros(1,j)]', 'fro')/cur_normb+1e-308);
0226 %         err = err+order_beta-order_normb;
0227         err = log(norm(H(1:j+1, 1:j)*y-[cur_beta zeros(1,j)]', 'fro')/cur_norm_f+1e-308);
0228         err = err+order_beta-order_norm_f;        
0229         err = exp(err);
0230         if (use_err_trunc==1)
0231 %             err_for_trunc=err;
0232             err_for_trunc = log(norm(H(1:j+1, 1:j)*y-[cur_beta zeros(1,j)]', 'fro')/cur_beta+1e-308);
0233             err_for_trunc = exp(err_for_trunc);
0234             err_for_trunc = (err_for_trunc*maxin*sigma_max_H/sigma_min_H)^(err_trunc_power);
0235             err_for_trunc = min(err_for_trunc, 1);
0236 %             err_for_trunc = err;
0237         end;
0238 %         err_to_f = log(norm(H(1:j+1, 1:j)*y-[cur_beta zeros(1,j)]', 'fro')/cur_norm_f+1e-308);
0239 %         err_to_f = err_to_f+order_beta-order_norm_f;
0240 %         err_to_f = exp(err_to_f);
0241         
0242         if (nargout>1)
0243             RESVEC(nitout,j)=err;
0244         end;
0245         
0246         
0247         x_new = x;
0248         max_x_rank = 0;
0249 %         dx = tt_scal2(v{j}, log(abs(y(j))+1e-308)+order_beta, sign(y(j)));
0250         for i=j:-1:1
0251 %             dx = tt_add(dx, tt_scal2(v{i}, log(abs(y(i))+1e-308)+order_beta, sign(y(i))));
0252 %             dx = tt_axpy2(0,1, dx, log(abs(y(i))+1e-308)+order_beta, sign(y(i)), v{i}, eps_x);
0253             x_new = tt_axpy2(0,1, x_new, log(abs(y(i))+1e-308)+order_beta, sign(y(i)), v{i}, eps_x);
0254 %             x_new = tt_axpy3(0,1, x_new, log(abs(y(i))+1e-308)+order_beta, sign(y(i)), v{i});
0255             max_x_rank = max([max_x_rank; tt_ranks(x_new)]);            
0256         end;
0257 %         x_new = tt_axpy2(0,1,x,0,1,dx,eps_x);
0258         if (nargout>3)
0259             rx(nitout,j)=max_x_rank;
0260         end;
0261 %         x_new = tt_wround([], x_new, eps_x, 'verb', 1);
0262 %         x_new = tt_add(x, dx);
0263 %         x_new = tt_compr2(x_new, eps_x);
0264 %         x_new = tt_axpy2(0,1,x,0,1,dx,eps_x);
0265        
0266         max_xrank = max(tt_ranks(x_new));
0267         max_rank = max_zrank_factor*max_xrank;
0268         
0269         if (compute_real_res==1)
0270             res = tt_iterapp('mtimes',afun,atype,afcnstr,x_new,eps_z, max_rank, max_swp,varargin{:});
0271             if (existsM1)
0272                 res = tt_iterapp('mldivide',m1fun,m1type,m1fcnstr,res,eps_z, max_rank, max_swp,varargin{:});
0273                 if (existsM2)
0274                     res = tt_iterapp('mldivide',m2fun,m2type,m2fcnstr,res,eps_z, max_rank, max_swp,varargin{:});
0275                 end;
0276                 if (existsM3)
0277                     res = tt_iterapp('mldivide',m3fun,m3type,m3fcnstr,res,eps_z, max_rank, max_swp,varargin{:});
0278                 end;
0279             end;
0280             res = tt_dist3(res, pre_b2)/sqrt(tt_dot(pre_b2,pre_b2));
0281         end;
0282         
0283         derr = old_err/err;
0284         old_err = err;
0285 %         dx = tt_add(x_new, tt_scal(x_old, -1));
0286 %         normx = sqrt(tt_dot(x_new,x_new));
0287 %         x_old=x_new;
0288 %         dx_norm = sqrt(tt_dot(dx, dx))/normx;
0289         if (derr<derr_tol_for_sp) stagpoints=stagpoints+1; end;       
0290         if (verbose==1)
0291 %             err_for_trunc
0292             if (compute_real_res==1)
0293                 fprintf('iter = [%d,%d], derr = %3.2f, resid=%3.2e, real_res=%3.2e, rank_w=%d, rank_x=%d, sp=%d, time=%g\n', nitout, j, derr, err, res, max_wrank, max_xrank, stagpoints, toc(t0));
0294             else                
0295                 fprintf('iter = [%d,%d], derr = %3.2f, resid=%3.2e, rank_w=%d, rank_x=%d, sp=%d, time=%g\n', nitout, j, derr, err, max_w_rank, max_x_rank, stagpoints, toc(t0));
0296             end;
0297         end;
0298         if (err<max_err)
0299             x_good = x_new;
0300             max_err=err;
0301         end;
0302         if (err<tol) break; end;
0303         if (stagpoints>=maxin*max_sp_factor) break; end;
0304     end;
0305     x = x_good;
0306 
0307     if (err<tol) break; end;
0308     if (stagpoints>=maxin*max_sp_factor) break; end;
0309 end;
0310 
0311 if (nargout>1)
0312     RESVEC=RESVEC(1:nitout,:);
0313     if (nitout==1)
0314         RESVEC=RESVEC(:,1:j);
0315     end;
0316 end;
0317 if (nargout>2)
0318     rw=rw(1:nitout,:);
0319     if (nitout==1)
0320         rw=rw(:,1:j);
0321     end;
0322 end;
0323 if (nargout>3)
0324     rx=rx(1:nitout,:);
0325     if (nitout==1)
0326         rx=rx(:,1:j);
0327     end;
0328 end;
0329 
0330 end

Generated on Wed 08-Feb-2012 18:20:24 by m2html © 2005