Computation of the approximate TT-matrix inverse using self-prec method [X]=TT_MINRES_SELFPREC2(A,EPS,OPTIONS) Computation of the approximate TT-matrix inverse using Saad self-prec method. Options are provided in form 'PropertyName1',PropertyValue1,'PropertyName2',PropertyValue2 and so on. The parameters are set to default (in brackets in the following) The list of option names and default values are: o matvec - type of the local matvec [ {mm+compr} | mmk2 ] o max_rank - maximal TT-rank bound [1000] o prec_type - left or right inversion [ {left} | right ] o maxit - maximal number of iterations [10] o tol - the requested inversion tolerance [EPS] 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 ---------------------------
0001 function [X]=tt_minres_selfprec2(A, eps, varargin) 0002 %Computation of the approximate TT-matrix inverse using self-prec method 0003 % [X]=TT_MINRES_SELFPREC2(A,EPS,OPTIONS) Computation of the approximate 0004 % TT-matrix inverse using Saad self-prec method. Options are provided 0005 % in form 'PropertyName1',PropertyValue1,'PropertyName2',PropertyValue2 0006 % and so on. The parameters are set to default (in brackets in the 0007 % following) The list of option names and default values are: 0008 % o matvec - type of the local matvec [ {mm+compr} | mmk2 ] 0009 % o max_rank - maximal TT-rank bound [1000] 0010 % o prec_type - left or right inversion [ {left} | right ] 0011 % o maxit - maximal number of iterations [10] 0012 % o tol - the requested inversion tolerance [EPS] 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 % matvec='mmk2'; 0026 matvec='mm+compr'; 0027 max_rank=1000; 0028 prec_type='left'; 0029 tol=eps; 0030 maxit=10; 0031 0032 for i=1:2:length(varargin)-1 0033 switch lower(varargin{i}) 0034 case 'matvec' 0035 matvec=varargin{i+1}; 0036 case 'max_rank' 0037 max_rank=lower(varargin{i+1}); 0038 case 'prec_type' 0039 prec_type=varargin{i+1}; 0040 case 'maxit' 0041 maxit=varargin{i+1}; 0042 otherwise 0043 error('Unrecognized option: %s\n',varargin{i}); 0044 end 0045 end 0046 0047 0048 %d=max(size(A)); 0049 d=ndims(A); 0050 ns=A.n; 0051 %ns=tt_size(A); 0052 0053 I=tt_eye(ns, d); I=tt_matrix(I); 0054 if (strcmp(prec_type, 'left')) 0055 % X=tt_scal(tt_transp(A), 1e-15); 0056 X=(1e-15)*conj(A'); 0057 else 0058 % X=tt_scal(I, 1e-150); 0059 X=1e-150*I; 0060 end; 0061 normf=norm(I); 0062 %normf=sqrt(tt_dot(tt_mat_to_vec(I), tt_mat_to_vec(I))); 0063 err=1; 0064 err_old=2; 0065 sp=0; 0066 0067 for iout=1:maxit 0068 if (strcmp(matvec, 'mmk2')) 0069 resid=tt_mmk2(A, X, eps, max_rank); 0070 else 0071 %resid=tt_mm(A,X); 0072 resid=round(A*X,eps,max_rank*2); 0073 %resid=tt_mat_compr(resid, eps,max_rank*2); 0074 end; 0075 %resid=ttm_add(I, tt_scal(resid, -1)); 0076 resid=I-resid; 0077 0078 if (iout>1)||(strcmp(prec_type, 'left')) 0079 if (strcmp(matvec, 'mmk2')) 0080 Xresid=tt_mmk2(X, resid, eps,max_rank); 0081 else 0082 %Xresid=tt_mm(X,resid); 0083 Xresid=round(X*resid,eps,max_rank*2); 0084 %Xresid=tt_mat_compr(Xresid, eps,max_rank*2); 0085 0086 end; 0087 else 0088 Xresid=resid; 0089 end; 0090 0091 max_xrrank=max(rank(Xresid)); 0092 0093 if (strcmp(prec_type, 'left')) % we use left preconditioner (igmres) 0094 if (strcmp(matvec, 'mmk2')) 0095 w=tt_mmk2(A, Xresid, eps,max_rank); 0096 w=tt_mmk2(X, w, eps,max_rank); 0097 else 0098 w=round(A*Xresid,eps,max_rank); 0099 w=round(X*w,eps); 0100 %w=tt_mm(A,Xresid); 0101 %w=tt_mat_compr(w, eps,max_rank); 0102 %w=tt_mm(X,w); 0103 % w=tt_mat_compr(w, eps); 0104 end; 0105 %w = tt_mat_to_vec(w); 0106 %beta=sqrt(tt_dot(tt_mat_to_vec(Xresid),tt_mat_to_vec(Xresid))); 0107 beta=norm(Xresid); 0108 wr=dot(w,Xresid); 0109 %wr = tt_dot(w, tt_mat_to_vec(Xresid)); 0110 ww=dot(w,w); 0111 0112 H=zeros(2,1); 0113 H(1,1)=wr/(beta^2); 0114 H(2,1)=sqrt(ww/(beta^2)-wr^2/(beta^4)); 0115 rhs=H(1,1)*beta; 0116 0117 y = rhs/(H'*H); 0118 err = err*norm(H*y-[beta; 0])/beta; 0119 y=y/beta; 0120 0121 % y = tt_dot(w,w); 0122 % if (y~=0) 0123 % y = tt_dot(w, tt_mat_to_vec(Xresid))/y; 0124 % else 0125 % y=1; 0126 % end; 0127 else % Use right preconditioner and fgmres 0128 if (strcmp(matvec, 'mmk2')) 0129 w=tt_mmk2(A, Xresid, eps,max_rank); 0130 else 0131 w=A*Xresid; 0132 w=round(w,eps); 0133 % w=tt_mat_compr(w, eps); 0134 end; 0135 w=tt_tensor(w); 0136 %w = tt_mat_to_vec(w); 0137 beta=norm(resid); 0138 %beta=sqrt(tt_dot(tt_mat_to_vec(resid),tt_mat_to_vec(resid))); 0139 wr=dot(w,tt_tensor(resid)); 0140 ww=dot(w,w); 0141 %wr = tt_dot(w, tt_mat_to_vec(resid)); 0142 %ww = tt_dot(w,w); 0143 0144 H=zeros(2,1); 0145 H(1,1)=wr/(beta^2); 0146 H(2,1)=sqrt(ww/(beta^2)-wr^2/(beta^4)); 0147 rhs=H(1,1)*beta; 0148 0149 y = rhs/(H'*H); 0150 err = norm(H*y-[beta; 0])/normf; 0151 y=y/beta; 0152 0153 % y2=ww; 0154 % if (y~=0) 0155 % y2 = wr/y2; 0156 % else 0157 % y2=1; 0158 % end; 0159 end; 0160 0161 max_wrank=max(rank(w)); 0162 0163 %Xresid=tt_mat_to_vec(Xresid); 0164 Xresid=tt_tensor(Xresid); 0165 X=tt_tensor(X); 0166 %X=tt_mat_to_vec(X); 0167 %if (err<err_old) 0168 %X = tt_axpy(1, X, y, Xresid, eps, max_rank); 0169 X=round(X+y*Xresid,eps,max_rank); 0170 0171 %err_old=err; 0172 %else 0173 % sp=sp+1; 0174 %end; 0175 %max_xrank=max(tt_ranks(X)); 0176 max_xrank=max(rank(X)); 0177 %X = tt_vec_to_mat(X, ns, ns); 0178 X=tt_matrix(X,ns,ns); 0179 % if (strcmp(matvec, 'mmk2')) 0180 % resid=tt_mmk2(A, X, eps, max_rank); 0181 % else 0182 % resid=tt_mm(A, X); 0183 % % resid=tt_mat_compr(resid, eps); 0184 % end; 0185 % resid=ttm_add(I, tt_scal(resid, -1)); 0186 % resid=tt_mat_to_vec(resid); 0187 % 0188 % err_real=sqrt(tt_dot(resid,resid))/normf; 0189 err_real=0; 0190 fprintf('iter [%d], res_real=%3.3e, x_rank=%d, Xr_rank=%d, w_rank=%d, err=%3.3e\n', iout, err_real, max_xrank, max_xrrank, max_wrank, err); 0191 if (err<tol)||(sp>0) 0192 break; 0193 end 0194 0195 end 0196 %keyboard;