Home > matlab > solve_two_boundaries.m

solve_two_boundaries

PURPOSE ^

Computes the deterministic simulation of a block of equation containing

SYNOPSIS ^

function y = solve_two_boundaries(fname, y, x, params, steady_state, y_index, nze, periods, y_kmin_l, y_kmax_l, is_linear, Block_Num, y_kmin, maxit_, solve_tolf, lambda, cutoff, stack_solve_algo, M, oo)

DESCRIPTION ^

 Computes the deterministic simulation of a block of equation containing
 both lead and lag variables using relaxation methods 

 INPUTS
   fname               [string]        name of the file containing the block
                                       to simulate
   y                   [matrix]        All the endogenous variables of the model
   x                   [matrix]        All the exogenous variables of the model
   params              [vector]        All the parameters of the model
   steady_state        [vector]        steady state of the model
   y_index             [vector of int] The index of the endogenous variables of
                                       the block
   nze                 [integer]       number of non-zero elements in the
                                       jacobian matrix
   periods             [integer]       number of simulation periods
   y_kmin_l            [integer]       maximum number of lag in the block
   y_kmax_l            [integer]       maximum number of lead in the block
   is_linear           [integer]       if is_linear=1 the block is linear
                                       if is_linear=0 the block is not linear
   Block_Num           [integer]       block number
   y_kmin              [integer]       maximum number of lag in the model
   maxit_              [integer]       maximum number of iteration in Newton
   solve_tolf          [double]        convergence criteria
   lambda              [double]        initial value of step size in
   Newton
   cutoff              [double]        cutoff to correct the direction in Newton in case
                                       of singular jacobian matrix
   stack_solve_algo    [integer]       linear solver method used in the
                                       Newton algorithm : 
                                            - 1 sprse LU
                                            - 2 GMRES
                                            - 3 BicGStab
                                            - 4 Optimal path length

 OUTPUTS
   y                   [matrix]        All endogenous variables of the model      
  
 ALGORITHM
   Newton with LU or GMRES or BicGstab
    
 SPECIAL REQUIREMENTS
   none.

CROSS-REFERENCE INFORMATION ^

This function calls: This function is called by:

SOURCE CODE ^

