Home > tt2 > exp > tt_minres_selfprec.m

tt_minres_selfprec

PURPOSE ^

Computation of the approximate TT-matrix inverse using self-prec method

SYNOPSIS ^

function [X]=tt_minres_selfprec(A, tol, eps, maxit, prec_type)

DESCRIPTION ^

Computation of the approximate TT-matrix inverse using self-prec method
   X=TT_MINRES_SELFPREC(A,TOL,EPS,MAXIT,PREC_TYPE) Computes the
   approximate matrix inverse of the TT-matrix A (given in TT1.0 format)


 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 function [X]=tt_minres_selfprec(A, tol, eps, maxit, prec_type)
0002 %Computation of the approximate TT-matrix inverse using self-prec method
0003 %   X=TT_MINRES_SELFPREC(A,TOL,EPS,MAXIT,PREC_TYPE) Computes the
0004 %   approximate matrix inverse of the TT-matrix A (given in TT1.0 format)
0005 %
0006 %
0007 % TT-Toolbox 2.2, 2009-2012
0008 %
0009 %This is TT Toolbox, written by Ivan Oseledets et al.
0010 %Institute of Numerical Mathematics, Moscow, Russia
0011 %webpage: http://spring.inm.ras.ru/osel
0012 %
0013 %For all questions, bugs and suggestions please mail
0014 %ivan.oseledets@gmail.com
0015 %---------------------------
0016 
0017 matvec='mm+compr';
0018 max_rank=10;
0019 
0020 d=max(size(A));
0021 ns=tt_size(A);
0022 
0023 I=tt_eye(ns, d); I=core(I); 
0024 if (strcmp(prec_type, 'left'))
0025     X=tt_scal(tt_transp(A), 1e-15);
0026 else
0027     X=tt_scal(I, 1e-150);
0028 end;
0029 
0030 normf=sqrt(tt_dot(tt_mat_to_vec(I), tt_mat_to_vec(I)));
0031 err=1;
0032 err_old=2;
0033 sp=0;
0034 
0035 for iout=1:maxit
0036     if (strcmp(matvec, 'mmk2'))
0037         resid=tt_mmk2(A, X, eps, max_rank);
0038     else
0039         resid=tt_mm(A,X);
0040         resid=tt_mat_compr(resid, eps,max_rank*2);
0041     end;
0042     resid=ttm_add(I, tt_scal(resid, -1));
0043     
0044     
0045     if (iout>1)||(strcmp(prec_type, 'left'))
0046         if (strcmp(matvec, 'mmk2'))
0047             Xresid=tt_mmk2(X, resid, eps,max_rank);
0048         else
0049             Xresid=tt_mm(X,resid);
0050             Xresid=tt_mat_compr(Xresid, eps,max_rank*2);
0051         end;
0052     else
0053         Xresid=resid;
0054     end;
0055     
0056     max_xrrank=max(tt_ranks(tt_mat_to_vec(Xresid)));
0057     
0058     if (strcmp(prec_type, 'left')) % we use left preconditioner (igmres)
0059         if (strcmp(matvec, 'mmk2'))
0060             w=tt_mmk2(A, Xresid, eps,max_rank);
0061             w=tt_mmk2(X, w, eps,max_rank);
0062         else
0063             w=tt_mm(A,Xresid);
0064             w=tt_mat_compr(w, eps,max_rank);
0065             w=tt_mm(X,w);
0066 %             w=tt_mat_compr(w, eps);
0067         end;
0068         
0069         w = tt_mat_to_vec(w);        
0070         beta=sqrt(tt_dot(tt_mat_to_vec(Xresid),tt_mat_to_vec(Xresid)));
0071         wr = tt_dot(w, tt_mat_to_vec(Xresid));
0072         ww = tt_dot(w,w); 
0073         
0074         H=zeros(2,1);
0075         H(1,1)=wr/(beta^2);
0076         H(2,1)=sqrt(ww/(beta^2)-wr^2/(beta^4));
0077         rhs=H(1,1)*beta;
0078         
0079         y = rhs/(H'*H);
0080         err = err*norm(H*y-[beta; 0])/beta; 
0081         y=y/beta;
0082         
0083 %         y = tt_dot(w,w);
0084 %         if (y~=0)
0085 %             y = tt_dot(w, tt_mat_to_vec(Xresid))/y;
0086 %         else
0087 %             y=1;
0088 %         end;
0089     else % Use right preconditioner and fgmres
0090         if (strcmp(matvec, 'mmk2'))
0091             w=tt_mmk2(A, Xresid, eps,max_rank);
0092         else
0093             w=tt_mm(A,Xresid);
0094 %             w=tt_mat_compr(w, eps);
0095         end;
0096         
0097         w = tt_mat_to_vec(w);        
0098         
0099         beta=sqrt(tt_dot(tt_mat_to_vec(resid),tt_mat_to_vec(resid)));
0100         wr = tt_dot(w, tt_mat_to_vec(resid));
0101         ww = tt_dot(w,w); 
0102         
0103         H=zeros(2,1);
0104         H(1,1)=wr/(beta^2);
0105         H(2,1)=sqrt(ww/(beta^2)-wr^2/(beta^4));
0106         rhs=H(1,1)*beta;
0107         
0108         y = rhs/(H'*H);
0109         err = norm(H*y-[beta; 0])/normf; 
0110         y=y/beta;
0111         
0112 %         y2=ww;
0113 %         if (y~=0)
0114 %             y2 = wr/y2;
0115 %         else
0116 %             y2=1;
0117 %         end;
0118     end;
0119     
0120     max_wrank=max(tt_ranks(w));
0121                 
0122     Xresid=tt_mat_to_vec(Xresid);    
0123     
0124     X=tt_mat_to_vec(X);
0125     if (err<err_old)
0126         X = tt_axpy(1, X, y, Xresid, eps, max_rank);
0127         err_old=err;
0128     else
0129         sp=sp+1;
0130     end;
0131     max_xrank=max(tt_ranks(X));
0132     
0133     X = tt_vec_to_mat(X, ns, ns);
0134     
0135 %     if (strcmp(matvec, 'mmk2'))
0136 %         resid=tt_mmk2(A, X, eps, max_rank);
0137 %     else
0138 %         resid=tt_mm(A, X);
0139 % %         resid=tt_mat_compr(resid, eps);
0140 %     end;
0141 %     resid=ttm_add(I, tt_scal(resid, -1));
0142 %     resid=tt_mat_to_vec(resid);
0143 %
0144 %     err_real=sqrt(tt_dot(resid,resid))/normf;
0145     err_real=0;
0146     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);    
0147     if (err<tol)||(sp>0)
0148         break;
0149     end;
0150     
0151 end;
0152 
0153 end

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