function x0 = stab_map_(OutputDirectoryName,opt_gsa) % % function x0 = stab_map_(OutputDirectoryName) % % Mapping of stability regions in the prior ranges applying % Monte Carlo filtering techniques. % % INPUTS (from opt_gsa structure) % Nsam = MC sample size % fload = 0 to run new MC; 1 to load prevoiusly generated analysis % alpha2 = significance level for bivariate sensitivity analysis % [abs(corrcoef) > alpha2] % prepSA = 1: save transition matrices for mapping reduced form % = 0: no transition matrix saved (default) % pprior = 1: sample from prior ranges (default): sample saved in % _prior.mat file % = 0: sample from posterior ranges: sample saved in % _mc.mat file % OUTPUT: % x0: one parameter vector for which the model is stable. % % GRAPHS % 1) Pdf's of marginal distributions under the stability (dotted % lines) and unstability (solid lines) regions % 2) Cumulative distributions of: % - stable subset (dotted lines) % - unacceptable subset (solid lines) % 3) Bivariate plots of significant correlation patterns % ( abs(corrcoef) > alpha2) under the stable and unacceptable subsets % % USES qmc_sequence, stab_map_1, stab_map_2 % % Written by Marco Ratto % Joint Research Centre, The European Commission, % (http://eemc.jrc.ec.europa.eu/), % marco.ratto@jrc.it % % Reference: % M. Ratto, Global Sensitivity Analysis for Macroeconomic models, MIMEO, 2006. % Copyright (C) 2012-2013 Dynare Team % % This file is part of Dynare. % % Dynare is free software: you can redistribute it and/or modify % it under the terms of the GNU General Public License as published by % the Free Software Foundation, either version 3 of the License, or % (at your option) any later version. % % Dynare is distributed in the hope that it will be useful, % but WITHOUT ANY WARRANTY; without even the implied warranty of % MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the % GNU General Public License for more details. % % You should have received a copy of the GNU General Public License % along with Dynare. If not, see . %global bayestopt_ estim_params_ dr_ options_ ys_ fname_ global bayestopt_ estim_params_ options_ oo_ M_ % opt_gsa=options_.opt_gsa; Nsam = opt_gsa.Nsam; fload = opt_gsa.load_stab; ksstat = opt_gsa.ksstat; alpha2 = opt_gsa.alpha2_stab; pvalue_ks = opt_gsa.pvalue_ks; pvalue_corr = opt_gsa.pvalue_corr; prepSA = (opt_gsa.redform | opt_gsa.identification); pprior = opt_gsa.pprior; neighborhood_width = opt_gsa.neighborhood_width; ilptau = opt_gsa.ilptau; nliv = opt_gsa.morris_nliv; ntra = opt_gsa.morris_ntra; dr_ = oo_.dr; %if isfield(dr_,'ghx'), ys_ = oo_.dr.ys; nspred = M_.nspred; %size(dr_.ghx,2); nboth = M_.nboth; nfwrd = M_.nfwrd; %end fname_ = M_.fname; np = estim_params_.np; nshock = estim_params_.nvx; nshock = nshock + estim_params_.nvn; nshock = nshock + estim_params_.ncx; nshock = nshock + estim_params_.ncn; lpmat0=[]; xparam1=[]; pshape = bayestopt_.pshape(nshock+1:end); p1 = bayestopt_.p1(nshock+1:end); p2 = bayestopt_.p2(nshock+1:end); p3 = bayestopt_.p3(nshock+1:end); p4 = bayestopt_.p4(nshock+1:end); if nargin==0, OutputDirectoryName=''; end opt=options_; options_.periods=0; options_.nomoments=1; options_.irf=0; options_.noprint=1; options_.simul=0; if fload==0, % if prepSA % T=zeros(size(dr_.ghx,1),size(dr_.ghx,2)+size(dr_.ghu,2),Nsam/2); % end if isfield(dr_,'ghx'), egg=zeros(length(dr_.eigval),Nsam); end yys=zeros(length(dr_.ys),Nsam); if opt_gsa.morris == 1 [lpmat, OutFact] = Sampling_Function_2(nliv, np+nshock, ntra, ones(np+nshock, 1), zeros(np+nshock,1), []); lpmat = lpmat.*(nliv-1)/nliv+1/nliv/2; Nsam=size(lpmat,1); lpmat0 = lpmat(:,1:nshock); lpmat = lpmat(:,nshock+1:end); % elseif opt_gsa.morris==3, % lpmat = prep_ide(Nsam,np,5); % Nsam=size(lpmat,1); else if np<52 && ilptau>0, [lpmat] = qmc_sequence(np, int64(1), 0, Nsam)'; if np>30 || ilptau==2, % scrambled lptau for j=1:np, lpmat(:,j)=lpmat(randperm(Nsam),j); end end else %ilptau==0 %[lpmat] = rand(Nsam,np); for j=1:np, lpmat(:,j) = randperm(Nsam)'./(Nsam+1); %latin hypercube end end end % try dummy=prior_draw_gsa(1); % catch % if pprior, % if opt_gsa.prior_range==0; % error('Some unknown prior is specified or ML estimation,: use prior_range=1 option!!'); % end % end % % end if pprior, for j=1:nshock, if opt_gsa.morris~=1, lpmat0(:,j) = randperm(Nsam)'./(Nsam+1); %latin hypercube end if opt_gsa.prior_range lpmat0(:,j)=lpmat0(:,j).*(bayestopt_.ub(j)-bayestopt_.lb(j))+bayestopt_.lb(j); end end if opt_gsa.prior_range % if opt_gsa.identification, % deltx=min(0.001, 1/Nsam/2); % for j=1:np, % xdelt(:,:,j)=prior_draw_gsa(0,[lpmat0 lpmat]+deltx); % end % end for j=1:np, lpmat(:,j)=lpmat(:,j).*(bayestopt_.ub(j+nshock)-bayestopt_.lb(j+nshock))+bayestopt_.lb(j+nshock); end else xx=prior_draw_gsa(0,[lpmat0 lpmat]); % if opt_gsa.identification, % deltx=min(0.001, 1/Nsam/2); % ldum=[lpmat0 lpmat]; % ldum = prior_draw_gsa(0,ldum+deltx); % for j=1:nshock+np, % xdelt(:,:,j)=xx; % xdelt(:,j,j)=ldum(:,j); % end % clear ldum % end lpmat0=xx(:,1:nshock); lpmat=xx(:,nshock+1:end); clear xx; end else % for j=1:nshock, % xparam1(j) = oo_.posterior_mode.shocks_std.(bayestopt_.name{j}); % sd(j) = oo_.posterior_std.shocks_std.(bayestopt_.name{j}); % lpmat0(:,j) = randperm(Nsam)'./(Nsam+1); %latin hypercube % lb = max(bayestopt_.lb(j), xparam1(j)-2*sd(j)); % ub1=xparam1(j)+(xparam1(j) - lb); % define symmetric range around the mode! % ub = min(bayestopt_.ub(j),ub1); % if ub30 & np<52 % lpmat(:,j) = lpmat(randperm(Nsam),j).*(ub-lb)+lb; % else % lpmat(:,j) = lpmat(:,j).*(ub-lb)+lb; % end % end %load([fname_,'_mode']) eval(['load ' options_.mode_file '.mat;']); if neighborhood_width>0, for j=1:nshock, lpmat0(:,j) = randperm(Nsam)'./(Nsam+1); %latin hypercube ub=min([bayestopt_.ub(j) xparam1(j)*(1+neighborhood_width)]); lb=max([bayestopt_.lb(j) xparam1(j)*(1-neighborhood_width)]); lpmat0(:,j)=lpmat0(:,j).*(ub-lb)+lb; end for j=1:np, ub=min([bayestopt_.ub(j+nshock) xparam1(j+nshock)*(1+neighborhood_width)]); lb=max([bayestopt_.lb(j+nshock) xparam1(j+nshock)*(1-neighborhood_width)]); lpmat(:,j)=lpmat(:,j).*(ub-lb)+lb; end else d = chol(inv(hh)); lp=randn(Nsam*2,nshock+np)*d+kron(ones(Nsam*2,1),xparam1'); for j=1:Nsam*2, lnprior(j) = any(lp(j,:)'<=bayestopt_.lb | lp(j,:)'>=bayestopt_.ub); end ireal=[1:2*Nsam]; ireal=ireal(find(lnprior==0)); lp=lp(ireal,:); Nsam=min(Nsam, length(ireal)); lpmat0=lp(1:Nsam,1:nshock); lpmat=lp(1:Nsam,nshock+1:end); clear lp lnprior ireal; end end % h = dyn_waitbar(0,'Please wait...'); istable=[1:Nsam]; jstab=0; iunstable=[1:Nsam]; iindeterm=zeros(1,Nsam); iwrong=zeros(1,Nsam); for j=1:Nsam, M_.params(estim_params_.param_vals(:,1)) = lpmat(j,:)'; %try stoch_simul([]); try [Tt,Rr,SteadyState,info,M_,options_,oo_] = dynare_resolve(M_,options_,oo_,'restrict'); infox(j,1)=info(1); if infox(j,1)==0 && ~exist('T'), dr_=oo_.dr; if prepSA, try T=zeros(size(dr_.ghx,1),size(dr_.ghx,2)+size(dr_.ghu,2),Nsam); catch ME if strcmp('MATLAB:nomem',ME.identifier), prepSA=0; disp('The model is too large for storing state space matrices ...') disp('for mapping reduced form or for identification') end T=[]; end else T=[]; end egg=zeros(length(dr_.eigval),Nsam); end if infox(j,1), % disp('no solution'), if isfield(oo_.dr,'ghx'), oo_.dr=rmfield(oo_.dr,'ghx'); end if (infox(j,1)<3 || infox(j,1)>5) && isfield(oo_.dr,'eigval'), oo_.dr=rmfield(oo_.dr,'eigval'); end end catch ME if isfield(oo_.dr,'eigval'), oo_.dr=rmfield(oo_.dr,'eigval'); end if isfield(oo_.dr,'ghx'), oo_.dr=rmfield(oo_.dr,'ghx'); end disp('No solution could be found'), end dr_ = oo_.dr; if isfield(dr_,'ghx'), egg(:,j) = sort(dr_.eigval); iunstable(j)=0; if prepSA jstab=jstab+1; T(:,:,jstab) = [dr_.ghx dr_.ghu]; % [A,B] = ghx2transition(squeeze(T(:,:,jstab)), ... % bayestopt_.restrict_var_list, ... % bayestopt_.restrict_columns, ... % bayestopt_.restrict_aux); end if ~exist('nspred'), nspred = dr_.nspred; %size(dr_.ghx,2); nboth = dr_.nboth; nfwrd = dr_.nfwrd; end else istable(j)=0; if isfield(dr_,'eigval') egg(:,j) = sort(dr_.eigval); if exist('nspred') if any(isnan(egg(1:nspred,j))) iwrong(j)=j; else if (nboth || nfwrd) && abs(egg(nspred+1,j))<=options_.qz_criterium, iindeterm(j)=j; end end end else if exist('egg'), egg(:,j)=ones(size(egg,1),1).*NaN; end iwrong(j)=j; end end ys_=real(dr_.ys); yys(:,j) = ys_; ys_=yys(:,1); dyn_waitbar(j/Nsam,h,['MC iteration ',int2str(j),'/',int2str(Nsam)]) end dyn_waitbar_close(h); if prepSA && jstab, T=T(:,:,1:jstab); else T=[]; end istable=istable(find(istable)); % stable params iunstable=iunstable(find(iunstable)); % unstable params iindeterm=iindeterm(find(iindeterm)); % indeterminacy iwrong=iwrong(find(iwrong)); % dynare could not find solution % % map stable samples % istable=[1:Nsam]; % for j=1:Nsam, % if any(isnan(egg(1:nspred,j))) % istable(j)=0; % else % if abs(egg(nspred,j))>=options_.qz_criterium; %(1-(options_.qz_criterium-1)); %1-1.e-5; % istable(j)=0; % %elseif (dr_.nboth | dr_.nfwrd) & abs(egg(nspred+1,j))<=options_.qz_criterium; %1+1.e-5; % elseif (nboth | nfwrd) & abs(egg(nspred+1,j))<=options_.qz_criterium; %1+1.e-5; % istable(j)=0; % end % end % end % istable=istable(find(istable)); % stable params % % % map unstable samples % iunstable=[1:Nsam]; % for j=1:Nsam, % %if abs(egg(dr_.npred+1,j))>1+1.e-5 & abs(egg(dr_.npred,j))<1-1.e-5; % %if (dr_.nboth | dr_.nfwrd), % if ~any(isnan(egg(1:5,j))) % if (nboth | nfwrd), % if abs(egg(nspred+1,j))>options_.qz_criterium & abs(egg(nspred,j))0 && length(iunstable)ksstat); indstab=find(probaksstat); indindet=find(probaksstat); indunst=find(probaksstat); indwrong=find(proba10, stab_map_2(lpmat(iunstable,:),alpha2, pvalue_corr, auname, OutputDirectoryName,xparam1,autitle); end if length(iindeterm)>10, stab_map_2(lpmat(iindeterm,:),alpha2, pvalue_corr, aindname, OutputDirectoryName,xparam1,aindtitle); end if length(ixun)>10, stab_map_2(lpmat(ixun,:),alpha2, pvalue_corr, aunstname, OutputDirectoryName,xparam1,aunsttitle); end if length(iwrong)>10, stab_map_2(lpmat(iwrong,:),alpha2, pvalue_corr, awrongname, OutputDirectoryName,xparam1,awrongtitle); end x0=0.5.*(bayestopt_.ub(1:nshock)-bayestopt_.lb(1:nshock))+bayestopt_.lb(1:nshock); x0 = [x0; lpmat(istable(1),:)']; if istable(end)~=Nsam M_.params(estim_params_.param_vals(:,1)) = lpmat(istable(1),:)'; [oo_.dr,info,M_,options_,oo_] = resol(0,M_,options_,oo_); % stoch_simul([]); end else if length(iunstable)==0, disp('All parameter values in the specified ranges give unique saddle-path solution!') x0=0.5.*(bayestopt_.ub(1:nshock)-bayestopt_.lb(1:nshock))+bayestopt_.lb(1:nshock); x0 = [x0; lpmat(istable(1),:)']; else disp('All parameter values in the specified ranges are not acceptable!') x0=[]; end end xparam1=x0; save prior_ok xparam1; options_.periods=opt.periods; if isfield(opt,'nomoments'), options_.nomoments=opt.nomoments; end options_.irf=opt.irf; options_.noprint=opt.noprint; if isfield(opt,'simul'), options_.simul=opt.simul; end