Home > tt2 > @tt_matrix > mvk.m

mvk

PURPOSE ^

DMRG fast matrix-by-vector product, Method 1 (less accurate)

SYNOPSIS ^

function [y]=mvk(a,x,eps,nswp,z,rmax)

DESCRIPTION ^

DMRG fast matrix-by-vector product, Method 1 (less accurate)
   [Y]=MVK(A,X,EPS,NSWP,Y,RMAX) Matrix-by-vector product of a TT-matrix A
   by a TT-tensor X with accuracy EPS. Also, one can specify the number of
   sweeps NSWP, initial approximation Z and the maximal TT-rank RMAX (if
   they become too large)


 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 [y]=mvk(a,x,eps,nswp,z,rmax)
0002 %DMRG fast matrix-by-vector product, Method 1 (less accurate)
0003 %   [Y]=MVK(A,X,EPS,NSWP,Y,RMAX) Matrix-by-vector product of a TT-matrix A
0004 %   by a TT-tensor X with accuracy EPS. Also, one can specify the number of
0005 %   sweeps NSWP, initial approximation Z and the maximal TT-rank RMAX (if
0006 %   they become too large)
0007 %
0008 %
0009 % TT-Toolbox 2.2, 2009-2012
0010 %
0011 %This is TT Toolbox, written by Ivan Oseledets et al.
0012 %Institute of Numerical Mathematics, Moscow, Russia
0013 %webpage: http://spring.inm.ras.ru/osel
0014 %
0015 %For all questions, bugs and suggestions please mail
0016 %ivan.oseledets@gmail.com
0017 %---------------------------
0018 n=a.n;
0019 m=a.m;
0020 att=a.tt;
0021 corea=att.core;
0022 psa=att.ps;
0023 ra=att.r;
0024 d=att.d;
0025 if ( nargin <= 5 || isempty(rmax) )
0026    rmax=1000;
0027 end
0028 if ( nargin <= 4 || isempty(z) )
0029   %z=tt_random(n,ndims(x),5);
0030   %z=tt_tensor(z);
0031   z=x;
0032 end
0033 if ( nargin <= 3 || isempty(nswp) )
0034   nswp = 40;
0035 end
0036 y=z;
0037 
0038 %Warmup is to orthogonalize Y from right-to-left and compute psi-matrices
0039 %for Ax
0040 psi=cell(d+1,1); %Psi-matrices
0041 psi{d+1}=1; psi{1}=1;
0042 
0043 
0044 %Parameters section
0045 
0046 %Warmup: right-to-left QR + computation of psi matrices
0047 
0048   corex=x.core;
0049   psx=x.ps;
0050   rx=x.r;
0051 
0052     psy=y.ps;
0053   ry=y.r;
0054   corey=y.core;
0055    pos1=psy(d+1);
0056 cr1=corey(psy(d):psy(d+1)-1);
0057 for i=d:-1:2  
0058    cr2=corey(psy(i-1):psy(i)-1);
0059    cr1=reshape(cr1,[ry(i),n(i)*ry(i+1)]);
0060    cr2=reshape(cr2,[ry(i-1)*n(i-1),ry(i)]);
0061    cr1=cr1.';
0062    [q,rm]=qr(cr1,0); rn=size(q,2); rm=rm.';
0063    q=q.'; 
0064    ry(i)=rn;
0065    corey(pos1-ry(i+1)*n(i)*ry(i):pos1-1)=q(:);
0066    %Convolution is now performed for psi(i) using psi(i+1) and corea, and
0067    %(new) core q
0068    cra=corea(psa(i):psa(i+1)-1); cra=reshape(cra,[ra(i),n(i),m(i),ra(i+1)]);
0069    cry=reshape(conj(q),[ry(i),n(i),ry(i+1)]);
0070    crx=corex(psx(i):psx(i+1)-1); crx=reshape(crx,[rx(i),m(i),rx(i+1)]);
0071    pscur=psi{i+1}; pscur=reshape(pscur,[ra(i+1),rx(i+1),ry(i+1)]); %ra,rx,ry
0072    %First, convolve over rx(i+1)
0073    crx=reshape(crx,[rx(i)*m(i),rx(i+1)]);
0074    pscur=permute(pscur,[2,1,3]); pscur=reshape(pscur,[rx(i+1),ra(i+1)*ry(i+1)]);
0075    pscur=crx*pscur; %pscur is now rx(i)*m(i)*ra(i+1)*ry(i+1)
0076    %Convolve over m(i),ra(i+1),n(i),ry(i+1)
0077     pscur=reshape(pscur,[rx(i),m(i)*ra(i+1),ry(i+1)]);
0078     pscur=permute(pscur,[1,3,2]); 
0079     pscur=reshape(pscur,[rx(i)*ry(i+1),m(i)*ra(i+1)]);
0080     cra=reshape(cra,[ra(i)*n(i),m(i)*ra(i+1)]); cra=cra.';  
0081     pscur=pscur*cra; 
0082     %pscur is now rx(i)*ry(i+1)*ra(i)*n(i), it is left to convolve over
0083     %n(i)*ry(i+1)
0084     pscur=reshape(pscur,[rx(i),ry(i+1),ra(i),n(i)]);
0085     pscur=permute(pscur,[3,1,4,2]);
0086     pscur=reshape(pscur,[rx(i)*ra(i),n(i)*ry(i+1)]);
0087     cry=reshape(cry,[ry(i),n(i)*ry(i+1)]); cry=cry.';
0088     pscur=pscur*cry;
0089     psi{i}=pscur;
0090    %End of psi-block
0091    pos1=pos1-ry(i+1)*n(i)*ry(i);
0092    cr1=cr2*rm;
0093    
0094 end
0095 corey(pos1-ry(2)*n(1)*ry(1):pos1-1)=cr1(:);
0096 pos1=pos1-ry(2)*n(1)*ry(1);
0097 corey=corey(pos1:numel(corey)); %Truncate unused elements
0098 
0099   
0100 psy=cumsum([1;n.*ry(1:d).*ry(2:d+1)]);
0101 y.core=corey;
0102 y.r=ry;
0103 y.ps=psy;
0104 swp=1;
0105    converged=false;
0106 while (swp <= nswp && ~converged) 
0107   %left-to-right dmrg sweep
0108 psy=y.ps;
0109   ry=y.r;
0110   corey=y.core;
0111   
0112   pos1=1;
0113   cry_old=corey;
0114   converged=true;
0115   ermax=0;
0116   
0117   for i=1:d-1
0118      %We care for two cores, with number i & number i+1, and use
0119      %psi(i) and psi(i+2) as a basis; also we will need to recompute
0120      %psi(i+1)
0121      ps1=psi{i}; ps2=psi{i+2};
0122      cra1=corea(psa(i):psa(i+1)-1); cra2=corea(psa(i+1):psa(i+2)-1);
0123      crx1=corex(psx(i):psx(i+1)-1); crx2=corex(psx(i+1):psx(i+2)-1);
0124      %our convolution is
0125      %ps1(ra(i),rx(i),ry(i))*cra1(ra(i),n(i),m(i),ra(i+1))*
0126      %cra2(ra(i+1),n(i+1),m(i+1),ra(i+2))*
0127      %*ps2(ra(i+2),rx(i+2),ry(i+2))
0128      %*crx1(rx(i),m(i),rx(i+1))*cr2x(rx(i+1)*m(i+1)*rx(i+2))
0129      %Scheme ps1*crx1 over rx(i)
0130      
0131      ps1=reshape(ps1,[ra(i),rx(i),ry(i)]); 
0132      ps1=permute(ps1,[1,3,2]); ps1=reshape(ps1,[ra(i)*ry(i),rx(i)]);
0133      crx1=reshape(crx1,[rx(i),m(i)*rx(i+1)]);
0134      ps1=ps1*crx1; %ps1 is now ra(i)*ry(i)*m(i)*rx(i+1)
0135      %Now convolve with matrix A over ra(i)*m(i)
0136      ps1=reshape(ps1,[ra(i),ry(i),m(i),rx(i+1)]);
0137      ps1=permute(ps1,[2,4,1,3]);
0138      ps1=reshape(ps1,[ry(i)*rx(i+1),ra(i)*m(i)]);
0139      cra1=reshape(cra1,[ra(i),n(i),m(i),ra(i+1)]);
0140      cra1=permute(cra1,[1,3,2,4]); 
0141      cra1=reshape(cra1,[ra(i)*m(i),n(i)*ra(i+1)]);
0142      ps1=ps1*cra1; %ps1 is now ry(i)*rx(i+1)*n(i)*ra(i+1)
0143      %Then the ``same'' convolution is carried over for second pair of
0144      %cores
0145      ps2=reshape(ps2,[ra(i+2),rx(i+2),ry(i+2)]);
0146      ps2=permute(ps2,[2,1,3]);
0147      ps2=reshape(ps2,[rx(i+2),ra(i+2)*ry(i+2)]);
0148      crx2=reshape(crx2,[rx(i+1)*m(i+1),rx(i+2)]);
0149      ps2=crx2*ps2; %ps2 is now rx*(i+1)*m(i+1)*ra(i+2)*ry(i+2)
0150      %Convolve over m(i+1)*ra(i+2)
0151      ps2=reshape(ps2,[rx(i+1),m(i+1)*ra(i+2),ry(i+2)]);
0152      ps2=permute(ps2,[2,1,3]);
0153      ps2=reshape(ps2,[m(i+1)*ra(i+2),rx(i+1)*ry(i+2)]);
0154      cra2=reshape(cra2,[ra(i+1)*n(i+1),m(i+1)*ra(i+2)]);
0155      ps2=cra2*ps2; 
0156      %ps2 is now ra(i+1)*n(i+1)*rx(i+1)*ry(i+2)
0157      %Now form superblock by contraction ps1 & ps2
0158      %over ra(i+1)*rx(i+1)
0159      ps0=reshape(ps1,[ry(i),rx(i+1),n(i),ra(i+1)]);
0160      ps0=permute(ps0,[1,3,2,4]);
0161      ps0=reshape(ps0,[ry(i)*n(i),rx(i+1)*ra(i+1)]);
0162      ps2=reshape(ps2,[ra(i+1),n(i+1),rx(i+1),ry(i+2)]);
0163      ps2=permute(ps2,[3,1,2,4]); 
0164      ps2=reshape(ps2,[rx(i+1)*ra(i+1),n(i+1)*ry(i+2)]);
0165      super_core=ps0*ps2; %super_core is ry(i)*n(i)*n(i+1)*ry(i+2)
0166      %Compute previous supercore
0167      cr1=corey(pos1:pos1+ry(i)*n(i)*ry(i+1)-1);
0168      cr2=cry_old(psy(i+1):psy(i+2)-1); 
0169      cr1=reshape(cr1,[ry(i)*n(i),ry(i+1)]);
0170      cr2=reshape(cr2,[ry(i+1),n(i+1)*ry(i+2)]);
0171      super_core_old=cr1*cr2; 
0172      er=norm(super_core_old(:)-super_core(:))/norm(super_core(:));
0173      if ( er > eps) 
0174         converged=false;
0175      end
0176      ermax=max(er,ermax);
0177      %fprintf('i=%d er=%3.2e \n',i,er);
0178      [u,s,v]=svd(super_core,'econ');
0179      s=diag(s); 
0180      r=my_chop2(s,eps/sqrt(d-1)*norm(s)); r=min(r,rmax);
0181      u=u(:,1:r); s=s(1:r); v=v(:,1:r); v=v*diag(s);
0182      ry(i+1)=r;
0183      %u is ry(i)*n(i)*ry(i+1)
0184      %core_new(pos1:pos1+ry(i)*n(i)*ry(i+1)-1)=u(:);
0185      corey(pos1:pos1+ry(i)*n(i)*ry(i+1)-1)=u(:); 
0186      
0187      u=reshape(u,[ry(i),n(i),ry(i+1)]);
0188 
0189      %Compute new psi
0190      %ps1 is ry(i)*rx(i+1)*n(i)*ra(i+1) with u over ry(i)*n(i)
0191      u=conj(u);
0192      ps1=reshape(ps1,[ry(i),rx(i+1),n(i),ra(i+1)]);
0193      ps1=permute(ps1,[4,2,1,3]);
0194      ps1=reshape(ps1,[ra(i+1)*rx(i+1),ry(i)*n(i)]);
0195      u=reshape(u,[ry(i)*n(i),ry(i+1)]);
0196      ps1=ps1*u;
0197      psi{i+1}=ps1;
0198      %Compute (?) new v
0199      pos1=pos1+ry(i)*n(i)*ry(i+1);
0200      v=v';
0201      %v=reshape(v,[ry(i+1),n(i+1),ry(i+2)]);
0202      corey(pos1:pos1+ry(i+1)*n(i+1)*ry(i+2)-1)=v(:);
0203   end
0204   psy=cumsum([1;n.*ry(1:d).*ry(2:d+1)]);
0205 y.core=corey;
0206 y.r=ry;
0207 y.ps=psy;
0208 cry_old=corey;
0209 %m=numel(cry);
0210 %cry=[zeros(size(cry)),cry];
0211 %pos1=pos1+m;
0212 %cry2=cry_old(psy(i+1):psy(i+2)-1);
0213 corey=corey(psy(d):psy(d+1)-1); %Start--only two cores are left!
0214 
0215  
0216   for i=d-1:-1:1
0217      %We care for two cores, with number i & number i+1, and use
0218      %psi(i) and psi(i+2) as a basis; also we will need to recompute
0219      %psi(i)
0220      ps1=psi{i}; ps2=psi{i+2};
0221      cra1=corea(psa(i):psa(i+1)-1); cra2=corea(psa(i+1):psa(i+2)-1);
0222      crx1=corex(psx(i):psx(i+1)-1); crx2=corex(psx(i+1):psx(i+2)-1);
0223      %our convolution is
0224      %ps1(ra(i),rx(i),ry(i))*cra1(ra(i),n(i),m(i),ra(i+1))*
0225      %cra2(ra(i+1),n(i+1),m(i+1),ra(i+2))*
0226      %*ps2(ra(i+2),rx(i+2),ry(i+2))
0227      %*crx1(rx(i),m(i),rx(i+1))*cr2x(rx(i+1)*m(i+1)*rx(i+2))
0228      %Scheme ps1*crx1 over rx(i)
0229      
0230      ps1=reshape(ps1,[ra(i),rx(i),ry(i)]); 
0231      ps1=permute(ps1,[1,3,2]); ps1=reshape(ps1,[ra(i)*ry(i),rx(i)]);
0232      crx1=reshape(crx1,[rx(i),m(i)*rx(i+1)]);
0233      ps1=ps1*crx1; %ps1 is now ra(i)*ry(i)*m(i)*rx(i+1)
0234      %Now convolve with matrix A over ra(i)*m(i)
0235      ps1=reshape(ps1,[ra(i),ry(i),m(i),rx(i+1)]);
0236      ps1=permute(ps1,[2,4,1,3]);
0237      ps1=reshape(ps1,[ry(i)*rx(i+1),ra(i)*m(i)]);
0238      cra1=reshape(cra1,[ra(i),n(i),m(i),ra(i+1)]);
0239      cra1=permute(cra1,[1,3,2,4]); 
0240      cra1=reshape(cra1,[ra(i)*m(i),n(i)*ra(i+1)]);
0241      ps1=ps1*cra1; %ps1 is now ry(i)*rx(i+1)*n(i)*ra(i+1)
0242      %Then the ``same'' convolution is carried over for second pair of
0243      %cores
0244      ps2=reshape(ps2,[ra(i+2),rx(i+2),ry(i+2)]);
0245      ps2=permute(ps2,[2,1,3]);
0246      ps2=reshape(ps2,[rx(i+2),ra(i+2)*ry(i+2)]);
0247      crx2=reshape(crx2,[rx(i+1)*m(i+1),rx(i+2)]);
0248      ps2=crx2*ps2; %ps2 is now rx*(i+1)*m(i+1)*ra(i+2)*ry(i+2)
0249      %Convolve over m(i+1)*ra(i+2)
0250      ps2=reshape(ps2,[rx(i+1),m(i+1)*ra(i+2),ry(i+2)]);
0251      ps2=permute(ps2,[2,1,3]);
0252      ps2=reshape(ps2,[m(i+1)*ra(i+2),rx(i+1)*ry(i+2)]);
0253      cra2=reshape(cra2,[ra(i+1)*n(i+1),m(i+1)*ra(i+2)]);
0254      ps2=cra2*ps2; 
0255      %ps2 is now ra(i+1)*n(i+1)*rx(i+1)*ry(i+2)
0256      %Now form superblock by contraction ps1 & ps2
0257      %over ra(i+1)*rx(i+1)
0258      ps0=reshape(ps1,[ry(i),rx(i+1),n(i),ra(i+1)]);
0259      ps0=permute(ps0,[1,3,2,4]);
0260      ps0=reshape(ps0,[ry(i)*n(i),rx(i+1)*ra(i+1)]);
0261      ps2=reshape(ps2,[ra(i+1),n(i+1),rx(i+1),ry(i+2)]);
0262      ps2=permute(ps2,[3,1,2,4]); 
0263      ps2=reshape(ps2,[rx(i+1)*ra(i+1),n(i+1)*ry(i+2)]);
0264      super_core=ps0*ps2; %super_core is ry(i)*n(i)*n(i+1)*ry(i+2)
0265      %Compute previous supercore
0266      cry1=cry_old(psy(i):psy(i+1)-1);
0267      %cry2=cry(ry(:ry(i+1)*n(i+1)*ry(i+2));
0268      %cry1=cry(1:ry(i)*n(i)*ry(i+1));
0269      %cry2=cry(ry(i)*n(i)*ry(i+1)+1:ry(i)*n(i)*ry(i+1)+ry(i+1)*n(i+1)*ry(i+2));
0270      cry2=corey(1:ry(i+1)*n(i+1)*ry(i+2));
0271      corey(1:ry(i+1)*n(i+1)*ry(i+2))=[];
0272      cry1=reshape(cry1,[ry(i)*n(i),ry(i+1)]);
0273      cry2=reshape(cry2,[ry(i+1),n(i+1)*ry(i+2)]);
0274      super_core_old=cry1*cry2;
0275      er=norm(super_core_old(:)-super_core(:))/norm(super_core(:));
0276      if ( er > eps) 
0277         converged=false;
0278      end
0279      er=max(er,ermax);
0280      [u,s,v]=svd(super_core,'econ');
0281      s=diag(s); 
0282      r=my_chop2(s,eps/sqrt(d-1)*norm(s)); r=min(r,rmax);
0283      u=u(:,1:r); s=s(1:r); u=u*diag(s); v=v(:,1:r);
0284      ry(i+1)=r;
0285      v=v.';
0286        u1=u(:); u1=u1.'; v1=v(:); v1=v1.';
0287      corey=[u1,v1,corey];
0288       %Compute new psi
0289       %ps2 is rx(i+1)*ra(i+1)*n(i+1)*ry(i+2)
0290       %v is ry(i+1)*n(i+1)*ry(i+2)
0291       %convolve over n(i+1)*ry(i+2);
0292       v = v.';
0293       ps2=reshape(ps2,rx(i+1)*ra(i+1),n(i+1)*ry(i+2));
0294       ps2=ps2*v; %ps2 is now rx(i+1)*ra(i+1)*ry(i+1)
0295       ps2=reshape(ps2,[rx(i+1),ra(i+1),ry(i+1)]);
0296       ps2=permute(ps2,[2,1,3]);
0297       psi{i+1}=ps2;
0298   end  
0299   swp=swp+1;
0300   psy=cumsum([1;n.*ry(1:d).*ry(2:d+1)]);
0301   %fprintf('swp=%d er=%3.2e \n',swp,ermax);
0302   y.core=corey;
0303   y.r=ry;
0304   y.ps=psy;
0305 end
0306 %Compute alpha*y0 = Ax - y
0307 %p=p/tt_dot(y0,y0);
0308 %keyboard;%
0309 %end
0310 %p1=tt_mvdot(core(a),core(x),y0);
0311 %p2=dot(y,tt_tensor(y0));
0312 %abs(p1-p2)/abs(p1)
0313 %keyboard;
0314 
0315 
0316

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