0001 function y = solve_two_boundaries(fname, y, x, params, steady_state, y_index, nze, periods, y_kmin_l, y_kmax_l, is_linear, Block_Num, y_kmin, maxit_, solve_tolf, lambda, cutoff, stack_solve_algo, M, oo)
0002 % Computes the deterministic simulation of a block of equation containing
0003 % both lead and lag variables using relaxation methods
0004 %
0005 % INPUTS
0006 %   fname               [string]        name of the file containing the block
0007 %                                       to simulate
0008 %   y                   [matrix]        All the endogenous variables of the model
0009 %   x                   [matrix]        All the exogenous variables of the model
0010 %   params              [vector]        All the parameters of the model
0011 %   steady_state        [vector]        steady state of the model
0012 %   y_index             [vector of int] The index of the endogenous variables of
0013 %                                       the block
0014 %   nze                 [integer]       number of non-zero elements in the
0015 %                                       jacobian matrix
0016 %   periods             [integer]       number of simulation periods
0017 %   y_kmin_l            [integer]       maximum number of lag in the block
0018 %   y_kmax_l            [integer]       maximum number of lead in the block
0019 %   is_linear           [integer]       if is_linear=1 the block is linear
0020 %                                       if is_linear=0 the block is not linear
0021 %   Block_Num           [integer]       block number
0022 %   y_kmin              [integer]       maximum number of lag in the model
0023 %   maxit_              [integer]       maximum number of iteration in Newton
0024 %   solve_tolf          [double]        convergence criteria
0025 %   lambda              [double]        initial value of step size in
0026 %   Newton
0027 %   cutoff              [double]        cutoff to correct the direction in Newton in case
0028 %                                       of singular jacobian matrix
0029 %   stack_solve_algo    [integer]       linear solver method used in the
0030 %                                       Newton algorithm :
0031 %                                            - 1 sprse LU
0032 %                                            - 2 GMRES
0033 %                                            - 3 BicGStab
0034 %                                            - 4 Optimal path length
0035 %
0036 % OUTPUTS
0037 %   y                   [matrix]        All endogenous variables of the model
0038 %
0039 % ALGORITHM
0040 %   Newton with LU or GMRES or BicGstab
0041 %
0042 % SPECIAL REQUIREMENTS
0043 %   none.
0044 %
0045 
0046 % Copyright (C) 1996-2011 Dynare Team
0047 %
0048 % This file is part of Dynare.
0049 %
0050 % Dynare is free software: you can redistribute it and/or modify
0051 % it under the terms of the GNU General Public License as published by
0052 % the Free Software Foundation, either version 3 of the License, or
0053 % (at your option) any later version.
0054 %
0055 % Dynare is distributed in the hope that it will be useful,
0056 % but WITHOUT ANY WARRANTY; without even the implied warranty of
0057 % MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
0058 % GNU General Public License for more details.
0059 %
0060 % You should have received a copy of the GNU General Public License
0061 % along with Dynare.  If not, see <http://www.gnu.org/licenses/>.
0062 
0063 cvg=0;
0064 iter=0;
0065 Per_u_=0;
0066 g2 = [];
0067 g3 = [];
0068 Blck_size=size(y_index,2);
0069 correcting_factor=0.01;
0070 luinc_tol=1e-10;
0071 max_resa=1e100;
0072 Jacobian_Size=Blck_size*(y_kmin+y_kmax_l +periods);
0073 g1=spalloc( Blck_size*periods, Jacobian_Size, nze*periods);
0074 reduced = 0;
0075 while ~(cvg==1 || iter>maxit_),
0076     [r, y, g1, g2, g3, b]=feval(fname, y, x, params, steady_state, periods, 0, y_kmin, Blck_size);
0077     %     fjac = zeros(Blck_size, Blck_size*(y_kmin_l+1+y_kmax_l));
0078     %     disp(['Blck_size=' int2str(Blck_size) ' size(y_index)=' int2str(size(y_index,2))]);
0079     %     dh = max(abs(y(y_kmin+1-y_kmin_l:y_kmin+1+y_kmax_l, y_index)),options_.gstep*ones(y_kmin_l+1+y_kmax_l, Blck_size))*eps^(1/3);
0080     %     fvec = r;
0081     %     %for i = y_kmin+1-y_kmin_l:y_kmin+1+y_kmax_l
0082     %     i = y_kmin+1;
0083     %       i
0084     %       for j = 1:Blck_size
0085     %             ydh = y ;
0086     %           ydh(i, y_index(j)) = ydh(i, y_index(j)) + dh(i, j)  ;
0087     %           if(j==11 && i==2)
0088     %               disp(['y(i,y_index(11)=' int2str(y_index(11)) ')= ' num2str(y(i,y_index(11))) ' ydh(i, y_index(j))=' num2str(ydh(i, y_index(j))) ' dh(i,j)= ' num2str(dh(i,j))]);
0089     %           end;
0090     %           [t, y1, g11, g21, g31, b1]=feval(fname, ydh, x, params, periods, 0, y_kmin, Blck_size);
0091     %           fjac(:,j+(i-(y_kmin+1-y_kmin_l))*Blck_size) = (t(:, 1+y_kmin) - fvec(:, 1+y_kmin))./dh(i, j) ;
0092     %           if(j==11 && i==2)
0093     %                disp(['fjac(:,' int2str(j+(i-(y_kmin+1-y_kmin_l))*Blck_size) ')=']);
0094     %                disp([num2str(fjac(:,j+(i-(y_kmin+1-y_kmin_l))*Blck_size))]);
0095     %           end;
0096     %       end;
0097     % %    end
0098     %     %diff = g1(1:Blck_size, 1:Blck_size*(y_kmin_l+1+y_kmax_l)) -fjac;
0099     %     diff = g1(1:Blck_size, y_kmin_l*Blck_size+1:(y_kmin_l+1)*Blck_size) -fjac(1:Blck_size, y_kmin_l*Blck_size+1:(y_kmin_l+1)*Blck_size);
0100     %     disp(diff);
0101     %     [c_max, i_c_max] = max(abs(diff));
0102     %     [l_c_max, i_r_max] = max(c_max);
0103     %     disp(['maximum element row=' int2str(i_c_max(i_r_max)) ' and column=' int2str(i_r_max) ' value = ' num2str(l_c_max)]);
0104     %     equation = i_c_max(i_r_max);
0105     %     variable = i_r_max;
0106     %     variable
0107     %     disp(['equation ' int2str(equation) ' and variable ' int2str(y_index(mod(variable, Blck_size))) ' ' M_.endo_names(y_index(mod(variable, Blck_size)), :)]);
0108     %     disp(['g1(' int2str(equation) ', ' int2str(variable) ')=' num2str(g1(equation, y_kmin_l*Blck_size+variable),'%3.10f') ' fjac(' int2str(equation) ', ' int2str(variable) ')=' num2str(fjac(equation, y_kmin_l*Blck_size+variable), '%3.10f')]);
0109     %     return;
0110 
0111 
0112 
0113     %     for i=1:periods;
0114     %       disp([sprintf('%5.14f ',[T9025(i) T1149(i) T11905(i)])]);
0115     %     end;
0116     %     return;
0117     %residual = r(:,y_kmin+1:y_kmin+1+y_kmax_l);
0118     %num2str(residual,' %1.6f')
0119     %jac_ = g1(1:(y_kmin)*Blck_size, 1:(y_kmin+1+y_kmax_l)*Blck_size);
0120     %jac_
0121     
0122     g1a=g1(:, y_kmin*Blck_size+1:(periods+y_kmin)*Blck_size);
0123     term1 = g1(:, 1:y_kmin_l*Blck_size)*reshape(y(1+y_kmin-y_kmin_l:y_kmin,y_index)',1,y_kmin_l*Blck_size)';
0124     term2 = g1(:, (periods+y_kmin_l)*Blck_size+1:(periods+y_kmin_l+y_kmax_l)*Blck_size)*reshape(y(periods+y_kmin+1:periods+y_kmin+y_kmax_l,y_index)',1,y_kmax_l*Blck_size)';
0125     b = b - term1 - term2;
0126     
0127     %      fid = fopen(['result' num2str(iter)],'w');
0128     %      fg1a = full(g1a);
0129     %      fprintf(fid,'%d\n',size(fg1a,1));
0130     %      fprintf(fid,'%d\n',size(fg1a,2));
0131     %      fprintf(fid,'%5.14f\n',fg1a);
0132     %      fprintf(fid,'%d\n',size(b,1));
0133     %      fprintf(fid,'%5.14f\n',b);
0134     %      fclose(fid);
0135     %      return;
0136     %ipconfigb_ = b(1:(1+y_kmin)*Blck_size);
0137     %b_
0138     
0139     
0140     [max_res, max_indx]=max(max(abs(r')));
0141     if(~isreal(r))
0142         max_res = (-max_res^2)^0.5;
0143     end;
0144     %     if(~isreal(r))
0145     %       max_res=(-(max(max(abs(r))))^2)^0.5;
0146     %     else
0147     %       max_res=max(max(abs(r)));
0148     %     end;
0149     if(~isreal(max_res) || isnan(max_res))
0150         cvg = 0;
0151     elseif(is_linear && iter>0)
0152         cvg = 1;
0153     else
0154         cvg=(max_res<solve_tolf);
0155     end;
0156     if(~cvg)
0157         if(iter>0)
0158             if(~isreal(max_res) || isnan(max_res) || (max_resa<max_res && iter>1))
0159                 if(~isreal(max_res))
0160                     disp(['Variable ' M.endo_names(max_indx,:) ' (' int2str(max_indx) ') returns an undefined value']);
0161                 end;
0162                 if(isnan(max_res))
0163                     detJ=det(g1aa);
0164                     if(abs(detJ)<1e-7)
0165                         max_factor=max(max(abs(g1aa)));
0166                         ze_elem=sum(diag(g1aa)<cutoff);
0167                         disp([num2str(full(ze_elem),'%d') ' elements on the Jacobian diagonal are below the cutoff (' num2str(cutoff,'%f') ')']);
0168                         if(correcting_factor<max_factor)
0169                             correcting_factor=correcting_factor*4;
0170                             disp(['The Jacobain matrix is singular, det(Jacobian)=' num2str(detJ,'%f') '.']);
0171                             disp(['    trying to correct the Jacobian matrix:']);
0172                             disp(['    correcting_factor=' num2str(correcting_factor,'%f') ' max(Jacobian)=' num2str(full(max_factor),'%f')]);
0173                             dx = (g1aa+correcting_factor*speye(periods*Blck_size))\ba- ya;
0174                             y(1+y_kmin:periods+y_kmin,y_index)=reshape((ya_save+lambda*dx)',length(y_index),periods)';
0175                             continue;
0176                         else
0177                             disp('The singularity of the jacobian matrix could not be corrected');
0178                             return;
0179                         end;
0180                     end;
0181                 elseif(lambda>1e-8)
0182                     lambda=lambda/2;
0183                     reduced = 1;
0184                     disp(['reducing the path length: lambda=' num2str(lambda,'%f')]);
0185                     y(1+y_kmin:periods+y_kmin,y_index)=reshape((ya_save+lambda*dx)',length(y_index),periods)';
0186                     continue;
0187                 else
0188                     if(cutoff == 0)
0189                         fprintf('Error in simul: Convergence not achieved in block %d, after %d iterations.\n Increase "options_.maxit_".\n',Block_Num, iter);
0190                     else
0191                         fprintf('Error in simul: Convergence not achieved in block %d, after %d iterations.\n Increase "options_.maxit_" or set "cutoff=0" in model options.\n',Block_Num, iter);
0192                     end;
0193                     oo.deterministic_simulation.status = 0;
0194                     oo.deterministic_simulation.error = max_res;
0195                     oo.deterministic_simulation.iterations = iter;
0196                     oo.deterministic_simulation.block(Block_Num).status = 0;% Convergency failed.
0197                     oo.deterministic_simulation.block(Block_Num).error = max_res;
0198                     oo.deterministic_simulation.block(Block_Num).iterations = iter;
0199                     return;
0200                 end;
0201             else
0202                 if(lambda<1)
0203                     lambda=max(lambda*2, 1);
0204                 end;
0205             end;
0206         end;
0207         ya = reshape(y(y_kmin+1:y_kmin+periods,y_index)',1,periods*Blck_size)';
0208         ya_save=ya;
0209         g1aa=g1a;
0210         ba=b;
0211         max_resa=max_res;
0212         if(stack_solve_algo==0),
0213             dx = g1a\b- ya;
0214             ya = ya + lambda*dx;
0215             y(1+y_kmin:periods+y_kmin,y_index)=reshape(ya',length(y_index),periods)';
0216         elseif(stack_solve_algo==1),
0217             for t=1:periods;
0218                 first_elem = (t-1)*Blck_size+1;
0219                 last_elem = t*Blck_size;
0220                 next_elem = (t+1)*Blck_size;
0221                 Elem = first_elem:last_elem;
0222                 Elem_1 = last_elem+1:next_elem;
0223                 B1_inv = inv(g1a(Elem, Elem));
0224                 if (t < periods)
0225                     S1 = B1_inv * g1a(Elem, Elem_1);
0226                 end;
0227                 g1a(Elem, Elem_1) = S1;
0228                 b(Elem) = B1_inv * b(Elem);
0229                 g1a(Elem, Elem) = ones(Blck_size, Blck_size);
0230                 if (t < periods)
0231                     g1a(Elem_1, Elem_1) = g1a(Elem_1, Elem_1) - g1a(Elem_1, Elem) * S1;
0232                     b(Elem_1) = b(Elem_1) - g1a(Elem_1, Elem) * b(Elem);
0233                     g1a(Elem_1, Elem) = zeros(Blck_size, Blck_size);
0234                 end;
0235             end;
0236             za = b(Elem);
0237             zaa = za;
0238             y_Elem = (periods - 1) * Blck_size + 1:(periods) * Blck_size;
0239             dx = ya;
0240             dx(y_Elem) = za - ya(y_Elem);
0241             ya(y_Elem) = ya(y_Elem) + lambda*dx(y_Elem);
0242             for t=periods-1:-1:1;
0243                 first_elem = (t-1)*Blck_size+1;
0244                 last_elem = t*Blck_size;
0245                 next_elem = (t+1)*Blck_size;
0246                 Elem_1 = last_elem+1:next_elem;
0247                 Elem = first_elem:last_elem;
0248                 za = b(Elem) - g1a(Elem, Elem_1) * zaa;
0249                 zaa = za;
0250                 y_Elem = Blck_size * (t-1)+1:Blck_size * (t);
0251                 dx(y_Elem) = za - ya(y_Elem);
0252                 ya(y_Elem) = ya(y_Elem) + lambda*dx(y_Elem);
0253                 y(y_kmin + t, y_index) = ya(y_Elem);
0254             end;
0255         elseif(stack_solve_algo==2),
0256             flag1=1;
0257             while(flag1>0)
0258                 [L1, U1]=luinc(g1a,luinc_tol);
0259                 [za,flag1] = gmres(g1a,b,Blck_size,1e-6,Blck_size*periods,L1,U1);
0260                 if (flag1>0 || reduced)
0261                     if(flag1==1)
0262                         disp(['Error in simul: No convergence inside GMRES after ' num2str(periods*10,'%6d') ' iterations, in block' num2str(Block_Size,'%3d')]);
0263                     elseif(flag1==2)
0264                         disp(['Error in simul: Preconditioner is ill-conditioned, in block' num2str(Block_Size,'%3d')]);
0265                     elseif(flag1==3)
0266                         disp(['Error in simul: GMRES stagnated (Two consecutive iterates were the same.), in block' num2str(Block_Size,'%3d')]);
0267                     end;
0268                     luinc_tol = luinc_tol/10;
0269                     reduced = 0;
0270                 else
0271                     dx = za - ya;
0272                     ya = ya + lambda*dx;
0273                     y(1+y_kmin:periods+y_kmin,y_index)=reshape(ya',length(y_index),periods)';
0274                 end;
0275             end;
0276         elseif(stack_solve_algo==3),
0277             flag1=1;
0278             while(flag1>0)
0279                 [L1, U1]=luinc(g1a,luinc_tol);
0280                 [za,flag1] = bicgstab(g1a,b,1e-7,Blck_size*periods,L1,U1);
0281                 if (flag1>0 || reduced)
0282                     if(flag1==1)
0283                         disp(['Error in simul: No convergence inside BICGSTAB after ' num2str(periods*10,'%6d') ' iterations, in block' num2str(Block_Size,'%3d')]);
0284                     elseif(flag1==2)
0285                         disp(['Error in simul: Preconditioner is ill-conditioned, in block' num2str(Block_Size,'%3d')]);
0286                     elseif(flag1==3)
0287                         disp(['Error in simul: GMRES stagnated (Two consecutive iterates were the same.), in block' num2str(Block_Size,'%3d')]);
0288                     end;
0289                     luinc_tol = luinc_tol/10;
0290                     reduced = 0;
0291                 else
0292                     dx = za - ya;
0293                     ya = ya + lambda*dx;
0294                     y(1+y_kmin:periods+y_kmin,y_index)=reshape(ya',length(y_index),periods)';
0295                 end;
0296             end;
0297         elseif(stack_solve_algo==4),
0298             ra = reshape(r(:, y_kmin+1:periods+y_kmin),periods*Blck_size, 1);
0299             stpmx = 100 ;
0300             stpmax = stpmx*max([sqrt(ya'*ya);size(y_index,2)]);
0301             nn=1:size(ra,1);
0302             g = (ra'*g1a)';
0303             f = 0.5*ra'*ra;
0304             p = -g1a\ra;
0305             [yn,f,ra,check]=lnsrch1(ya,f,g,p,stpmax, ...
0306                                     'lnsrch1_wrapper_two_boundaries',nn,nn,  fname, y, y_index, x, params, steady_state, periods, y_kmin, Blck_size);
0307             dx = ya - yn;
0308             y(1+y_kmin:periods+y_kmin,y_index)=reshape(yn',length(y_index),periods)';
0309         end
0310     end
0311     iter=iter+1;
0312     disp(['iteration: ' num2str(iter,'%d') ' error: ' num2str(max_res,'%e')]);
0313 end;
0314 if (iter>maxit_)
0315     disp(['No convergence after ' num2str(iter,'%4d') ' iterations in Block ' num2str(Block_Num,'%d')]);
0316     oo.deterministic_simulation.status = 0;
0317     oo.deterministic_simulation.error = max_res;
0318     oo.deterministic_simulation.iterations = iter;
0319     oo.deterministic_simulation.block(Block_Num).status = 0;% Convergency failed.
0320     oo.deterministic_simulation.block(Block_Num).error = max_res;
0321     oo.deterministic_simulation.block(Block_Num).iterations = iter;
0322     return;
0323 end
0324 oo.deterministic_simulation.status = 1;
0325 oo.deterministic_simulation.error = max_res;
0326 oo.deterministic_simulation.iterations = iter;
0327 oo.deterministic_simulation.block(Block_Num).status = 1;% Convergency obtained.
0328 oo.deterministic_simulation.block(Block_Num).error = max_res;
0329 oo.deterministic_simulation.block(Block_Num).iterations = iter;
0330 return;

Generated on Tue 22-May-2012 02:40:23 by m2html © 2005