Home > tt2 > solve > dmrg_rake_solve2.m

dmrg_rake_solve2

PURPOSE ^

DMRG-type method for the solution of linear systems in QTT-Tucker format

SYNOPSIS ^

function [x]=dmrg_rake_solve2(A, y, tol, varargin)

DESCRIPTION ^

DMRG-type method for the solution of linear systems in QTT-Tucker format
   [X]=DMRG_RAKE_SOLVE2(A,Y,TOL,VARARGIN) Attempts to solve the linear
   system A*X = Y with accuracy EPS using the two-sided DMRG iteration.
   Matrix A has to be given in the QTT-Tucker, right-hand side Y should be
   given in the QTT-Tucker format also. 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 x0 - initial approximation [random rank-2 tensor] 
       o nswp - maximal number of DMRG sweeps [10]
       o rmax - maximal TT-rank of the solution [1000]
       o verb - verbosity level, 0-silent, 1-sweep info, 2-block info [1]
       o kick_rank - stabilization parameter [2]
       o max_full_size - maximal size of the local matrix to full solver 
       [2500]
       o local_prec: Local preconditioner, 'als' - ALS-Richardson
       iteration, 'selfprec' (Saad selfpreconditioner) ['als']
       o gmres_iters - number of local gmres restarts [2]
       o nrestart - dimension of local gmres [25]


 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:

SUBFUNCTIONS ^

SOURCE CODE ^

0001 function [x]=dmrg_rake_solve2(A, y, tol, varargin)
0002 %DMRG-type method for the solution of linear systems in QTT-Tucker format
0003 %   [X]=DMRG_RAKE_SOLVE2(A,Y,TOL,VARARGIN) Attempts to solve the linear
0004 %   system A*X = Y with accuracy EPS using the two-sided DMRG iteration.
0005 %   Matrix A has to be given in the QTT-Tucker, right-hand side Y should be
0006 %   given in the QTT-Tucker format also. Options are provided in form
0007 %   'PropertyName1',PropertyValue1,'PropertyName2',PropertyValue2 and so
0008 %   on. The parameters are set to default (in brackets in the following)
0009 %   The list of option names and default values are:
0010 %       o x0 - initial approximation [random rank-2 tensor]
0011 %       o nswp - maximal number of DMRG sweeps [10]
0012 %       o rmax - maximal TT-rank of the solution [1000]
0013 %       o verb - verbosity level, 0-silent, 1-sweep info, 2-block info [1]
0014 %       o kick_rank - stabilization parameter [2]
0015 %       o max_full_size - maximal size of the local matrix to full solver
0016 %       [2500]
0017 %       o local_prec: Local preconditioner, 'als' - ALS-Richardson
0018 %       iteration, 'selfprec' (Saad selfpreconditioner) ['als']
0019 %       o gmres_iters - number of local gmres restarts [2]
0020 %       o nrestart - dimension of local gmres [25]
0021 %
0022 %
0023 % TT-Toolbox 2.2, 2009-2012
0024 %
0025 %This is TT Toolbox, written by Ivan Oseledets et al.
0026 %Institute of Numerical Mathematics, Moscow, Russia
0027 %webpage: http://spring.inm.ras.ru/osel
0028 %
0029 %For all questions, bugs and suggestions please mail
0030 %ivan.oseledets@gmail.com
0031 %---------------------------
0032 
0033 nswp = 20;
0034 local_format = 'full';
0035 % local_format = 'tt';
0036 max_full_size = 1500;
0037 max_full_size2 = Inf;
0038 nrestart = 25;
0039 gmres_iters = 2;
0040 verb = 1;
0041 kickrank = 2;
0042 checkrank = 1;
0043 resid_damp_glob = 1.1;
0044 resid_damp_loc = 1.1;
0045 rmax = Inf;
0046 
0047 x = [];
0048 
0049 for i=1:2:length(varargin)-1
0050     switch lower(varargin{i})
0051         case 'nswp'
0052             nswp=varargin{i+1};
0053         case 'rmax'
0054            rmax=lower(varargin{i+1});
0055         case 'x0'
0056             x=varargin{i+1};
0057         case 'verb'
0058             verb=varargin{i+1};
0059 %         case 'local_prec'
0060 %             local_prec=varargin{i+1};
0061         case 'nrestart'
0062             nrestart=varargin{i+1};
0063         case 'gmres_iters'
0064             gmres_iters=varargin{i+1};
0065         case 'kickrank'
0066             kickrank=varargin{i+1};
0067         case  'max_full_size'
0068             max_full_size=varargin{i+1};
0069         case  'resid_damp'
0070             resid_damp_loc=varargin{i+1};            
0071 %         case 'prec_compr'
0072 %             prec_compr=varargin{i+1};
0073 %         case 'prec_tol'
0074 %             prec_tol=varargin{i+1};
0075 %         case 'prec_iters'
0076 %             prec_iters=varargin{i+1};
0077 %         case 'use_self_prec'
0078 %             use_self_prec=varargin{i+1};
0079 %         case 'ddpow'
0080 %             ddpow=varargin{i+1};
0081 %         case 'ddrank'
0082 %             ddrank=varargin{i+1};
0083 %         case 'd_pow_check'
0084 %             d_pow_check=varargin{i+1};
0085 %         case 'bot_conv'
0086 %             bot_conv=varargin{i+1};
0087 %         case 'top_conv'
0088 %             top_conv=varargin{i+1};
0089 %         case 'min_dpow'
0090 %             min_dpow=varargin{i+1};
0091 %         case 'min_drank'
0092 %             min_drank=varargin{i+1};
0093         otherwise
0094             error('Unrecognized option: %s\n',varargin{i});
0095     end
0096 end
0097 
0098 d = y.dphys; % Physical dim.
0099 yc = y.core;
0100 yf = y.tuck;
0101 Af = A.tuck;
0102 Ac = A.core;
0103 
0104 L = zeros(1,d); % Quantics dims
0105 n = zeros(max(L), d); % Physical mode sizes
0106 for i=1:d
0107     L(i) = yf{i}.d;
0108     n(1:L(i), i) = yf{i}.n;
0109 end;
0110 
0111 if (isempty(x))
0112     xc = tt_rand(2,d,2);
0113     xf = cell(d,1);
0114     for i=1:d
0115         xf{i} = tt_rand(n(1:L(i),i), L(i), [1;2*ones(L(i),1)]);
0116     end;
0117 else
0118     xc = x.core;
0119     xf = x.tuck;
0120 end;
0121 
0122 
0123 % Extract ranks; Note that rf(L(i)+1,i) = r_tuck(i)
0124 rcy = yc.r;
0125 rfy = zeros(max(L)+1, d);
0126 for i=1:d
0127     rfy(1:L(i)+1, i) = yf{i}.r;
0128 end;
0129 rcA = Ac.r;
0130 rfA = zeros(max(L)+1, d);
0131 for i=1:d
0132     rfA(1:L(i)+1, i) = Af{i}.r;
0133 end;
0134 rcx = xc.r;
0135 rfx = zeros(max(L)+1, d);
0136 for i=1:d
0137     rfx(1:L(i)+1, i) = xf{i}.r;
0138 end;
0139 
0140 % Init phis. Thousands of them... (c)
0141 phcA = cell(d+1,1); phcA{1} = 1; phcA{d+1}=1; % Core, matrix
0142 phcy = cell(d+1,1); phcy{1} = 1; phcy{d+1}=1; % Core, rhs
0143 phfA = cell(d,1); % factors, matrix
0144 phfy = cell(d,1); % factors, rhs
0145 phAfc = cell(d,1); % factor-core, matrix
0146 phyfc = cell(d,1); % factor-core, rhs
0147 for i=1:d
0148     phfA{i} = cell(L(i)+1,1);
0149     phfA{i}{1} = 1; phfA{i}{L(i)+1} = 1;
0150     phfy{i} = cell(L(i)+1,1);
0151     phfy{i}{1} = 1; phfy{i}{L(i)+1} = 1;
0152 end;
0153 
0154 % For random check
0155 cphcA = cell(d+1,1); cphcA{1} = 1; cphcA{d+1}=1; % Core, matrix
0156 cphcy = cell(d+1,1); cphcy{1} = 1; cphcy{d+1}=1; % Core, rhs
0157 cphfA = cell(d,1); % factors, matrix
0158 cphfy = cell(d,1); % factors, rhs
0159 cphAfc = cell(d,1); % factor-core, matrix
0160 cphyfc = cell(d,1); % factor-core, rhs
0161 for i=1:d
0162     cphfA{i} = cell(L(i)+1,1);
0163     cphfA{i}{1} = 1; cphfA{i}{L(i)+1} = 1;
0164     cphfy{i} = cell(L(i)+1,1);
0165     cphfy{i}{1} = 1; cphfy{i}{L(i)+1} = 1;
0166 end;
0167 
0168 
0169 last_sweep = false;
0170 
0171 for swp=1:nswp
0172     % init check vector
0173 %      chkc = tt_rand(kickrank, d, kickrank);
0174     rcchk = [1; checkrank*ones(d-1,1); 1];
0175 %      chkf = cell(d,1);
0176     rfchk = zeros(max(L)+1, d);
0177     for i=1:d
0178 %          chkf{i} = tt_rand(n(1:L(i),i), L(i), [1;kickrank*ones(L(i),1)]);
0179         rfchk(1:L(i)+1, i) = [1; checkrank*ones(L(i),1)];
0180     end;
0181 
0182     dx_max = 0;
0183     res_max = 0;
0184     r_max = 0;
0185     chk_res_max = 0;
0186     % bottom-to-top QR and phis
0187     for i=d:-1:1 % physical dims/core
0188         for j=1:L(i) % quantics dims
0189             cr = xf{i}{j};
0190             cr = reshape(cr, rfx(j,i)*n(j,i), rfx(j+1,i));
0191             [cr, rv] = qr(cr, 0);
0192             % What is our next core?
0193             if (j<L(i))
0194                 % we are still on a "tooth"
0195                 cr2 = xf{i}{j+1};
0196                 cr2 = reshape(cr2, rfx(j+1,i), n(j+1,i)*rfx(j+2,i));
0197                 cr2 = rv*cr2;
0198                 rfx(j+1,i) = size(cr, 2);
0199                 xf{i}{j} = reshape(cr, rfx(j,i), n(j,i), rfx(j+1,i));
0200                 xf{i}{j+1} = reshape(cr2, rfx(j+1,i), n(j+1,i), rfx(j+2,i));
0201             else
0202                 % We have to convlove rv to the tucker core
0203                 cr2 = xc{i};
0204                 cr2 = permute(cr2, [2, 1, 3]);
0205                 cr2 = reshape(cr2, rfx(j+1,i), rcx(i)*rcx(i+1));
0206                 cr2 = rv*cr2;
0207                 rfx(j+1,i) = size(cr, 2);
0208                 cr2 = reshape(cr2, rfx(j+1,i), rcx(i), rcx(i+1));
0209                 cr2 = permute(cr2, [2, 1, 3]);
0210                 xf{i}{j} = reshape(cr, rfx(j,i), n(j,i), rfx(j+1,i));
0211                 xc{i} = cr2;
0212             end;
0213             % Update bottom phis
0214             cr = reshape(cr, rfx(j,i), n(j,i), rfx(j+1,i));
0215             if (j<L(i))
0216                 phfA{i}{j+1} = compute_next_Phi(phfA{i}{j}, cr, Af{i}{j}, cr, 'lr');
0217                 phfy{i}{j+1} = compute_next_Phi(phfy{i}{j}, cr, [], yf{i}{j}, 'lr');
0218             else
0219                 phAfc{i} = compute_next_Phi(phfA{i}{j}, cr, Af{i}{j}, cr, 'lr');
0220                 phyfc{i} = compute_next_Phi(phfy{i}{j}, cr, [], yf{i}{j}, 'lr');
0221             end;
0222 
0223             % check vector
0224 %              ccr = ones(rfchk(j,i)*n(j,i), checkrank);
0225 %  %              [ccr, rv] = qr(ccr, 0);
0226 %              rfchk(j+1,i) = size(ccr, 2);
0227 %              % Update bottom phis
0228 %              ccr = reshape(ccr, rfchk(j,i), n(j,i), rfchk(j+1,i));
0229 %              chkf{i}{j} = ccr;
0230             ccr = ones(1, n(j,i), 1);
0231             if (j<L(i))
0232                 cphfA{i}{j+1} = compute_next_Phi(cphfA{i}{j}, ccr, Af{i}{j}, cr, 'lr');
0233                 cphfy{i}{j+1} = compute_next_Phi(cphfy{i}{j}, ccr, [], yf{i}{j}, 'lr');
0234             else
0235                 cphAfc{i} = compute_next_Phi(cphfA{i}{j}, ccr, Af{i}{j}, cr, 'lr');
0236                 cphyfc{i} = compute_next_Phi(cphfy{i}{j}, ccr, [], yf{i}{j}, 'lr');
0237             end;
0238         end;
0239     end;
0240 
0241     % QRs and phis over the core
0242     % Project the system on the factors
0243     Acr = tt_matrix(Ac, Ac.n, ones(d,1));
0244     cAcr = tt_matrix(Ac, Ac.n, ones(d,1));
0245     ycr = yc;
0246     cycr = yc;
0247     for i=1:d
0248         Acr{i} = core_matrix(Ac{i}, phAfc{i});
0249         ycr{i} = core_vector(yc{i}, phyfc{i});
0250         cAcr{i} = core_matrix(Ac{i}, cphAfc{i});
0251         cycr{i} = core_vector(yc{i}, cphyfc{i});
0252     end;
0253     for i=d:-1:2
0254         rtx = rfx(L(i)+1, i);
0255         cr = xc{i};
0256         cr = reshape(cr, rcx(i), rtx*rcx(i+1));
0257         [cr, rv] = qr(cr.', 0);
0258         cr2 = xc{i-1};
0259         rtx2 = rfx(L(i-1)+1, i-1);
0260         cr2 =  reshape(cr2, rcx(i-1)*rtx2, rcx(i));
0261         cr2 = cr2*(rv.');
0262         rcx(i) = size(cr, 2);
0263         cr = reshape(cr.', rcx(i), rtx, rcx(i+1));
0264         xc{i-1} = reshape(cr2, rcx(i-1), rtx2, rcx(i));
0265         xc{i} = cr;
0266         % Update right phi
0267         phcA{i} = compute_next_Phi(phcA{i+1}, cr, Acr{i}, cr, 'rl');
0268         phcy{i} = compute_next_Phi(phcy{i+1}, cr, [], ycr{i}, 'rl');
0269 
0270         % check vector
0271 %          rtchk = rfchk(L(i)+1, i);
0272 %          ccr = ones(rtchk*rcchk(i+1), checkrank);
0273 %          [ccr, rv] = qr(ccr.', 0);
0274 %          rcchk(i) = size(ccr, 2);
0275 %          ccr = reshape(ccr.', rcchk(i), rtchk, rcchk(i+1));
0276         ccr = 1;
0277 %          chkc{i} = ccr;
0278         % Update right phi
0279         cphcA{i} = compute_next_Phi(cphcA{i+1}, ccr, cAcr{i}, cr, 'rl');
0280         cphcy{i} = compute_next_Phi(cphcy{i+1}, ccr, [], cycr{i}, 'rl');
0281     end;
0282 
0283 %     % DMRG over the core
0284 %     % Compute the reduced matrix and rhs;
0285 %     rta = zeros(d,1);
0286 %     for i=1:d
0287 %         rta(i) = rfA(L(i)+1,i);
0288 %     end;
0289 %     Acr = tt_matrix(Ac, rta, ones(d,1));
0290 %     ycr = yc;
0291 %     for i=1:d
0292 %         a1 = Ac{i};
0293 %         a1 = permute(a1, [2, 1, 3]);
0294 %         a1 = reshape(a1, rfA(L(i)+1,i), rcA(i)*rcA(i+1));
0295 %         ph2 = phfAb{i};
0296 %         ph2 = permute(ph2, [1,3,2]);
0297 %         ph2 = reshape(ph2, rfx(L(i)+1,i)*rfx(L(i)+1,i), rfA(L(i)+1,i));
0298 %         a1 = ph2*a1;
0299 %         a1 = reshape(a1, rfx(L(i)+1,i), rfx(L(i)+1,i), rcA(i), rcA(i+1));
0300 %         a1 = permute(a1, [3, 1, 2, 4]);
0301 %         Acr{i} = a1;
0302 %
0303 %         y1 = yc{i};
0304 %         y1 = permute(y1, [2, 1, 3]);
0305 %         y1 = reshape(y1, rfy(L(i)+1,i), rcy(i)*rcy(i+1));
0306 %         ph2 = phfyb{i};
0307 %         y1 = ph2*y1;
0308 %         y1 = reshape(y1, rfx(L(i)+1,i), rcy(i), rcy(i+1));
0309 %         y1 = permute(y1, [2, 1, 3]);
0310 %         ycr{i} = y1;
0311 %     end;
0312 %     fprintf('=rake_solve========= core processing, sweep %d =======\n', swp);
0313 %     xcold = xc;
0314 %     xc = dmrg_solve2(Acr, ycr, tol, 'x0', xcold, 'max_full_size', max_full_size, 'nrestart', nrestart, 'gmres_iters', gmres_iters, 'nswp', 1);
0315 %     res = norm(Acr*xcold-ycr)/norm(ycr);
0316 %     dx = norm(xcold-xc)/norm(xc);
0317 %     dx_max = max(dx_max, dx);
0318 %     r_max = max([r_max; rank(xc)]);
0319 %     fprintf('=rake_solve========= core res_prev: %3.3e, dx: %3.3e, rmax: %d\n', res, dx, max(rank(xc)));
0320 %     res_max = max(res_max, res);
0321 
0322 %     % Horisontal phis
0323 %     rcx = xc.r;
0324 %     for i=d:-1:2
0325 %         % n here is in fact a tucker rank
0326 %         cr2 = xc{i};
0327 %         cr2 = reshape(cr2, rcx(i), rfx(L(i)+1,i)*rcx(i+1));
0328 %         [cr2, rv]=qr(cr2.', 0);
0329 %         cr3 = xc{i-1};
0330 %         cr3 = reshape(cr3, rcx(i-1)*rfx(L(i-1)+1,i-1), rcx(i));
0331 %         cr3 = cr3*(rv.');
0332 %         rcx(i) = size(cr2, 2);
0333 %         cr2 = cr2.'; % sizes rcx1, rtuck*rcx2
0334 %         xc{i} = reshape(cr2, rcx(i), rfx(L(i)+1,i), rcx(i+1));
0335 %         xc{i-1} = reshape(cr3, rcx(i-1), rfx(L(i-1)+1,i-1), rcx(i));
0336 %
0337 %         % Update right phis
0338 %         cr2 = reshape(cr2, rcx(i), rfx(L(i)+1,i), rcx(i+1));
0339 %         phcAr{i} = compute_next_Phi(phcAr{i+1}, cr2, Acr{i}, cr2, 'rl');
0340 %         % New size: rcx1, rcA1, rcx1
0341 %         phcyr{i} = compute_next_Phi(phcyr{i+1}, cr2, [], ycr{i}, 'rl');
0342 %         % New size: rcx1, rcy1
0343 %     end;
0344 
0345     % DMRG over the factors
0346     for i=1:d
0347         % Convolve the tucker block to the last physical
0348         % Preparing the matrix and rhs in this case smells like
0349         % manure.
0350         % phcAl - Ac{i} - phcAr     <- this has to be convolved
0351         %          |                <- this has to be convolved
0352         %         Af{i}{j}
0353         %          |
0354         %         Af{i}{j-1}
0355         %          |                <- this has to be convolved
0356         %         phfAb{j-1}
0357         rta = rfA(L(i)+1,i); rtx = rfx(L(i)+1,i); rty = rfy(L(i)+1,i); rtchk = rfchk(L(i)+1,i);
0358         a1left = phcA{i};
0359         a1left = permute(a1left, [1,3,2]);
0360         a1left = reshape(a1left, rcx(i)*rcx(i), rcA(i));
0361         a1left = a1left*reshape(Ac{i}, rcA(i), rta*rcA(i+1));
0362 
0363         ca1left = cphcA{i};
0364         ca1left = permute(ca1left, [1,3,2]);
0365         ca1left = reshape(ca1left, rcchk(i)*rcx(i), rcA(i));
0366         ca1left = ca1left*reshape(Ac{i}, rcA(i), rta*rcA(i+1));
0367         % old
0368 %         % a1left is to use later to compute phcAl
0369 %         ph2 = phcAr{i+1};
0370 %         ph2 = permute(ph2, [2, 1, 3]);
0371 %         ph2 = reshape(ph2, rcA(i+1), rcx(i+1)*rcx(i+1));
0372 %         a1top = reshape(a1left, rcx(i)*rcx(i)*rta, rcA(i+1))*ph2;
0373 %         a1top = reshape(a1top, rcx(i), rcx(i), rta, rcx(i+1), rcx(i+1));
0374 %         a1top = permute(a1top, [1, 4, 2, 5, 3]);
0375 %         % we need a1 as well, to check the  residual in tucker block
0376 %         % splitting
0377 %         a1top = reshape(a1top, rcx(i)*rcx(i+1)*rcx(i)*rcx(i+1), rta);
0378         % new
0379         a1left = reshape(a1left, rcx(i)*rcx(i), rta, rcA(i+1));
0380         a1left = permute(a1left, [2, 1, 3]);
0381         a1left = reshape(a1left, rta, rcx(i)*rcx(i)*rcA(i+1));
0382 
0383         ca1left = reshape(ca1left, rcchk(i)*rcx(i), rta, rcA(i+1));
0384         ca1left = permute(ca1left, [2, 1, 3]);
0385         ca1left = reshape(ca1left, rta, rcchk(i)*rcx(i)*rcA(i+1));
0386 
0387 
0388         a2 = Af{i}{L(i)};
0389         a2 = reshape(a2, rfA(L(i),i)*n(L(i),i)*n(L(i),i), rta);
0390         % old
0391 %         a2 = a2*a1top.';
0392 %         a2 = reshape(a2, rfA(L(i),i), n(L(i),i), n(L(i),i), rcx(i)*rcx(i+1), rcx(i)*rcx(i+1));
0393 %         a2 = permute(a2, [1, 2, 4, 3, 5]);
0394 %         a2 = reshape(a2, rfA(L(i),i), n(L(i),i)*rcx(i)*rcx(i+1), n(L(i),i)*rcx(i)*rcx(i+1), 1);
0395         % new
0396         a2 = a2*a1left;
0397         a2 = reshape(a2, rfA(L(i),i), n(L(i),i), n(L(i),i), rcx(i), rcx(i), rcA(i+1));
0398         a2 = permute(a2, [1, 2, 4, 3, 5, 6]);
0399         a2 = reshape(a2, rfA(L(i),i), n(L(i),i)*rcx(i), n(L(i),i)*rcx(i), rcA(i+1));
0400 
0401         ca2 = Af{i}{L(i)};
0402         ca2 = reshape(ca2, rfA(L(i),i)*n(L(i),i)*n(L(i),i), rta);
0403         ca2 = ca2*ca1left;
0404         ca2 = reshape(ca2, rfA(L(i),i), n(L(i),i), n(L(i),i), rcchk(i), rcx(i), rcA(i+1));
0405         ca2 = permute(ca2, [1, 2, 4, 3, 5, 6]);
0406         ca2 = reshape(ca2, rfA(L(i),i), n(L(i),i)*rcchk(i), n(L(i),i)*rcx(i), rcA(i+1));
0407 
0408         Afr = Af{i};
0409         Afr{L(i)} = a2;
0410 
0411         cAfr = Af{i};
0412         cAfr{L(i)} = ca2;
0413 
0414         y1left = phcy{i};
0415         y1left = y1left*reshape(yc{i}, rcy(i), rty*rcy(i+1));
0416         % y1left is to use later to compute phcyl
0417         ph2 = phcy{i+1};
0418         ph2 = ph2.';
0419         y1top = reshape(y1left, rcx(i)*rty, rcy(i+1))*ph2;
0420         y1top = reshape(y1top, rcx(i), rty, rcx(i+1));
0421         y1top = permute(y1top, [1, 3, 2]);
0422         y1top = reshape(y1top, rcx(i)*rcx(i+1), rty);
0423         y2 = yf{i}{L(i)};
0424         y2 = reshape(y2, rfy(L(i),i)*n(L(i),i), rty);
0425         y2 = y2*y1top.';
0426         y2 = reshape(y2, rfy(L(i),i), n(L(i),i)*rcx(i)*rcx(i+1), 1);
0427         yfr = yf{i};
0428         yfr{L(i)} = y2;
0429 
0430         cy1left = cphcy{i};
0431         cy1left = cy1left*reshape(yc{i}, rcy(i), rty*rcy(i+1));
0432         % y1left is to use later to compute phcyl
0433         ph2 = cphcy{i+1};
0434         ph2 = ph2.';
0435         cy1top = reshape(cy1left, rcchk(i)*rty, rcy(i+1))*ph2;
0436         cy1top = reshape(cy1top, rcchk(i), rty, rcchk(i+1));
0437         cy1top = permute(cy1top, [1, 3, 2]);
0438         cy1top = reshape(cy1top, rcchk(i)*rcchk(i+1), rty);
0439         cy2 = yf{i}{L(i)};
0440         cy2 = reshape(cy2, rfy(L(i),i)*n(L(i),i), rty);
0441         cy2 = cy2*cy1top.';
0442         cy2 = reshape(cy2, rfy(L(i),i), n(L(i),i)*rcchk(i)*rcchk(i+1), 1);
0443         cyfr = yf{i};
0444         cyfr{L(i)} = cy2;
0445 
0446         x1 = xc{i};
0447         x1 = permute(x1, [2, 1, 3]);
0448         x1 = reshape(x1, rtx, rcx(i)*rcx(i+1));
0449         x2 = xf{i}{L(i)};
0450         x2 = reshape(x2, rfx(L(i),i)*n(L(i),i), rtx);
0451         x2 = x2*x1;
0452         x2 = reshape(x2, rfx(L(i),i), n(L(i),i)*rcx(i)*rcx(i+1), 1);
0453         xfr = xf{i};
0454         xfr{L(i)} = x2;
0455 
0456 %          chk1 = chkc{i};
0457 %          chk1 = permute(chk1, [2, 1, 3]);
0458 %          chk1 = reshape(chk1, rtchk, rcchk(i)*rcchk(i+1));
0459 %          chk2 = chkf{i}{L(i)};
0460 %          chk2 = reshape(chk2, rfchk(L(i),i)*n(L(i),i), rtchk);
0461 %          chk2 = chk2*chk1;
0462 %          chk2 = reshape(chk2, rfchk(L(i),i), n(L(i),i)*rcchk(i)*rcchk(i+1), 1);
0463 %          chkfr = chkf{i};
0464 %          chkfr{L(i)} = chk2;
0465 
0466         curn = xfr.n;
0467         curm = xfr.n;
0468         curm(L(i)) = n(L(i),i); % *rcchk(i)*rcchk(i+1);
0469         curra = Afr.r;
0470         curry = yfr.r;
0471         currx = xfr.r;
0472         currchk = [rfchk(1:L(i),i); 1];
0473 %          currchk = ones(L(i)+1,1);
0474 %          currchk = chkfr.r;
0475         % DMRG over the factor from L(i) to 1
0476         for j=L(i):-1:2
0477             % new
0478             if (j<L(i))
0479                 Phi2 = phfA{i}{j+1};
0480                 cPhi2 = cphfA{i}{j+1};
0481             else
0482                 Phi2 = phcA{i+1};
0483                 cPhi2 = cphcA{i+1};
0484             end;
0485             a2 = Afr{j};
0486             a1 = Afr{j-1};
0487             Phi1 = phfA{i}{j-1};
0488             ca2 = cAfr{j};
0489             ca1 = cAfr{j-1};
0490             cPhi1 = cphfA{i}{j-1};
0491             % old
0492 %             a2 = phfAt{i}{j+1};
0493 %             a2 = permute(a2, [2, 1, 3]);
0494 %             a2 = reshape(a2, curra(j+1), currx(j+1)*currx(j+1));
0495 %             a2 = reshape(Afr{j}, curra(j)*curn(j)*curn(j), curra(j+1))*a2;
0496 %             a2 = reshape(a2, curra(j), curn(j), curn(j), currx(j+1), currx(j+1));
0497 %             a2 = permute(a2, [2, 4, 3, 5, 1]);
0498 %             a2 = reshape(a2, curn(j)*currx(j+1), curn(j)*currx(j+1), curra(j));
0499 %             a1 = phfAb{i}{j-1};
0500 %             a1 = permute(a1, [1, 3, 2]);
0501 %             a1 = reshape(a1, currx(j-1)*currx(j-1), curra(j-1));
0502 %             a1 = a1*reshape(Afr{j-1}, curra(j-1), curn(j-1)*curn(j-1)*curra(j));
0503 %             a1 = reshape(a1, currx(j-1), currx(j-1), curn(j-1), curn(j-1), curra(j));
0504 %             a1 = permute(a1, [1, 3, 2, 4, 5]);
0505 %             a1 = reshape(a1, currx(j-1)*curn(j-1), currx(j-1)*curn(j-1), curra(j));
0506 
0507             y2 = phfy{i}{j+1}.';
0508             y2 = reshape(yfr{j}, curry(j)*curn(j), curry(j+1))*y2;
0509             y2 = reshape(y2, curry(j), curn(j)*currx(j+1));
0510             y2 = y2.';
0511 
0512             y1 = phfy{i}{j-1};
0513             y1 = y1*reshape(yfr{j-1}, curry(j-1), curn(j-1)*curry(j));
0514             y1 = reshape(y1, currx(j-1)*curn(j-1), curry(j));
0515 
0516             x2 = reshape(xfr{j}, currx(j), curn(j)*currx(j+1));
0517             x2 = x2.';
0518             x1 = reshape(xfr{j-1}, currx(j-1)*curn(j-1), currx(j));
0519 
0520             % check
0521             cy2 = cphfy{i}{j+1}.';
0522             cy2 = reshape(cyfr{j}, curry(j)*curm(j), curry(j+1))*cy2;
0523             cy2 = reshape(cy2, curry(j), curm(j)*currchk(j+1));
0524             cy1 = cphfy{i}{j-1};
0525             cy1 = cy1*reshape(cyfr{j-1}, curry(j-1), curm(j-1)*curry(j));
0526             cy1 = reshape(cy1, currchk(j-1)*curm(j-1), curry(j));
0527             cy = cy1*cy2;
0528 
0529             % new
0530             if (j==L(i))
0531                 currx(j+1)=rcx(i+1);
0532                 currchk(j+1)=rcchk(i+1);
0533                 curn(j) = n(L(i),i)*rcx(i);
0534                 curm(j) = n(L(i),i)*rcchk(i);
0535             end;
0536 
0537             if (verb>1)
0538                 fprintf('=rake_solve2= swp %d, factor {%d}{%d}, ', swp, i, j);
0539             end;
0540             local_format = 'full';
0541             if (currx(j-1)*curn(j-1)*curn(j)*currx(j+1)>max_full_size2)
0542                 local_format = 'tt';
0543             end;
0544             % old
0545 %             [u,s,v,r,dx_max,res_max]=local_solve(a1, a2, y1, y2, x1, x2, ...
0546 %                 currx(j-1), curn(j-1), curn(j), currx(j+1), curra(j), ...
0547 %                 tol/sqrt(L(i))/sqrt(d)/2, res_max, dx_max, ...
0548 %                 local_format, max_full_size, nrestart, gmres_iters, verb);
0549             % new
0550             [u,s,v,r,dx_max,res_max]=local_solve(Phi1,a1, a2, Phi2, y1, y2, x1, x2, ...
0551                 currx(j-1), curn(j-1), curn(j), currx(j+1), curra(j), ...
0552                 tol/sqrt(L(i))/sqrt(d), res_max, dx_max, resid_damp_loc, ...
0553                 local_format, max_full_size, nrestart, gmres_iters, verb);
0554         r = min(rmax, r);
0555         u = u(:,1:r); s = s(1:r,1:r); v = v(:,1:r);
0556             u = u*s;
0557 
0558             % check
0559             Asol = bfun3(cPhi1,ca1,ca2,cPhi2, u*(v.'));
0560             chk_res_max = max(chk_res_max, norm(Asol-cy(:))/norm(cy(:)));
0561 
0562             % kick
0563             if (~last_sweep)
0564 %                 Axprev = bfun3(Phi1,a1, a2, Phi2, x1*x2.');
0565 %                 Axprev = reshape(Axprev, currx(j-1)*curn(j-1), curn(j)*currx(j+1));
0566 %                 [unew,snew,vnew]=svd(Axprev, 'econ');
0567 %                 v = reort(v, vnew(:,1:min(kickrank, size(vnew,2))));
0568                 v = reort(v, randn(curn(j)*currx(j+1), kickrank));
0569             end;
0570             radd = size(v,2)-r;
0571             u = [u, zeros(currx(j-1)*curn(j-1), radd)];
0572             r = r+radd;
0573             xfr{j} = reshape(v.', r, curn(j), currx(j+1));
0574             xfr{j-1} = reshape(u, currx(j-1), curn(j-1), r);
0575             currx(j) = r;
0576             rfx(j,i)=r;
0577             r_max = max(r_max, r);
0578             % old
0579             % new phis
0580 %             a2 = permute(a2, [1, 3, 2]);
0581 %             a2 = reshape(a2, curn(j)*currx(j+1)*curra(j), curn(j)*currx(j+1));
0582 %             a2 = a2*v;
0583 %             a2 = reshape(a2, curn(j)*currx(j+1), curra(j)*r);
0584 %             a2 = (v')*a2;
0585 %             phfAt{i}{j} = reshape(a2, r, curra(j), r);
0586             phfy{i}{j} = (v')*y2;
0587             %new
0588             phfA{i}{j} = compute_next_Phi(Phi2, xfr{j}, a2, xfr{j}, 'rl');
0589 
0590             % check vector
0591             ccr = ones(n(j,i), 1);
0592 %              [ccr,rv]=qr(ccr.', 0);
0593 %              rfchk(j,i) = size(ccr,2);
0594             cphfy{i}{j} = (ccr')*(cy2.');
0595 %              ccr = reshape(ccr.', rfchk(j,i), curm(j), currchk(j+1));
0596             cphfA{i}{j} = compute_next_Phi(cPhi2, ccr.', ca2, xfr{j}, 'rl');
0597         end;
0598 
0599 %         fprintf('=rake_solve========= factor {%d} processing, Sweep %d =======\n', i, swp);
0600 %         xfrold = xfr;
0601 %         xfr = dmrg_solve2(Afr, yfr, tol, 'x0', xfrold, 'max_full_size', max_full_size, 'nrestart', nrestart, 'gmres_iters', gmres_iters, 'nswp', 1);
0602 %         res = norm(Afr*xfrold-yfr)/norm(yfr);
0603 %         dx = norm(xfrold-xfr)/norm(xfr);
0604 %         dx_max = max(dx_max, dx);
0605 %         r_max = max([r_max; rank(xfr)]);
0606 %         fprintf('=rake_solve========= factor {%d} res_prev: %3.3e, dx: %3.3e, rmax: %d\n', i, res, dx, max(rank(xfr)));
0607 %         res_max = max(res_max, res);
0608 %         rfx(1:L(i)+1,i) = xfr.r;
0609 
0610         % We have to split the tucker block, and compute new phf*b
0611         for j=1:L(i)
0612             cr = xfr{j};
0613             % What is our next core?
0614             if (j<L(i))
0615                 % we are still on a "tooth"
0616                 cr = reshape(cr, rfx(j,i)*n(j,i), rfx(j+1,i));
0617                 [cr, rv] = qr(cr, 0);
0618                 cr2 = xfr{j+1};
0619                 ncur = size(cr2, 2);
0620                 r3cur = size(cr2, 3);
0621                 cr2 = reshape(cr2, rfx(j+1,i), ncur*r3cur);
0622                 cr2 = rv*cr2;
0623                 rfx(j+1,i) = size(cr, 2);
0624                 xfr{j} = reshape(cr, rfx(j,i), n(j,i), rfx(j+1,i));
0625                 xfr{j+1} = reshape(cr2, rfx(j+1,i), ncur, r3cur);
0626             else
0627                 % We have to split the tucker core and update the tucker
0628                 % rank
0629                 cr = reshape(cr, rfx(j,i)*n(j,i), rcx(i)*rcx(i+1));
0630                 [u,s,v]=svd(cr, 'econ');
0631                 % Prepare the local matrix and rhs for residue check
0632                 % new
0633                 Phi2 = phcA{i+1};
0634                 curA2 = reshape(a1left, rta, rcx(i), rcx(i), rcA(i+1));
0635                 curA1 = Af{i}{L(i)};
0636                 Phi1 = phfA{i}{j};
0637                 % old
0638 %                 curA = cell(2,1);
0639 %                 curA{2} = reshape(a1top, rcx(i)*rcx(i+1), rcx(i)*rcx(i+1), rta);
0640 %                 curA{1} = Af{i}{L(i)};
0641 %                 curA{1} = reshape(curA{1}, rfA(j,i), n(j,i)*n(j,i)*rta);
0642 %                 ph2 = phfAb{i}{j};
0643 %                 ph2 = permute(ph2, [1,3,2]);
0644 %                 ph2 = reshape(ph2, rfx(j,i)*rfx(j,i), rfA(j,i));
0645 %                 curA{1} = ph2*curA{1};
0646 %                 curA{1} = reshape(curA{1}, rfx(j,i), rfx(j,i), n(j,i), n(j,i), rta);
0647 %                 curA{1} = permute(curA{1}, [1, 3, 2, 4, 5]);
0648 %                 curA{1} = reshape(curA{1}, rfx(j,i)*n(j,i), rfx(j,i)*n(j,i), rta);
0649 
0650                 rhs = reshape(yf{i}{j}, rfy(j,i), n(j,i)*rty);
0651                 rhs = phfy{i}{j}*rhs;
0652                 rhs = reshape(rhs, rfx(j,i)*n(j,i), rty);
0653                 rhs = rhs*(y1top.');
0654                 rhs = reshape(rhs, rfx(j,i)*n(j,i)*rcx(i)*rcx(i+1),1);
0655                 r = 1;
0656                 normy = norm(rhs);
0657                 % old
0658 %                 res_true = norm(bfun2(curA, cr, rfx(j,i), n(j,i), rcx(i), rcx(i+1), rfx(j,i), n(j,i), rcx(i), rcx(i+1))-rhs)/normy;
0659                 % new
0660                 res_true = norm(bfun3(Phi1, curA1, curA2, Phi2, cr)-rhs)/normy;
0661                 while (r<=size(s,1))
0662                     cursol = u(:,1:r)*s(1:r,1:r)*(v(:,1:r)');
0663                     % old
0664 %                     res = norm(bfun2(curA, cursol, rfx(j,i), n(j,i), rcx(i), rcx(i+1), rfx(j,i), n(j,i), rcx(i), rcx(i+1))-rhs)/normy;
0665                     % new
0666                     res = norm(bfun3(Phi1, curA1, curA2, Phi2, cursol)-rhs)/normy;
0667                     if (res<max(tol/sqrt(L(i)), res_true*2))
0668                         break;
0669                     end;
0670                     r = r+1;
0671                 end;
0672                 if (verb>1)
0673                     fprintf('=rake_solve2= swp %d, tuckerrank {%d}, res: %3.3e, r: %d\n', swp, i, res, r);
0674                 end;
0675         r = min(rmax, r);
0676                 u = u(:,1:r);
0677                 v = conj(v(:,1:r));
0678                 s = s(1:r,1:r);
0679                 v = v*s;
0680                 if (~last_sweep)
0681 %                     Axprev = bfun3(Phi1, curA1, curA2, Phi2, cr);
0682 %                     Axprev = reshape(Axprev, rfx(j,i)*n(j,i), rcx(i)*rcx(i+1));
0683 %                     [unew,snew,vnew]=svd(Axprev, 'econ');
0684 %                     u = reort(u, unew(:,1:min(kickrank, size(unew,2))));
0685                     u = reort(u, randn(rfx(j,i)*n(j,i), kickrank));
0686                 end;
0687                 radd = size(u,2)-r;
0688                 v = [v, zeros(rcx(i)*rcx(i+1), radd)];
0689                 r = r+radd;
0690                 rfx(j+1,i) = r;
0691                 cr = u;
0692                 xfr{j} = reshape(cr, rfx(j,i), n(j,i), r);
0693                 v = reshape(v.', r, rcx(i), rcx(i+1));
0694                 v = permute(v, [2, 1, 3]);
0695                 xc{i} = v;
0696                 r_max = max(r_max, r);
0697             end;
0698             % Update bottom phis
0699             cr = reshape(cr, rfx(j,i), n(j,i), rfx(j+1,i));
0700             if (j<L(i))
0701                 phfA{i}{j+1} = compute_next_Phi(phfA{i}{j}, cr, Af{i}{j}, cr, 'lr');
0702                 phfy{i}{j+1} = compute_next_Phi(phfy{i}{j}, cr, [], yf{i}{j}, 'lr');
0703             else
0704                 phAfc{i} = compute_next_Phi(phfA{i}{j}, cr, Af{i}{j}, cr, 'lr');
0705                 phyfc{i} = compute_next_Phi(phfy{i}{j}, cr, [], yf{i}{j}, 'lr');
0706             end;
0707 
0708             % check vector
0709 %              ccr = ones(rfchk(j,i)*n(j,i), checkrank);
0710 %              [ccr, rv] = qr(ccr, 0);
0711 %              rfchk(j+1,i) = size(ccr, 2);
0712 %              ccr = reshape(ccr, rfchk(j,i), n(j,i), rfchk(j+1,i));
0713             ccr = ones(1, n(j,i), 1);
0714             if (j<L(i))
0715                 cphfA{i}{j+1} = compute_next_Phi(cphfA{i}{j}, ccr, Af{i}{j}, cr, 'lr');
0716                 cphfy{i}{j+1} = compute_next_Phi(cphfy{i}{j}, ccr, [], yf{i}{j}, 'lr');
0717             else
0718                 cphAfc{i} = compute_next_Phi(cphfA{i}{j}, ccr, Af{i}{j}, cr, 'lr');
0719                 cphyfc{i} = compute_next_Phi(cphfy{i}{j}, ccr, [], yf{i}{j}, 'lr');
0720             end;
0721         end;
0722         xf{i} = xfr;
0723 
0724         % Now, perform the DMRG over the core to the next factor, update phc*l as well
0725         if (i<d)
0726             rtx = rfx(L(i)+1,i); rtx2 = rfx(L(i+1)+1,i+1);
0727             rx1 = rcx(i); rx2 = rcx(i+1); rx3 = rcx(i+2);
0728             ra2 = rcA(i+1); ra3 = rcA(i+2);
0729             ry2 = rcy(i+1); ry3 = rcy(i+2);
0730 
0731             rtchk = rfchk(L(i)+1,i); rtchk2 = rfchk(L(i+1)+1,i+1);
0732             rchk1 = rcchk(i); rchk3 = rcchk(i+2);
0733             % Acr and ycr are okey, except the i-th core
0734             % new
0735             Phi1 = phcA{i};
0736             a1 = core_matrix(Ac{i}, phAfc{i});
0737             a2 = Acr{i+1};
0738             Phi2 = phcA{i+2};
0739 
0740             cPhi1 = cphcA{i};
0741             ca1 = core_matrix(Ac{i}, cphAfc{i});
0742             ca2 = cAcr{i+1};
0743             cPhi2 = cphcA{i+2};
0744             % old
0745 %             a1 = reshape(a1left, rx1*rx1, rta, ra2);
0746 %             a1 = core_matrix(a1, phAfc{i});
0747 %             a1 = reshape(a1, rx1, rx1, rtx, rtx, ra2);
0748 %             a1 = permute(a1, [1, 3, 2, 4, 5]);
0749 %             a1 = reshape(a1, rx1*rtx, rx1*rtx, ra2);
0750 %             a2 = phcAr{i+2};
0751 %             a2 = permute(a2, [2, 1, 3]);
0752 %             a2 = reshape(a2, ra3, rx3*rx3);
0753 %             a2 = reshape(Acr{i+1}, ra2*rtx2*rtx2, ra3)*a2;
0754 %             a2 = reshape(a2, ra2, rtx2, rtx2, rx3, rx3);
0755 %             a2 = permute(a2, [2, 4, 3, 5, 1]);
0756 %             a2 = reshape(a2, rtx2*rx3, rtx2*rx3, ra2);
0757 
0758             y1 = reshape(y1left, rx1, rty, ry2);
0759             y1 = core_vector(y1, phyfc{i});
0760             y1 = reshape(y1, rx1*rtx, ry2);
0761 
0762             y2 = phcy{i+2}.';
0763             y2 = reshape(ycr{i+1}, ry2*rtx2, ry3)*y2;
0764             y2 = reshape(y2, ry2, rtx2*rx3);
0765             y2 = y2.';
0766 
0767             x1 = reshape(xc{i}, rx1*rtx, rx2);
0768             x2 = reshape(xc{i+1}, rx2, rtx2*rx3);
0769             x2 = x2.';
0770 
0771             if (verb>1)
0772                 fprintf('=rake_solve2= swp %d, core {%d}, ', swp, i);
0773             end;
0774             local_format = 'full';
0775             if (rx1*rtx*rtx2*rx3>max_full_size2)
0776                 local_format = 'tt';
0777             end;
0778             % new
0779             [u,s,v,r,dx_max,res_max]=local_solve(Phi1, a1, a2, Phi2, y1, y2, x1, x2, ...
0780                 rx1, rtx, rtx2, rx3, ra2, ...
0781                 tol/sqrt(d), res_max, dx_max, resid_damp_loc,...
0782                 local_format, max_full_size, nrestart, gmres_iters, verb);
0783             % old
0784 %             [u,s,v,r,dx_max,res_max]=local_solve(a1, a2, y1, y2, x1, x2, ...
0785 %                 rx1, rtx, rtx2, rx3, ra2, ...
0786 %                 tol/sqrt(d)/2, res_max, dx_max, ...
0787 %                 local_format, max_full_size, nrestart, gmres_iters, verb);
0788             r = min(rmax, r);
0789             u = u(:,1:r); s = s(1:r,1:r); v = v(:,1:r);
0790             v = v*s;
0791 
0792             % check
0793             Asol = bfun3(cPhi1, ca1, ca2, cPhi2, u*(v.'));
0794             cy1 = reshape(cy1left, rchk1, rty, ry2);
0795             cy1 = core_vector(cy1, cphyfc{i});
0796             cy1 = reshape(cy1, rchk1*rtchk, ry2);
0797             cy2 = cphcy{i+2}.';
0798             cy2 = reshape(cycr{i+1}, ry2*rtchk2, ry3)*cy2;
0799             cy2 = reshape(cy2, ry2, rtchk2*rchk3);
0800             cy = cy1*cy2;
0801             chk_res_max = max(chk_res_max, norm(Asol-cy(:))/norm(cy(:)));
0802 
0803             % kick
0804             if (~last_sweep)
0805 %                 Axprev = bfun3(Phi1, a1, a2, Phi2, x1*x2.');
0806 %                 Axprev = reshape(Axprev, rx1*rtx, rtx2*rx3);
0807 %                 [unew,snew,vnew]=svd(Axprev, 'econ');
0808 %                 u = reort(u, unew(:,1:min(kickrank, size(unew,2))));
0809                 u = reort(u, randn(rx1*rtx, kickrank));
0810             end;
0811             radd = size(u,2)-r;
0812             v = [v, zeros(rtx2*rx3, radd)];
0813             r = r+radd;
0814             xc{i} = reshape(u, rx1, rtx, r);
0815             xc{i+1} = reshape(v.', r, rtx2, rx3);
0816             rcx(i+1)=r;
0817             r_max = max(r_max, r);
0818             % new phis
0819             % old
0820 %             a1 = permute(a1, [1, 3, 2]);
0821 %             a1 = reshape(a1, rx1*rtx*ra2, rx1*rtx);
0822 %             a1 = a1*u;
0823 %             a1 = reshape(a1, rx1*rtx, ra2*r);
0824 %             a1 = (u')*a1;
0825 %             phcAl{i+1} = reshape(a1, r, ra2, r);
0826             phcy{i+1} = (u')*y1;
0827             % new
0828             phcA{i+1} = compute_next_Phi(phcA{i}, xc{i}, a1, xc{i}, 'lr');
0829 
0830 %              ccr = ones(rchk1*rtchk, checkrank);
0831 %              [ccr,rv]=qr(ccr,0);
0832 %              rcchk(i+1) = size(ccr, 2);
0833             ccr = 1;
0834             cphcy{i+1} = (ccr')*cy1;
0835             ccr = reshape(ccr, rchk1, rtchk, rcchk(i+1));
0836             cphcA{i+1} = compute_next_Phi(cphcA{i}, ccr, ca1, xc{i}, 'lr');
0837 
0838 
0839 %             cr = xc{i};
0840 %             cr = reshape(cr, rcx(i)*rtx, rcx(i+1));
0841 %             [cr, rv]=qr(cr, 0);
0842 %             cr2 = xc{i+1};
0843 %             cr2 = reshape(cr2, rcx(i+1), rfx(L(i+1)+1,i+1)*rcx(i+2));
0844 %             cr2 = rv*cr2;
0845 %             rcx(i+1) = size(cr, 2);
0846 %             xc{i+1} = reshape(cr2, rcx(i+1), rfx(L(i+1)+1,i+1), rcx(i+2));
0847 %             xc{i} = reshape(cr, rcx(i), rtx, rcx(i+1));
0848 %
0849 %             % We have a1left, y1left of sizes rcx(i)(^2), rtuck, rc*(i+1)
0850 %             curph = reshape(a1left, rcx(i)*rcx(i), rta, rcA(i+1));
0851 %             curph = permute(curph, [2, 1, 3]);
0852 %             curph = reshape(curph, rta, rcx(i)*rcx(i)*rcA(i+1));
0853 %             ph2 = phfAb{i};
0854 %             ph2 = permute(ph2, [1, 3, 2]);
0855 %             ph2 = reshape(ph2, rtx*rtx, rta);
0856 %             curph = ph2*curph;
0857 %             curph = reshape(curph, rtx, rtx, rcx(i), rcx(i), rcA(i+1));
0858 %             curph = permute(curph, [3, 1, 5, 4, 2]);
0859 %             curph = reshape(curph, rcx(i)*rtx*rcA(i+1), rcx(i)*rtx);
0860 %             curph = curph*cr;
0861 %             curph = reshape(curph, rcx(i)*rtx, rcA(i+1)*rcx(i+1));
0862 %             curph = (cr')*curph;
0863 %             phcAl{i+1} = reshape(curph, rcx(i+1), rcA(i+1), rcx(i+1));
0864 %
0865 %             curph = reshape(y1left, rcx(i), rty, rcy(i+1));
0866 %             curph = permute(curph, [2, 1, 3]);
0867 %             curph = reshape(curph, rty, rcx(i)*rcy(i+1));
0868 %             ph2 = phfyb{i};
0869 %             curph = ph2*curph;
0870 %             curph = reshape(curph, rtx, rcx(i), rcy(i+1));
0871 %             curph = permute(curph, [2, 1, 3]);
0872 %             curph = reshape(curph, rcx(i)*rtx, rcy(i+1));
0873 %             curph = (cr')*curph;
0874 %             phcyl{i+1} = curph;
0875         end;
0876     end;
0877 
0878     if (verb>0)
0879         real_res = NaN;
0880 %      x = qtt_tucker;
0881 %      x.dphys = d;
0882 %      x.tuck = xf;
0883 %      x.core = xc;
0884 %
0885 %      real_res = mvrk(A, x, tol, 'verb', 0);
0886 %  %          real_res = A*x;
0887 %      real_res = norm(real_res-y)/norm(y);
0888         fprintf('\n=rake_solve2= swp %d, dx_max: %3.3e, res_max: %3.3e, r_max: %d, real_res: %3.3e, chk_res: %3.3e\n\n', swp, dx_max, res_max, r_max, real_res, chk_res_max);
0889     end;
0890     if (last_sweep)
0891         break;
0892     end;
0893     if (chk_res_max<tol)||(swp==nswp-1)
0894 %      if (res_max<tol*resid_damp_glob)||(swp==nswp-1)
0895 %          last_sweep = true;
0896         break;
0897     end;
0898 end;
0899 
0900 x = qtt_tucker;
0901 x.dphys = d;
0902 x.tuck = xf;
0903 x.core = xc;
0904 end
0905 
0906 
0907 function [Phi] = compute_next_Phi(Phi_prev, x, A, y, direction)
0908 % Performs the recurrent Phi (or Psi) matrix computation
0909 % Phi = Phi_prev * (x'Ay).
0910 % If direction is 'lr', computes Psi
0911 % if direction is 'rl', computes Phi
0912 % A can be empty, then only x'y is computed.
0913 
0914 if (strcmp(direction, 'rl'))
0915   % Revert ranks to perform the right-to-left recursion
0916   x = permute(x, [3, 2, 1]);
0917   y = permute(y, [3, 2, 1]);
0918   if (~isempty(A))
0919     A = permute(A, [4, 2, 3, 1]);
0920   end
0921 end
0922 
0923 rx1 = size(x,1); n = size(x,2); rx2 = size(x,3);
0924 ry1 = size(y,1); m = size(y,2); ry2 = size(y,3);
0925 if (~isempty(A))
0926   ra1 = size(A,1); ra2 = size(A,4);
0927 else
0928   ra1 = 1; ra2 = 1;
0929 end
0930 
0931 Phi = reshape(Phi_prev, [rx1*ra1, ry1]);
0932 y = reshape(y, [ry1, m*ry2]);
0933 Phi = Phi*y;    % complexity §\mcommentfont$\mathcal{O}(n  r_x r_A r_y^2)$§
0934 Phi = reshape(Phi, [rx1, ra1, m, ry2]);
0935 Phi = permute(Phi, [2, 3, 1, 4]);
0936 if (~isempty(A))
0937   Phi = reshape(Phi, [ra1*m, rx1*ry2]);
0938   A = permute(A, [4, 2, 1, 3]);
0939   A = reshape(A, [ra2*n, ra1*m]);
0940   Phi = A*Phi;    % complexity §\mcommentfont$\mathcal{O}(n^2  r_x r_A^2 r_y)$§
0941   Phi = reshape(Phi, [ra2, n, rx1, ry2]);
0942 end
0943 Phi = permute(Phi, [3, 2, 1, 4]);
0944 Phi = reshape(Phi, [rx1*n, ra2*ry2]);
0945 x = reshape(x, [rx1*n, rx2]);
0946 Phi = (x')*Phi;    % complexity §\mcommentfont$\mathcal{O}(n  r_x^2 r_A r_y)$§
0947 if (~isempty(A))
0948   Phi = reshape(Phi, [rx2, ra2, ry2]);
0949 end
0950 end
0951 
0952 function [A] = core_matrix(core_block, Phi_factor)
0953 % Computes the matrix to use in DMRG by the core, by convolving the core of
0954 % the Tucker block and the last factor Phi matrix over the Tucker rank
0955 % core_block is rc1 x rtuck x rc2, and
0956 % Phi_factor is n1 x rtuck x n2
0957 % A is then rc1 x n1 x n2 x rc2
0958 
0959 r1 = size(core_block, 1); rtuck = size(core_block, 2); r2 = size(core_block, 3);
0960 n1 = size(Phi_factor, 1); n2 = size(Phi_factor, 3);
0961 
0962 Phi_factor = permute(Phi_factor, [1, 3, 2]);
0963 Phi_factor = reshape(Phi_factor, n1*n2, rtuck);
0964 A = permute(core_block, [2, 1, 3]);
0965 A = reshape(A, rtuck, r1*r2);
0966 A = Phi_factor*A;
0967 A = reshape(A, n1, n2, r1, r2);
0968 A = permute(A, [3, 1, 2, 4]);
0969 
0970 end
0971 
0972 
0973 function [A] = core_vector(core_block, Phi_factor)
0974 % Computes the vector to use in DMRG by the core, by convolving the core of
0975 % the Tucker block and the last factor Phi matrix over the Tucker rank
0976 % core_block is rc1 x rtuck x rc2, and
0977 % Phi_factor is n1 x rtuck
0978 % A is then rc1 x n1 x rc2
0979 
0980 r1 = size(core_block, 1); rtuck = size(core_block, 2); r2 = size(core_block, 3);
0981 n1 = size(Phi_factor, 1);
0982 
0983 A = permute(core_block, [2, 1, 3]);
0984 A = reshape(A, rtuck, r1*r2);
0985 A = Phi_factor*A;
0986 A = reshape(A, n1, r1, r2);
0987 A = permute(A, [2, 1, 3]);
0988 
0989 end
0990 
0991 
0992 function [y]=bfun2(B, x, rxm1, m1, m2, rxm3, rxn1, k1, k2, rxn3)
0993 % Computes (B{1} \otimes B{2})x
0994 % B{1} is of sizes rxn1*k1, rxm1*m1, rB
0995 % B{2} is of sizes k2*rxn3, m2*rxm3, rB
0996 rB=size(B{1},3);
0997 x = reshape(x, rxm1*m1, m2*rxm3);
0998 B1 = permute(B{1}, [3 1 2]);
0999 B1 = reshape(B1, rB*rxn1*k1, rxm1*m1);
1000 y = B1*x; % size rB*rxn1*k1,m2*rxm3  % cplx rB*rx^3*n^3
1001 y = reshape(y, rB, rxn1*k1, m2*rxm3);
1002 y = permute(y, [3 1 2]);
1003 y = reshape(y, m2*rxm3*rB, rxn1*k1);
1004 B2 = reshape(B{2}, k2*rxn3, m2*rxm3*rB);
1005 y = B2*y; % size k2*rxn3,rxn1*k1 % cplx rB*rx^3*n^3
1006 y = reshape(y.', rxn1*k1*k2*rxn3, 1);
1007 end
1008 
1009 
1010 function [y]=bfun3(Phi1,B1,B2,Phi2, x)
1011 % Computes (Phi1 * B1 * B2 * Phi2)*x
1012 % Phi1 is of sizes ry1, rB1, rx1
1013 % B1 is of sizes rB1, k1, m1, rB2
1014 % B2 is of sizes rB2, k2, m2, rB3
1015 % Phi2 is of sizes ry3, rB3, rx3
1016 ry1 = size(Phi1,1); ry3 = size(Phi2,1);
1017 rx1 = size(Phi1,3); rx3 = size(Phi2,3);
1018 rb1=size(B1,1); rb2=size(B1,4); rb3 = size(B2, 4);
1019 m1 = size(B1,3); m2 = size(B2,3);
1020 k1 = size(B1,2); k2 = size(B2,2);
1021 
1022 y = reshape(x, rx1, m1*m2*rx3);
1023 Phi1 = reshape(Phi1, ry1*rb1, rx1);
1024 y = Phi1*y; % size ry1*rb1,m1*m2*rx3 % cplx rb*rx^3*m^2
1025 y = reshape(y, ry1, rb1*m1, m2, rx3);
1026 y = permute(y, [2, 1, 3, 4]);
1027 y = reshape(y, rb1*m1, ry1*m2*rx3);
1028 B1 = permute(B1, [2, 4, 1, 3]);
1029 B1 = reshape(B1, k1*rb2, rb1*m1);
1030 y = B1*y; % size k1*rb2, ry1*m2*rx3 % cplx rb^2*rx^2*n^3
1031 y = reshape(y, k1, rb2, ry1, m2, rx3);
1032 y = permute(y, [2, 4, 3, 1, 5]);
1033 y = reshape(y, rb2*m2, ry1*k1*rx3);
1034 B2 = permute(B2, [2, 4, 1, 3]);
1035 B2 = reshape(B2, k2*rb3, rb2*m2);
1036 y = B2*y; % size k2*rb3, ry1*k1*rx3 % cplx rb^2*rx^2*n^3
1037 y = reshape(y, k2, rb3, ry1*k1, rx3);
1038 y = permute(y, [2, 4, 3, 1]);
1039 y = reshape(y, rb3*rx3, ry1*k1*k2);
1040 Phi2 = reshape(Phi2, ry3, rb3*rx3);
1041 y = Phi2*y; % size ry3, ry1*k1*k2 % cplx rb*rx^3*n^2
1042 y = y.';
1043 y = reshape(y, ry1*k1*k2*ry3, 1);
1044 end
1045 
1046 
1047 % old
1048 % function [u,s,v,r,dx_max,res_max]=local_solve(a1, a2, y1, y2, x1, x2, ...
1049 % new
1050 function [u,s,v,r,dx_max,res_max]=local_solve(Phi1,a1, a2, Phi2, y1, y2, x1, x2, ...
1051     rx1, n1, n2, rx3, ra2, ...
1052     real_tol, res_max, dx_max,  resid_damp, ...
1053     local_format, max_full_size, nrestart, gmres_iters, verb)
1054 
1055 if (strcmp(local_format, 'full'))
1056     sol_prev = x1*(x2.');
1057     sol_prev = sol_prev(:);
1058     rhs = y1*(y2.');
1059     rhs = rhs(:);
1060     normy = norm(rhs);
1061     if (rx1*n1*n2*rx3<max_full_size)
1062         % new
1063         B = permute(Phi1, [1,3,2]);
1064         B = reshape(B, rx1*rx1, size(a1,1));
1065         a1 = reshape(a1, size(a1,1), n1*n1*ra2);
1066         B = B*a1;
1067         B = reshape(B, rx1, rx1, n1, n1, ra2);
1068         B = permute(B, [1, 3, 2, 4, 5]);
1069         B = reshape(B, rx1*n1*rx1*n1, ra2);
1070         ra3 = size(a2,4);
1071         a2 = reshape(a2, ra2, n2*n2*ra3);
1072         B = B*a2;
1073         B = reshape(B, rx1*n1*rx1*n1*n2*n2, ra3);
1074         Phi2 = permute(Phi2, [2, 1, 3]);
1075         Phi2 = reshape(Phi2, ra3, rx3*rx3);
1076         B = B*Phi2;
1077         B = reshape(B, rx1*n1, rx1*n1, n2, n2, rx3, rx3);
1078         B = permute(B, [1, 3, 5, 2, 4, 6]);
1079         % old
1080 %         B = reshape(a1, rx1*n1*rx1*n1, ra2);
1081 %         B = B*(reshape(a2, n2*rx3*n2*rx3, ra2).');
1082 %         B = reshape(B, rx1*n1,  rx1*n1, n2*rx3, n2*rx3);
1083 %         B = permute(B, [1, 3, 2, 4]);
1084 
1085         B = reshape(B, rx1*n1*n2*rx3, rx1*n1*n2*rx3);
1086          sol = (B'*B) \ (B'*rhs);
1087 %          sol = B \ rhs;
1088          
1089          % If the direct solution sucked
1090         [sol,flg]=gmres(B, rhs, ...
1091             nrestart, real_tol/resid_damp, gmres_iters, [], [], sol);
1092         if (flg>0)
1093             fprintf('--warn-- gmres did not converge\n');
1094         end;         
1095         
1096 %          sol = (B'*B + 1e-16*norm(B'*B, 'fro')) \ (B'*rhs);
1097 %         sol = pinv(B)*rhs;
1098         res_true = norm(B*sol-rhs)/normy;
1099         res_prev = norm(B*sol_prev-rhs)/normy;
1100     else
1101         B = cell(2,1);
1102         B{1} = a1;
1103         B{2} = a2;
1104         % old
1105 %         res_prev = norm(bfun2(B, sol_prev, rx1, n1, n2, rx3, rx1, n1, n2, rx3)-rhs)/normy;
1106         % new
1107         res_prev = norm(bfun3(Phi1, a1, a2, Phi2, sol_prev)-rhs)/normy;
1108 
1109 %         sol = als_solve_rx_2(B, rhs, real_tol, 10, sol_prev, [], 3);
1110 %         [sol,flg]=bicgstab(@(v)bfun2(B, v, rx1, n1, n2, rx3, rx1, n1, n2, rx3), rhs, ...
1111 %             max(real_tol,res_prev*0.1), nrestart*gmres_iters, [], [], sol_prev);
1112         % old
1113 %         [sol,flg]=gmres(@(v)bfun2(B, v, rx1, n1, n2, rx3, rx1, n1, n2, rx3), rhs, ...
1114 %             nrestart, max(real_tol,res_prev*0.05), gmres_iters, [], [], sol_prev);
1115         % new
1116         [sol,flg]=gmres(@(v)bfun3(Phi1, a1, a2, Phi2, v), rhs, ...
1117             nrestart, max(real_tol/resid_damp,res_prev*0), gmres_iters, [], [], sol_prev);
1118         if (flg>0)
1119             fprintf('--warn-- gmres did not converge\n');
1120         end;
1121         % old
1122 %         res_true = norm(bfun2(B, sol, rx1, n1, n2, rx3, rx1, n1, n2, rx3)-rhs)/normy;
1123         % new
1124         res_true = norm(bfun3(Phi1, a1, a2, Phi2, sol)-rhs)/normy;
1125     end;
1126 
1127     if ((res_prev/res_true)<resid_damp)&&(res_true>real_tol/resid_damp)
1128         fprintf('--warn-- the residual damp by gmres was smaller than in the truncation\n');
1129     end;
1130 
1131     dx = norm(sol-sol_prev)/norm(sol);
1132     dx_max = max(dx_max, dx);
1133     if (rx1*n1*n2*rx3<max_full_size)
1134         Bx = B*sol; Bx_prev = B*sol_prev;
1135     else
1136         Bx = bfun3(Phi1, a1, a2, Phi2, sol);
1137         Bx_prev = bfun3(Phi1, a1, a2, Phi2, sol_prev);
1138     end;
1139     res_max = max(res_max, norm(Bx-Bx_prev)/norm(Bx));
1140 %      res_max = max(res_max, res_prev);
1141 
1142     sol = reshape(sol, rx1*n1, n2*rx3);
1143     [u,s,v]=svd(sol, 'econ');
1144     s = diag(s);
1145     r1 = 1; r2 = numel(s); r = round((r1+r2)/2);
1146     while (r2-r1>1)
1147         cursol = u(:,1:r)*diag(s(1:r))*(v(:,1:r)');
1148         if (rx1*n1*n2*rx3<max_full_size)
1149             res = norm(B*cursol(:)-rhs)/normy;
1150         else
1151             % old
1152 %             res = norm(bfun2(B, cursol, rx1, n1, n2, rx3, rx1, n1, n2, rx3)-rhs)/normy;
1153             % new
1154             res = norm(bfun3(Phi1, a1, a2, Phi2, cursol)-rhs)/normy;
1155         end;
1156         if (res<max(real_tol, res_true*resid_damp))
1157             r2 = r;
1158         else
1159             r1 = r;
1160         end;
1161         r = round((r1+r2)/2);
1162     end;
1163 %     r = 1;
1164     while (r<=numel(s))
1165         cursol = u(:,1:r)*diag(s(1:r))*(v(:,1:r)');
1166         if (rx1*n1*n2*rx3<max_full_size)
1167             res = norm(B*cursol(:)-rhs)/normy;
1168         else
1169             % old
1170 %             res = norm(bfun2(B, cursol, rx1, n1, n2, rx3, rx1, n1, n2, rx3)-rhs)/normy;
1171             % new
1172             res = norm(bfun3(Phi1, a1, a2, Phi2, cursol)-rhs)/normy;
1173         end;
1174         if (res<max(real_tol, res_true*resid_damp))
1175             break;
1176         end;
1177         r = r+1;
1178     end;
1179     r = min(r, numel(s));
1180     if (verb>1)
1181         fprintf('dx: %3.3e, res: %3.3e, res_prev: %3.3e, r: %d\n', dx, res, res_prev, r);
1182     end;
1183 
1184     s = diag(s(1:r));
1185     u = u(:,1:r);
1186     v = conj(v(:,1:r));
1187 else
1188     % implement tt-gmres here
1189     B = cell(2,1);
1190     B{1} = a1;
1191     B{2} = a2;
1192 %     iB = tt_minres_selfprec(B, 1e-1, 1e-2, 10, 'right');
1193     iB = [];
1194     sol_prev = cell(2,1);
1195     sol_prev{1} = x1;
1196     sol_prev{2} = x2;
1197     rhs = cell(2,1);
1198     rhs{1} = y1;
1199     rhs{2} = y2;
1200     normy = tt_dist3(rhs, tt_scal(rhs,0));
1201     drhs = tt_mv(B, sol_prev);
1202     res_prev = tt_dist3(drhs, rhs)/normy;
1203     drhs = tt_add(rhs, tt_scal(drhs, -1));
1204     drhs = tt_compr2(drhs, real_tol);
1205     dsol = tt_gmres(B, drhs, real_tol*resid_damp_loc/res_prev, gmres_iters, nrestart, real_tol, real_tol, iB, [], [], [], 1);
1206 %     sol = tt_gmres(B, rhs, real_tol*2, gmres_iters*5, nrestart/5, real_tol, real_tol, iB, [], [], sol_prev, 1);
1207     sol = tt_add(sol_prev, dsol);
1208     sol = tt_compr2(sol, real_tol);
1209     normsol = tt_dist3(sol, tt_scal(sol,0));
1210     dx = tt_dist3(sol, sol_prev)/normsol;
1211     res = tt_dist3(tt_mv(B, sol), rhs)/normy;
1212 
1213     dx_max = max(dx_max, dx);
1214     res_max = max(res_max, res_prev);
1215     [v, s]=qr(sol{2}, 0);
1216     sol{1} = sol{1}*(s.');
1217     [u, s]=qr(sol{1}, 0);
1218     r = size(sol{1},2);
1219     if (verb>1)
1220         fprintf('dx: %3.3e, res: %3.3e, res_prev: %3.3e, r: %d\n', dx, res, res_prev, r);
1221     end;
1222 end;
1223 end

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