dynare/matlab/estimation/dsge_likelihood.m

913 lines
36 KiB
Matlab
Raw Normal View History

function [fval,info,exit_flag,DLIK,Hess,SteadyState,trend_coeff,M_,options_,bayestopt_,dr] = dsge_likelihood(xparam1,dataset_,dataset_info,options_,M_,estim_params_,bayestopt_,BoundsInfo,dr, endo_steady_state, exo_steady_state, exo_det_steady_state,derivatives_info)
% [fval,info,exit_flag,DLIK,Hess,SteadyState,trend_coeff,M_,options_,bayestopt_,oo_] = dsge_likelihood(xparam1,dataset_,dataset_info,options_,M_,estim_params_,bayestopt_,BoundsInfo,oo_,derivatives_info)
% Evaluates the posterior kernel of a DSGE model using the specified
% kalman_algo; the resulting posterior includes the 2*pi constant of the
% likelihood function
%
% INPUTS
% - xparam1 [double] current values for the estimated parameters.
% - dataset_ [structure] dataset after transformations
% - dataset_info [structure] storing informations about the
% sample; not used but required for interface
% - options_ [structure] Matlab's structure describing the current options
% - M_ [structure] Matlab's structure describing the model
% - estim_params_ [structure] characterizing parameters to be estimated
% - bayestopt_ [structure] describing the priors
% - BoundsInfo [structure] containing prior bounds
% - dr [structure] Reduced form model.
% - endo_steady_state [vector] steady state value for endogenous variables
% - exo_steady_state [vector] steady state value for exogenous variables
% - exo_det_steady_state [vector] steady state value for exogenous deterministic variables
% - derivatives_info [structure] derivative info for identification
%
% OUTPUTS
% - fval [double] scalar, value of the likelihood or posterior kernel.
% - info [integer] 4×1 vector, informations resolution of the model and evaluation of the likelihood.
% - exit_flag [integer] scalar, equal to 1 (no issues when evaluating the likelihood) or 0 (not able to evaluate the likelihood).
% - DLIK [double] Vector with score of the likelihood
% - Hess [double] asymptotic hessian matrix.
% - SteadyState [double] steady state level for the endogenous variables
% - trend_coeff [double] Matrix of doubles, coefficients of the deterministic trend in the measurement equation.
% - M_ [struct] Updated M_ structure described in INPUTS section.
% - options_ [struct] Updated options_ structure described in INPUTS section.
% - bayestopt_ [struct] See INPUTS section.
% - dr [structure] Reduced form model.
%
% This function is called by: dynare_estimation_1, mode_check
% This function calls: dynare_resolve, lyapunov_symm, lyapunov_solver, compute_Pinf_Pstar, kalman_filter_d, missing_observations_kalman_filter_d,
% univariate_kalman_filter_d, kalman_steady_state, get_perturbation_params_deriv, kalman_filter, missing_observations_kalman_filter, univariate_kalman_filter, priordens
% Copyright © 2004-2023 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 <https://www.gnu.org/licenses/>.
2021-12-13 12:42:21 +01:00
% Initial author: stephane DOT adjemian AT univ DASH lemans DOT FR
% Initialization of the returned variables and others...
fval = [];
SteadyState = [];
trend_coeff = [];
exit_flag = 1;
2016-06-01 18:22:51 +02:00
info = zeros(4,1);
if options_.analytic_derivation
DLIK = NaN(1,length(xparam1));
else
DLIK = [];
end
Hess = [];
% Ensure that xparam1 is a column vector.
% (Don't do the transformation if xparam1 is empty, otherwise it would become a
% 0×1 matrix, which create issues with older MATLABs when comparing with [] in
% check_bounds_and_definiteness_estimation)
if ~isempty(xparam1)
xparam1 = xparam1(:);
end
% Set flag related to analytical derivatives.
analytic_derivation = options_.analytic_derivation;
if analytic_derivation
if options_.loglinear
error('The analytic_derivation and loglinear options are not compatible')
end
if options_.endogenous_prior
error('The analytic_derivation and endogenous_prior options are not compatible')
end
end
2017-05-16 12:42:01 +02:00
if nargout==1
analytic_derivation=0;
end
2017-05-16 12:42:01 +02:00
if analytic_derivation
kron_flag=options_.analytic_derivation_mode;
end
%------------------------------------------------------------------------------
% 1. Get the structural parameters & define penalties
%------------------------------------------------------------------------------
M_ = set_all_parameters(xparam1,estim_params_,M_);
[fval,info,exit_flag,Q,H]=check_bounds_and_definiteness_estimation(xparam1, M_, estim_params_, BoundsInfo);
if info(1)
return
end
%------------------------------------------------------------------------------
% 2. call model setup & reduction program
%------------------------------------------------------------------------------
is_restrict_state_space = true;
if options_.occbin.likelihood.status
occbin_options = set_occbin_options(options_);
if occbin_options.opts_simul.restrict_state_space
[T,R,SteadyState,info,dr, M_.params,TTx,RRx,CCx, T0, R0] = ...
occbin.dynare_resolve(M_,options_,dr, endo_steady_state, exo_steady_state, exo_det_steady_state,[],'restrict');
else
is_restrict_state_space = false;
oldoo.restrict_var_list = dr.restrict_var_list;
oldoo.restrict_columns = dr.restrict_columns;
dr.restrict_var_list = bayestopt_.smoother_var_list;
dr.restrict_columns = bayestopt_.smoother_restrict_columns;
% Linearize the model around the deterministic steady state and extract the matrices of the state equation (T and R).
[T,R,SteadyState,info,M_,dr, M_.params,TTx,RRx,CCx, T0, R0] = ...
occbin.dynare_resolve(M_,options_,dr, endo_steady_state, exo_steady_state, exo_det_steady_state);
dr.restrict_var_list = oldoo.restrict_var_list;
dr.restrict_columns = oldoo.restrict_columns;
end
occbin_.status = true;
occbin_.info= {options_, dr,endo_steady_state,exo_steady_state,exo_det_steady_state, M_, occbin_options, TTx, RRx, CCx,T0,R0};
else
% Linearize the model around the deterministic steady state and extract the matrices of the state equation (T and R).
[T,R,SteadyState,info,dr, M_.params] = dynare_resolve(M_,options_,dr, endo_steady_state, exo_steady_state, exo_det_steady_state,'restrict');
occbin_.status = false;
end
% Return, with endogenous penalty when possible, if dynare_resolve issues an error code (defined in resol).
2016-06-01 18:22:51 +02:00
if info(1)
if info(1) == 3 || info(1) == 4 || info(1) == 5 || info(1)==6 ||info(1) == 19 ||...
2017-05-16 15:10:20 +02:00
info(1) == 20 || info(1) == 21 || info(1) == 23 || info(1) == 26 || ...
info(1) == 81 || info(1) == 84 || info(1) == 85 || info(1) == 86 || ...
info(1) == 401 || info(1) == 402 || info(1) == 403 || ... %cycle reduction
info(1) == 411 || info(1) == 412 || info(1) == 413 % logarithmic reduction
2016-06-01 18:22:51 +02:00
%meaningful second entry of output that can be used
fval = Inf;
if isnan(info(2))
info(4) = 0.1;
else
info(4) = info(2);
end
2016-06-01 18:22:51 +02:00
exit_flag = 0;
2017-05-16 12:42:01 +02:00
if analytic_derivation
2016-06-01 18:22:51 +02:00
DLIK=ones(length(xparam1),1);
end
return
else
fval = Inf;
info(4) = 0.1;
exit_flag = 0;
2017-05-16 12:42:01 +02:00
if analytic_derivation
2016-06-01 18:22:51 +02:00
DLIK=ones(length(xparam1),1);
end
return
end
end
% check endogenous prior restrictions
info=endogenous_prior_restrictions(T,R,M_,options_,dr,endo_steady_state,exo_steady_state,exo_det_steady_state);
2017-05-16 12:42:01 +02:00
if info(1)
2016-06-01 18:22:51 +02:00
fval = Inf;
info(4)=info(2);
exit_flag = 0;
2017-05-16 12:42:01 +02:00
if analytic_derivation
DLIK=ones(length(xparam1),1);
end
return
end
if is_restrict_state_space
%% Define a vector of indices for the observed variables. Is this really usefull?...
bayestopt_.mf = bayestopt_.mf1;
else
%get location of observed variables and requested smoothed variables in
%decision rules
bayestopt_.mf = bayestopt_.smoother_var_list(bayestopt_.smoother_mf);
end
% Define the constant vector of the measurement equation.
if options_.noconstant
constant = zeros(dataset_.vobs,1);
else
if options_.loglinear
constant = log(SteadyState(bayestopt_.mfys));
else
constant = SteadyState(bayestopt_.mfys);
end
end
% Define the deterministic linear trend of the measurement equation.
if bayestopt_.with_trend
[trend_addition, trend_coeff]=compute_trend_coefficients(M_,options_,dataset_.vobs,dataset_.nobs);
trend = repmat(constant,1,dataset_.nobs)+trend_addition;
else
trend_coeff = zeros(dataset_.vobs,1);
trend = repmat(constant,1,dataset_.nobs);
end
% Get needed informations for kalman filter routines.
start = options_.presample+1;
Z = bayestopt_.mf; %selector for observed variables
no_missing_data_flag = ~dataset_info.missing.state;
mm = length(T); %number of states
pp = dataset_.vobs; %number of observables
rr = length(Q); %number of shocks
kalman_tol = options_.kalman_tol;
diffuse_kalman_tol = options_.diffuse_kalman_tol;
riccati_tol = options_.riccati_tol;
Y = transpose(dataset_.data)-trend;
smpl = size(Y,2);
%------------------------------------------------------------------------------
% 3. Initial condition of the Kalman filter
%------------------------------------------------------------------------------
kalman_algo = options_.kalman_algo;
diffuse_periods = 0;
expanded_state_vector_for_univariate_filter=0;
singular_diffuse_filter = 0;
if options_.heteroskedastic_filter
Qvec=get_Qvec_heteroskedastic_filter(Q,smpl,M_);
end
switch options_.lik_init
case 1% Standard initialization with the steady state of the state equation.
if kalman_algo~=2
% Use standard kalman filter except if the univariate filter is explicitely choosen.
kalman_algo = 1;
end
Pstar=lyapunov_solver(T,R,Q,options_);
Pinf = [];
a = zeros(mm,1);
a=set_Kalman_starting_values(a,M_,dr,options_,bayestopt_);
a_0_given_tm1=T*a; %set state prediction for first Kalman step;
if options_.occbin.likelihood.status
Z =zeros(length(bayestopt_.mf),size(T,1));
for i = 1:length(bayestopt_.mf)
Z(i,bayestopt_.mf(i))=1;
end
Zflag = 1;
else
Zflag = 0;
end
case 2% Initialization with large numbers on the diagonal of the covariance matrix if the states (for non stationary models).
if kalman_algo ~= 2
% Use standard kalman filter except if the univariate filter is explicitely choosen.
kalman_algo = 1;
end
Pstar = options_.Harvey_scale_factor*eye(mm);
Pinf = [];
a = zeros(mm,1);
a = set_Kalman_starting_values(a,M_,dr,options_,bayestopt_);
a_0_given_tm1 = T*a; %set state prediction for first Kalman step;
if options_.occbin.likelihood.status
Z =zeros(length(bayestopt_.mf),size(T,1));
for i = 1:length(bayestopt_.mf)
Z(i,bayestopt_.mf(i))=1;
end
Zflag = 1;
else
Zflag = 0;
end
case 3% Diffuse Kalman filter (Durbin and Koopman)
% Use standard kalman filter except if the univariate filter is explicitely choosen.
if kalman_algo == 0
kalman_algo = 3;
elseif ~((kalman_algo == 3) || (kalman_algo == 4))
2017-05-16 15:10:20 +02:00
error(['The model requires Diffuse filter, but you specified a different Kalman filter. You must set options_.kalman_algo ' ...
'to 0 (default), 3 or 4'])
end
[Pstar,Pinf] = compute_Pinf_Pstar(Z,T,R,Q,options_.qz_criterium);
Z =zeros(length(bayestopt_.mf),size(T,1));
for i = 1:length(bayestopt_.mf)
Z(i,bayestopt_.mf(i))=1;
end
Zflag = 1;
if options_.heteroskedastic_filter
QQ=Qvec;
else
QQ=Q;
end
% Run diffuse kalman filter on first periods.
if (kalman_algo==3)
% Multivariate Diffuse Kalman Filter
a = zeros(mm,1);
a = set_Kalman_starting_values(a,M_,dr,options_,bayestopt_);
a_0_given_tm1 = T*a; %set state prediction for first Kalman step;
Pstar0 = Pstar; % store Pstar
if no_missing_data_flag
[dLIK,dlik,a_0_given_tm1,Pstar] = kalman_filter_d(Y, 1, size(Y,2), ...
a_0_given_tm1, Pinf, Pstar, ...
kalman_tol, diffuse_kalman_tol, riccati_tol, options_.presample, ...
T,R,QQ,H,Z,mm,pp,rr);
else
[dLIK,dlik,a_0_given_tm1,Pstar] = missing_observations_kalman_filter_d(dataset_info.missing.aindex,dataset_info.missing.number_of_observations,dataset_info.missing.no_more_missing_observations, ...
Y, 1, size(Y,2), ...
a_0_given_tm1, Pinf, Pstar, ...
kalman_tol, diffuse_kalman_tol, riccati_tol, options_.presample, ...
T,R,QQ,H,Z,mm,pp,rr);
end
diffuse_periods = length(dlik);
if isinf(dLIK)
% Go to univariate diffuse filter if singularity problem.
singular_diffuse_filter = 1;
Pstar = Pstar0;
end
end
if singular_diffuse_filter || (kalman_algo==4)
% Univariate Diffuse Kalman Filter
if isequal(H,0)
H1 = zeros(pp,1);
mmm = mm;
else
if all(all(abs(H-diag(diag(H)))<1e-14))% ie, the covariance matrix is diagonal...
H1 = diag(H);
mmm = mm;
else
%Augment state vector (follows Section 6.4.3 of DK (2012))
expanded_state_vector_for_univariate_filter=1;
if Zflag
Z1=Z;
else
Z1=zeros(pp,size(T,2));
for jz=1:length(Z)
Z1(jz,Z(jz))=1;
end
end
Z = [Z1, eye(pp)];
Zflag=1;
T = blkdiag(T,zeros(pp));
Q = blkdiag(Q,H);
R = blkdiag(R,eye(pp));
Pstar = blkdiag(Pstar,H);
Pinf = blkdiag(Pinf,zeros(pp));
H1 = zeros(pp,1);
mmm = mm+pp;
if options_.heteroskedastic_filter
clear QQ
for kv=1:size(Qvec,3)
QQ(:,:,kv) = blkdiag(Qvec(:,:,kv),H);
end
Qvec=QQ;
else
QQ = Q;
end
end
end
a = zeros(mmm,1);
a = set_Kalman_starting_values(a,M_,dr,options_,bayestopt_);
a_0_given_tm1 = T*a;
[dLIK,dlik,a_0_given_tm1,Pstar] = univariate_kalman_filter_d(dataset_info.missing.aindex,...
dataset_info.missing.number_of_observations,...
dataset_info.missing.no_more_missing_observations, ...
2017-05-16 15:10:20 +02:00
Y, 1, size(Y,2), ...
a_0_given_tm1, Pinf, Pstar, ...
kalman_tol, diffuse_kalman_tol, riccati_tol, options_.presample, ...
T,R,QQ,H1,Z,mmm,pp,rr);
diffuse_periods = size(dlik,1);
end
2017-05-16 12:42:01 +02:00
if isnan(dLIK)
2016-06-01 18:22:51 +02:00
fval = Inf;
info(1) = 45;
info(4) = 0.1;
exit_flag = 0;
return
end
2017-05-16 15:10:20 +02:00
case 4% Start from the solution of the Riccati equation.
if kalman_algo ~= 2
kalman_algo = 1;
end
try
if isequal(H,0)
Pstar = kalman_steady_state(transpose(T),R*Q*transpose(R),transpose(build_selection_matrix(Z,mm,length(Z))));
else
Pstar = kalman_steady_state(transpose(T),R*Q*transpose(R),transpose(build_selection_matrix(Z,mm,length(Z))),H);
end
catch ME
disp(ME.message)
disp(['dsge_likelihood:: I am not able to solve the Riccati equation, so I switch to lik_init=1!']);
options_.lik_init = 1;
Pstar=lyapunov_solver(T,R,Q,options_);
end
Pinf = [];
a = zeros(mm,1);
a = set_Kalman_starting_values(a,M_,dr,options_,bayestopt_);
a_0_given_tm1 = T*a;
if options_.occbin.likelihood.status
Z =zeros(length(bayestopt_.mf),size(T,1));
for i = 1:length(bayestopt_.mf)
Z(i,bayestopt_.mf(i))=1;
end
Zflag = 1;
else
Zflag = 0;
end
case 5 % Old diffuse Kalman filter only for the non stationary variables
[eigenvect, eigenv] = eig(T);
eigenv = diag(eigenv);
nstable = length(find(abs(abs(eigenv)-1) > 1e-7));
unstable = find(abs(abs(eigenv)-1) < 1e-7);
V = eigenvect(:,unstable);
indx_unstable = find(sum(abs(V),2)>1e-5);
stable = find(sum(abs(V),2)<1e-5);
nunit = length(eigenv) - nstable;
Pstar = options_.Harvey_scale_factor*eye(nunit);
if kalman_algo ~= 2
kalman_algo = 1;
end
R_tmp = R(stable, :);
T_tmp = T(stable,stable);
Pstar_tmp=lyapunov_solver(T_tmp,R_tmp,Q,options_);
Pstar(stable, stable) = Pstar_tmp;
Pinf = [];
a = zeros(mm,1);
a = set_Kalman_starting_values(a,M_,dr,options_,bayestopt_);
a_0_given_tm1 = T*a;
if options_.occbin.likelihood.status
Z =zeros(length(bayestopt_.mf),size(T,1));
for i = 1:length(bayestopt_.mf)
Z(i,bayestopt_.mf(i))=1;
end
Zflag = 1;
else
Zflag = 0;
end
otherwise
error('dsge_likelihood:: Unknown initialization approach for the Kalman filter!')
end
2017-05-16 12:42:01 +02:00
if analytic_derivation
offset = estim_params_.nvx;
offset = offset+estim_params_.nvn;
offset = offset+estim_params_.ncx;
offset = offset+estim_params_.ncn;
no_DLIK = 0;
full_Hess = analytic_derivation==2;
asy_Hess = analytic_derivation==-2;
outer_product_gradient = analytic_derivation==-1;
2017-05-16 12:42:01 +02:00
if asy_Hess
analytic_derivation=1;
end
2017-05-16 12:42:01 +02:00
if outer_product_gradient
analytic_derivation=1;
end
DLIK = [];
AHess = [];
iv = dr.restrict_var_list;
if nargin<13 || isempty(derivatives_info)
[A,B,nou,nou,dr, M_.params] = dynare_resolve(M_,options_,dr, endo_steady_state, exo_steady_state, exo_det_steady_state);
if ~isempty(estim_params_.var_exo)
indexo=estim_params_.var_exo(:,1);
else
indexo=[];
end
if ~isempty(estim_params_.param_vals)
indparam=estim_params_.param_vals(:,1);
else
indparam=[];
end
old_order = options_.order;
if options_.order > 1%not sure whether this check is necessary
options_.order = 1; fprintf('Reset order to 1 for analytical parameter derivatives.\n');
Added parameter derivatives of perturbation solution up to 3 order # Preliminary comments I finished the identification toolbox at orders two and three using the pruned state space system, but before I merge request this, I decided to first merge the new functionality to compute parameter derivatives of perturbation solution matrices at higher orders. So after this is approved, I merge the identification toolbox. I guess @rattoma, @sebastien, and @michel are best choices to review this. I outline the main idea first and then provide some more detailed changes I made to the functions. *** # Main idea This merge request is concerned with the *analytical*computation of the parameter derivatives of first, second and third order perturbation solution matrices, i.e. using _closed-form_ expressions to efficiently compute the derivative of $g_x$ , $g_u$, $g_{xx}$, $g_{xu}$, $g_{uu}$, $g_{\sigma\sigma}$, $g_{xxx}$, $g_{xxu}$, $g_{xuu}$, $g_{uuu}$, $g_{x\sigma\sigma}$, $g_{u\sigma\sigma}$ *with respect to model parameters* $\theta$. Note that $\theta$ contains model parameters, stderr and corr parameters of shocks. stderr and corr parameters of measurement errors are not yet supported, (they can easily be included as exogenous shocks). The availability of such derivatives is beneficial in terms of more reliable analysis of model sensitivity and parameter identifiability as well as more efficient estimation methods, in particular for models solved up to third order, as it is well-known that numerical derivatives are a tricky business, especially for large models. References for my approach are: * Iskrev (2008, 2010) and Schmitt-Grohé and Uribe (2012, Appendix) who were the first to compute the parameter derivatives analytically at first order, however, using inefficient (sparse) Kronecker products. * Mutschler (2015) who provides the expressions for a second-order, but again using inefficient (sparse) Kronecker products. * Ratto and Iskrev (2012) who show how the first-order system can be solved accurately, fast and efficiently using existing numerical algorithms for generalized Sylvester equations by taking the parameter derivative with respect to each parameter separately. * Julliard and Kamenik (2004) who provide the perturbation solution equation system in tensor notation at any order k. * Levintal (2017) who introduces permutation matrices to express the perturbation solution equation system in matrix notation up to fifth order. Note that @rattoma already implemented the parameter derivatives of $g_x$ and $g_u$ analytically (and numerically), and I rely heavily on his work in `get_first_order_solution_params_derivs.m` (previously `getH.m`). My additions are mainly to this function and thus it is renamed to `get_perturbation_params_derivs.m`. The basic idea of this merge request is to take the second- and third-order perturbation solution systems in Julliard and Kamenik (2004), unfold these into an equivalent matrix representation using permutation matrices as in Levintal (2017). Then extending Ratto and Iskrev (2012) one takes the derivative with respect to each parameter separately and gets a computational problem that is linear, albeit large, as it involves either solving generalized Sylvester equations or taking inverses of highly sparse matrices. I will now briefly summarize the perturbation solution system at third order and the system that results when taking the derivative with respect to parameters. ## Perturbation Solution The following systems arise at first, second, and third order: $(ghx): f_{x} z_{x} = f_{y_{-}^*} + f_{y_0} g_{x} + f_{y_{+}^{**}} g^{**}_{x} g^{*}_{x}= A g_{x} + f_{y_{-}^*}=0$ $(ghu): f_{z} z_{u} = f_{y_0} g_{u} + f_{y_{+}^{**}} g^{**}_{x} g^{*}_{u} + f_{u}= A g_u + f_u = 0$ $(ghxx) : A g_{xx} + B g_{xx} \left(g^{*}_{x} \otimes g^{*}_{x}\right) + f_{zz} \left( z_{x} \otimes z_{x} \right) = 0$ $(ghxu) : A g_{xu} + B g_{xx} \left(g^{*}_{x} \otimes g^{*}_{u}\right) + f_{zz} \left( z_{x} \otimes z_{u} \right) = 0$ $(ghuu) : A g_{uu} + B g_{xx} \left(g^{*}_{u} \otimes g^{*}_{u}\right) + f_{zz} \left( z_{u} \otimes z_{u} \right) = 0$ $(ghs2) : (A+B) g_{\sigma\sigma} + \left( f_{y^{**}_{+}y^{**}_{+}} \left(g^{**}_{u} \otimes g^{**}_{u}\right) + f_{y^{**}_{+}} g^{**}_{uu}\right)vec(\Sigma) = 0$ $(ghxxx) : A g_{xxx} + B g_{xxx} \left(g^{*}_{x} \otimes g^{*}_{x} \otimes g^{*}_{x}\right) + f_{y_{+}}g^{**}_{xx} \left(g^{*}_x \otimes g^{*}_{xx}\right)P_{x\_xx} + f_{zz} \left( z_{x} \otimes z_{xx} \right)P_{x\_xx} + f_{zzz} \left( z_{x} \otimes z_{x} \otimes z_{x} \right) = 0$ $(ghxxu) : A g_{xxu} + B g_{xxx} \left(g^{*}_{x} \otimes g^{*}_{x} \otimes g^{*}_{u}\right) + f_{zzz} \left( z_{x} \otimes z_{x} \otimes z_{u} \right) + f_{zz} \left( \left( z_{x} \otimes z_{xu} \right)P_{x\_xu} + \left(z_{xx} \otimes z_{u}\right) \right) + f_{y_{+}}g^{**}_{xx} \left( \left(g^{*}_{x} \otimes g^{*}_{xu}\right)P_{x\_xu} + \left(g^{*}_{xx} \otimes g^{*}_{u}\right) \right) = 0$ $(ghxuu) : A g_{xuu} + B g_{xxx} \left(g^{*}_{x} \otimes g^{*}_{u} \otimes g^{*}_{u}\right) + f_{zzz} \left( z_{x} \otimes z_{u} \otimes z_{u} \right)+ f_{zz} \left( \left( z_{xu} \otimes z_{u} \right)P_{xu\_u} + \left(z_{x} \otimes z_{uu}\right) \right) + f_{y_{+}}g^{**}_{xx} \left( \left(g^{*}_{xu} \otimes g^{*}_{u}\right)P_{xu\_u} + \left(g^{*}_{x} \otimes g^{*}_{uu}\right) \right) = 0$ $(ghuuu) : A g_{uuu} + B g_{xxx} \left(g^{*}_{u} \otimes g^{*}_{u} \otimes g^{*}_{u}\right) + f_{zzz} \left( z_{u} \otimes z_{u} \otimes z_{u} \right)+ f_{zz} \left( z_{u} \otimes z_{uu} \right)P_{u\_uu} + f_{y_{+}}g^{**}_{xx} \left(g^{*}_{u} \otimes g^{*}_{uu}\right)P_{u\_uu} = 0$ $(ghx\sigma\sigma) : A g_{x\sigma\sigma} + B g_{x\sigma\sigma} g^{*}_x + f_{y_{+}} g^{**}_{xx}\left(g^{*}_{x} \otimes g^{*}_{\sigma\sigma}\right) + f_{zz} \left(z_{x} \otimes z_{\sigma\sigma}\right) + F_{xu_{+}u_{+}}\left(I_{n_x} \otimes vec(\Sigma)\right) = 0$ $F_{xu_{+}u_{+}} = f_{y_{+}^{\ast\ast}} g_{xuu}^{\ast\ast} (g_x^{\ast} \otimes I_{n_u^2}) + f_{zz} \left( \left( z_{xu_{+}} \otimes z_{u_{+}} \right)P_{xu\_u} + \left(z_{x} \otimes z_{u_{+}u_{+}}\right) \right) + f_{zzz}\left(z_{x} \otimes z_{u_{+}} \otimes z_{u_{+}}\right)$ $(ghu\sigma\sigma) : A g_{u\sigma\sigma} + B g_{x\sigma\sigma} g^{*}_{u} + f_{y_{+}} g^{**}_{xx}\left(g^{*}_{u} \otimes g^{*}_{\sigma\sigma}\right) + f_{zz} \left(z_{u} \otimes z_{\sigma\sigma}\right) + F_{uu_{+}u_{+}}\left(I_{n_u} \otimes vec(\Sigma_u)\right) = 0$ $F_{uu_{+}u_{+}} = f_{y_{+}^{\ast\ast}} g_{xuu}^{\ast\ast} (g_u^{\ast} \otimes I_{n_u^2}) + f_{zz} \left( \left( z_{uu_{+}} \otimes z_{u_{+}} \right)P_{uu\_u} + \left(z_{u} \otimes z_{u_{+}u_{+}}\right) \right) + f_{zzz}\left(z_{u} \otimes z_{u_{+}} \otimes z_{u_{+}}\right)$ A and B are the common perturbation matrices: $A = f_{y_0} + \begin{pmatrix} \underbrace{0}_{n\times n_{static}} &\vdots& \underbrace{f_{y^{**}_{+}} \cdot g^{**}_{x}}_{n \times n_{spred}} &\vdots& \underbrace{0}_{n\times n_{frwd}} \end{pmatrix}$and $B = \begin{pmatrix} \underbrace{0}_{n \times n_{static}}&\vdots & \underbrace{0}_{n \times n_{pred}} & \vdots & \underbrace{f_{y^{**}_{+}}}_{n \times n_{sfwrd}} \end{pmatrix}$ and $z=(y_{-}^{\ast}; y; y_{+}^{\ast\ast}; u)$ denotes the dynamic model variables as in `M_.lead_lag_incidence`, $y^\ast$ denote state variables, $y^{\ast\ast}$ denote forward looking variables, $y_+$ denote the variables with a lead, $y_{-}$ denote variables with a lag, $y_0$ denote variables at period t, $f$ the model equations, and $f_z$ the first-order dynamic model derivatives, $f_{zz}$ the second-order dynamic derivatives, and $f_{zzz}$ the third-order dynamic model derivatives. Then: $z_{x} = \begin{pmatrix}I\\g_{x}\\g^{**}_{x} g^{*}_{x}\\0\end{pmatrix}$, $z_{u} =\begin{pmatrix}0\\g_{u}\\g^{**}_{x} \cdot g^{*}_{u}\\I\end{pmatrix}$, $z_{u_{+}} =\begin{pmatrix}0\\0\\g^{**}_{u}\\0\end{pmatrix}$ $z_{xx} = \begin{pmatrix} 0\\g_{xx}\\g^{**}_{x} \left( g^{*}_x \otimes g^{*}_{x} \right) + g^{**}_{x} g^{*}_{x}\\0\end{pmatrix}$, $z_{xu} =\begin{pmatrix}0\\g_{xu}\\g^{**}_{xx} \left( g^{*}_x \otimes g^{*}_{u} \right) + g^{**}_{x} g^{*}_{xu}\\0\end{pmatrix}$, $z_{uu} =\begin{pmatrix}0\\g_{uu}\\g^{**}_{xx} \left( g^{*}_u \otimes g^{*}_{u} \right) + g^{**}_{x} g^{*}_{uu}\\0\end{pmatrix}$, $z_{xu_{+}} =\begin{pmatrix}0\\0\\g^{**}_{xu} \left( g^{*}_x \otimes I \right)\\0\end{pmatrix}$, $z_{uu_{+}} =\begin{pmatrix}0\\0\\g^{**}_{xu} \left( g^{*}_{u} \otimes I \right)\\0\end{pmatrix}$, $z_{u_{+}u_{+}} =\begin{pmatrix}0\\0\\g^{\ast\ast}_{uu}\\0\end{pmatrix}$, $z_{\sigma\sigma} = \begin{pmatrix}0\\ g_{\sigma\sigma}\\ g^{\ast\ast}_{x}g^{\ast}_{\sigma\sigma} + g^{\ast\ast}_{\sigma\sigma}\\0 \end{pmatrix}$ $P$ are permutation matrices that can be computed using Matlab's `ipermute` function. ## Parameter derivatives of perturbation solutions First, we need the parameter derivatives of first, second, third, and fourth derivatives of the dynamic model (i.e. g1,g2,g3,g4 in dynamic files), I make use of the implicit function theorem: Let $f_{z^k}$ denote the kth derivative (wrt all dynamic variables) of the dynamic model, then let $df_{z^k}$ denote the first-derivative (wrt all model parameters) of $f_{z^k}$ evaluated at the steady state. Note that $f_{z^k}$ is a function of both the model parameters $\theta$ and of the steady state of all dynamic variables $\bar{z}$, which also depend on the parameters. Hence, implicitly $f_{z^k}=f_{z^k}(\theta,\bar{z}(\theta))$ and $df_{z^k}$ consists of two parts: 1. direct derivative wrt to all model parameters given by the preprocessor in the `_params_derivs.m` files 2. contribution of derivative of steady state of dynamic variables (wrt all model parameters): $f_{z^{k+1}} \cdot d\bar{z}$ Note that we already have functionality to compute $d\bar{z}$ analytically. Having this, the above perturbation systems are basically equations of the following types $AX +BXC = RHS$ or $AX = RHS$ Now when taking the derivative (wrt to one single parameter $\theta_j$), we get $A\mathrm{d}\{X\} + B\mathrm{d}\{X\}C = \mathrm{d}\{RHS\} - \mathrm{d}\{A\}X - \mathrm{d}\{B\}XC - BX\mathrm{d}\{C\}$ or $A\mathrm{d}\{X\} = \mathrm{d}\{RHS\} - \mathrm{d}\{A\}X$ The first one is a Sylvester type equation, the second one can be solved by taking the inverse of $A$. The only diffculty and tedious work arrises in computing (the highly sparse) derivatives of $RHS$. *** # New functions: ` ## get_perturbation_params_derivs.m`and `get_perturbation_params_derivs_numerical_objective.m` * The parameter derivatives up to third order are computed in the new function`get_perturbation_params_derivs.m` both analytically and numerically. For numerical derivatives `get_perturbation_params_derivs_numerical_objective.m` is the objective for `fjaco.m` or `hessian_sparse.m` or `hessian.m`. * `get_perturbation_params_derivs.m` is basically an extended version of the previous `get_first_order_solution_params_derivs.m` function. * * `get_perturbation_params_derivs_numerical_objective.m`builds upon `identification_numerical_objective.m`. It is used for numerical derivatives, whenever `analytic_derivation_mode=-1|-2`. It takes from `identification_numerical_objective.m` the parts that compute numerical parameter Jacobians of steady state, dynamic model equations, and perturbation solution matrices. Hence, these parts are removed in `identification_numerical_objective.m` and it only computes numerical parameter Jacobian of moments and spectrum which are needed for identification analysis in `get_identification_jacobians.m`, when `analytic_derivation_mode=-1` only. * Detailed changes: * Most important: notation of this function is now in accordance to the k_order_solver, i.e. we do not compute derivatives of Kalman transition matrices A and B, but rather the solution matrices ghx,ghu,ghxx,ghxu,ghuu,ghs2,ghxxx,ghxxu,ghxuu,ghuuu,ghxss,ghuss in accordance with notation used in `oo_.dr`. As a byproduct at first-order, focusing on ghx and ghu instead of Kalman transition matrices A and B makes the computations slightly faster for large models (e.g. for Quest the computations were faster by a couple of seconds, not much, but okay). * Removed use of `kstate`, see also Dynare/dynare#1653 and Dynare/dynare!1656 * Output arguments are stored in a structure `DERIVS`, there is also a flag `d2flag` that computes parameter hessians needed only in `dsge_likelihood.m`. * Removed `kronflag` as input. `options_.analytic_derivation_mode` is now used instead of `kronflag`. * Removed `indvar`, an index that was used to selected specific variables in the derivatives. This is not needed, as we always compute the parameter derivatives for all variables first and then select a subset of variables. The selection now takes place in other functions, like `dsge_likelihood.m`. * Introduced some checks: (i) deterministic exogenous variables are not supported, (ii) Kronecker method only compatible with first-order approximation so reset to sylvester method, (iii) for purely backward or forward models we need to be careful with the rows in `M_.lead_la g_incidence`, (iv) if `_params_derivs.m` files are missing an error is thrown. * For numerical derivatives, if mod file does not contain an `estimated_params_block`, a temporary one with the most important parameter information is created. ## `unfold_g4.m` * When evaluating g3 and g4 one needs to take into account that these do not contain symmetric elements, so one needs to use `unfold_g3.m` and the new function `unfold_g4.m`. This returns an unfolded version of the same matrix (i.e. with symmetric elements). *** # New test models `.gitignore` and `Makefile.am` are changed accordingly. Also now it is possible to run test suite on analytic_derivatives, i.e. run `make check m/analytic_derivatives` ## `analytic_derivatives/BrockMirman_PertParamsDerivs.mod` * This is the Brock Mirman model, where we know the exact policy function $g$ for capital and consumption. As this does not imply a nonzero $g_{\sigma\sigma}$, $g_{x\sigma\sigma}$, $g_{u\sigma\sigma}$ I added some artificial equations to get nonzero solution matrices with respect to $\sigma$. The true perturbation solution matrices $g_x$ , $g_u$, $g_{xx}$, $g_{xu}$, $g_{uu}$, $g_{\sigma\sigma}$, $g_{xxx}$, $g_{xxu}$, $g_{xuu}$, $g_{uuu}$, $g_{x\sigma\sigma}$, $g_{u\sigma\sigma}$ are then computed analytically with Matlab's symbolic toolbox and saved in `nBrockMirmanSYM.mat`. There is a preprocessor flag that recreates these analytical computations if changes are needed (and to check whether I made some errors here ;-) ) * Then solution matrices up to third order and their parameter Jacobians are then compared to the ones computed by Dynare's `k_order_solver` and by `get_perturbation_params_derivs` for all `analytic_derivation_mode`'s. There will be an error if the maximum absolute deviation is too large, i.e. for numerical derivatives (`analytic_derivation_mode=-1|-2`) the tolerance is choosen lower (around 1e-5); for analytical methods we are stricter: around 1e-13 for first-order, 1e-12 for second order, and 1e-11 for third-order. * As a side note, this mod file also checks Dynare's `k_order_solver` algorithm and throws an error if something is wrong. * This test model shows that the new functionality works well. And analytical derivatives perform way better and accurate than numerical ones, even for this small model. ## `analytic_derivatives/burnside_3_order_PertParamsDerivs.mod` * This builds upon `tests/k_order_perturbation/burnside_k_order.mod` and computes the true parameter derivatives analytically by hand. * This test model also shows that the new functionality works well. ## `analytic_derivatives/LindeTrabandt2019.mod` * Shows that the new functionality also works for medium-sized models, i.e. a SW type model solved at third order with 35 variables (11 states). 2 shocks and 20 parameters. * This mod file can be used to tweak the speed of the computations in the future. * Compares numerical versus analytical parameter derivatives (for first, second and third order). Note that this model clearly shows that numerical ones are quite different than analytical ones even at first order! ## `identification/LindeTrabandt2019_xfail.mod` * This model is a check for issue Dynare/dynare#1595, see fjaco.m below, and will fail. * Removed `analytic_derivatives/ls2003.mod` as this mod file is neither in the testsuite nor does it work. *** # Detailed changes in other functions ## `get_first_order_solution_params_derivs.m` * Deleted, or actually, renamed to `get_perturbation_params_derivs.m`, as this function now is able to compute the derivatives up to third order ## `identification_numerical_objective.m` * `get_perturbation_params_derivs_numerical_objective.m`builds upon `identification_numerical_objective.m`. It takes from `identification_numerical_objective.m` the parts that compute numerical parameter Jacobians of steady state, dynamic model equations, and perturbation solution matrices. Hence, these parts are removed in `identification_numerical_objective.m` and it only computes numerical parameter Jacobian of moments and spectrum which are needed for identification analysis in `get_identification_jacobians.m`, when `analytic_derivation_mode=-1` only. ## `dsge_likelihood.m` * As `get_first_order_solution_params_derivs.m`is renamed to `get_perturbation_params_derivs.m`, the call is adapted. That is,`get_perturbation_params_derivs` does not compute the derivatives of the Kalman transition `T`matrix anymore, but instead of the dynare solution matrix `ghx`. So we recreate `T` here as this amounts to adding some zeros and focusing on selected variables only. * Added some checks to make sure the first-order approximation is selected. * Removed `kron_flag` as input, as `get_perturbation_params_derivs` looks into `options_.analytic_derivation_mode` for `kron_flag`. ## `dynare_identification.m` * make sure that setting `analytic_derivation_mode` is set both in `options_ident` and `options_`. Note that at the end of the function we restore the `options_` structure, so all changes are local. In a next merge request, I will remove the global variables to make all variables local. ## `get_identification_jacobians.m` * As `get_first_order_solution_params_derivs.m`is renamed to `get_perturbation_params_derivs.m`, the call is adapted. That is,`get_perturbation_params_derivs` does not compute the derivatives of the Kalman transition `A` and `B` matrix anymore, but instead of the dynare solution matrix `ghx` and `ghu`. So we recreate these matrices here instead of in `get_perturbation_params_derivs.m`. * Added `str2func` for better function handles in `fjaco.m`. ## `fjaco.m` * make `tol`an option, which can be adjusted by changing `options_.dynatol.x`for identification and parameter derivatives purposes. * include a check and an informative error message, if numerical derivatives (two-sided finite difference method) yield errors in `resol.m` for identification and parameter derivatives purposes. This closes issue Dynare/dynare#1595. * Changed year of copyright to 2010-2017,2019 *** # Further suggestions and questions * Ones this is merged, I will merge request an improvement of the identification toolbox, which will work up to third order using the pruned state space. This will also remove some issues and bugs, and also I will remove global variables in this request. * The third-order derivatives can be further improved by taking sparsity into account and use mex versions for kronecker products etc. I leave this for further testing (and if anybody actually uses this ;-) )
2019-12-17 19:17:09 +01:00
end
old_analytic_derivation_mode = options_.analytic_derivation_mode;
options_.analytic_derivation_mode = kron_flag;
2017-05-16 12:42:01 +02:00
if full_Hess
DERIVS = get_perturbation_params_derivs(M_, options_, estim_params_, dr, endo_steady_state, exo_steady_state, exo_det_steady_state, indparam, indexo, [], true);
indD2T = reshape(1:M_.endo_nbr^2, M_.endo_nbr, M_.endo_nbr);
indD2Om = dyn_unvech(1:M_.endo_nbr*(M_.endo_nbr+1)/2);
Added parameter derivatives of perturbation solution up to 3 order # Preliminary comments I finished the identification toolbox at orders two and three using the pruned state space system, but before I merge request this, I decided to first merge the new functionality to compute parameter derivatives of perturbation solution matrices at higher orders. So after this is approved, I merge the identification toolbox. I guess @rattoma, @sebastien, and @michel are best choices to review this. I outline the main idea first and then provide some more detailed changes I made to the functions. *** # Main idea This merge request is concerned with the *analytical*computation of the parameter derivatives of first, second and third order perturbation solution matrices, i.e. using _closed-form_ expressions to efficiently compute the derivative of $g_x$ , $g_u$, $g_{xx}$, $g_{xu}$, $g_{uu}$, $g_{\sigma\sigma}$, $g_{xxx}$, $g_{xxu}$, $g_{xuu}$, $g_{uuu}$, $g_{x\sigma\sigma}$, $g_{u\sigma\sigma}$ *with respect to model parameters* $\theta$. Note that $\theta$ contains model parameters, stderr and corr parameters of shocks. stderr and corr parameters of measurement errors are not yet supported, (they can easily be included as exogenous shocks). The availability of such derivatives is beneficial in terms of more reliable analysis of model sensitivity and parameter identifiability as well as more efficient estimation methods, in particular for models solved up to third order, as it is well-known that numerical derivatives are a tricky business, especially for large models. References for my approach are: * Iskrev (2008, 2010) and Schmitt-Grohé and Uribe (2012, Appendix) who were the first to compute the parameter derivatives analytically at first order, however, using inefficient (sparse) Kronecker products. * Mutschler (2015) who provides the expressions for a second-order, but again using inefficient (sparse) Kronecker products. * Ratto and Iskrev (2012) who show how the first-order system can be solved accurately, fast and efficiently using existing numerical algorithms for generalized Sylvester equations by taking the parameter derivative with respect to each parameter separately. * Julliard and Kamenik (2004) who provide the perturbation solution equation system in tensor notation at any order k. * Levintal (2017) who introduces permutation matrices to express the perturbation solution equation system in matrix notation up to fifth order. Note that @rattoma already implemented the parameter derivatives of $g_x$ and $g_u$ analytically (and numerically), and I rely heavily on his work in `get_first_order_solution_params_derivs.m` (previously `getH.m`). My additions are mainly to this function and thus it is renamed to `get_perturbation_params_derivs.m`. The basic idea of this merge request is to take the second- and third-order perturbation solution systems in Julliard and Kamenik (2004), unfold these into an equivalent matrix representation using permutation matrices as in Levintal (2017). Then extending Ratto and Iskrev (2012) one takes the derivative with respect to each parameter separately and gets a computational problem that is linear, albeit large, as it involves either solving generalized Sylvester equations or taking inverses of highly sparse matrices. I will now briefly summarize the perturbation solution system at third order and the system that results when taking the derivative with respect to parameters. ## Perturbation Solution The following systems arise at first, second, and third order: $(ghx): f_{x} z_{x} = f_{y_{-}^*} + f_{y_0} g_{x} + f_{y_{+}^{**}} g^{**}_{x} g^{*}_{x}= A g_{x} + f_{y_{-}^*}=0$ $(ghu): f_{z} z_{u} = f_{y_0} g_{u} + f_{y_{+}^{**}} g^{**}_{x} g^{*}_{u} + f_{u}= A g_u + f_u = 0$ $(ghxx) : A g_{xx} + B g_{xx} \left(g^{*}_{x} \otimes g^{*}_{x}\right) + f_{zz} \left( z_{x} \otimes z_{x} \right) = 0$ $(ghxu) : A g_{xu} + B g_{xx} \left(g^{*}_{x} \otimes g^{*}_{u}\right) + f_{zz} \left( z_{x} \otimes z_{u} \right) = 0$ $(ghuu) : A g_{uu} + B g_{xx} \left(g^{*}_{u} \otimes g^{*}_{u}\right) + f_{zz} \left( z_{u} \otimes z_{u} \right) = 0$ $(ghs2) : (A+B) g_{\sigma\sigma} + \left( f_{y^{**}_{+}y^{**}_{+}} \left(g^{**}_{u} \otimes g^{**}_{u}\right) + f_{y^{**}_{+}} g^{**}_{uu}\right)vec(\Sigma) = 0$ $(ghxxx) : A g_{xxx} + B g_{xxx} \left(g^{*}_{x} \otimes g^{*}_{x} \otimes g^{*}_{x}\right) + f_{y_{+}}g^{**}_{xx} \left(g^{*}_x \otimes g^{*}_{xx}\right)P_{x\_xx} + f_{zz} \left( z_{x} \otimes z_{xx} \right)P_{x\_xx} + f_{zzz} \left( z_{x} \otimes z_{x} \otimes z_{x} \right) = 0$ $(ghxxu) : A g_{xxu} + B g_{xxx} \left(g^{*}_{x} \otimes g^{*}_{x} \otimes g^{*}_{u}\right) + f_{zzz} \left( z_{x} \otimes z_{x} \otimes z_{u} \right) + f_{zz} \left( \left( z_{x} \otimes z_{xu} \right)P_{x\_xu} + \left(z_{xx} \otimes z_{u}\right) \right) + f_{y_{+}}g^{**}_{xx} \left( \left(g^{*}_{x} \otimes g^{*}_{xu}\right)P_{x\_xu} + \left(g^{*}_{xx} \otimes g^{*}_{u}\right) \right) = 0$ $(ghxuu) : A g_{xuu} + B g_{xxx} \left(g^{*}_{x} \otimes g^{*}_{u} \otimes g^{*}_{u}\right) + f_{zzz} \left( z_{x} \otimes z_{u} \otimes z_{u} \right)+ f_{zz} \left( \left( z_{xu} \otimes z_{u} \right)P_{xu\_u} + \left(z_{x} \otimes z_{uu}\right) \right) + f_{y_{+}}g^{**}_{xx} \left( \left(g^{*}_{xu} \otimes g^{*}_{u}\right)P_{xu\_u} + \left(g^{*}_{x} \otimes g^{*}_{uu}\right) \right) = 0$ $(ghuuu) : A g_{uuu} + B g_{xxx} \left(g^{*}_{u} \otimes g^{*}_{u} \otimes g^{*}_{u}\right) + f_{zzz} \left( z_{u} \otimes z_{u} \otimes z_{u} \right)+ f_{zz} \left( z_{u} \otimes z_{uu} \right)P_{u\_uu} + f_{y_{+}}g^{**}_{xx} \left(g^{*}_{u} \otimes g^{*}_{uu}\right)P_{u\_uu} = 0$ $(ghx\sigma\sigma) : A g_{x\sigma\sigma} + B g_{x\sigma\sigma} g^{*}_x + f_{y_{+}} g^{**}_{xx}\left(g^{*}_{x} \otimes g^{*}_{\sigma\sigma}\right) + f_{zz} \left(z_{x} \otimes z_{\sigma\sigma}\right) + F_{xu_{+}u_{+}}\left(I_{n_x} \otimes vec(\Sigma)\right) = 0$ $F_{xu_{+}u_{+}} = f_{y_{+}^{\ast\ast}} g_{xuu}^{\ast\ast} (g_x^{\ast} \otimes I_{n_u^2}) + f_{zz} \left( \left( z_{xu_{+}} \otimes z_{u_{+}} \right)P_{xu\_u} + \left(z_{x} \otimes z_{u_{+}u_{+}}\right) \right) + f_{zzz}\left(z_{x} \otimes z_{u_{+}} \otimes z_{u_{+}}\right)$ $(ghu\sigma\sigma) : A g_{u\sigma\sigma} + B g_{x\sigma\sigma} g^{*}_{u} + f_{y_{+}} g^{**}_{xx}\left(g^{*}_{u} \otimes g^{*}_{\sigma\sigma}\right) + f_{zz} \left(z_{u} \otimes z_{\sigma\sigma}\right) + F_{uu_{+}u_{+}}\left(I_{n_u} \otimes vec(\Sigma_u)\right) = 0$ $F_{uu_{+}u_{+}} = f_{y_{+}^{\ast\ast}} g_{xuu}^{\ast\ast} (g_u^{\ast} \otimes I_{n_u^2}) + f_{zz} \left( \left( z_{uu_{+}} \otimes z_{u_{+}} \right)P_{uu\_u} + \left(z_{u} \otimes z_{u_{+}u_{+}}\right) \right) + f_{zzz}\left(z_{u} \otimes z_{u_{+}} \otimes z_{u_{+}}\right)$ A and B are the common perturbation matrices: $A = f_{y_0} + \begin{pmatrix} \underbrace{0}_{n\times n_{static}} &\vdots& \underbrace{f_{y^{**}_{+}} \cdot g^{**}_{x}}_{n \times n_{spred}} &\vdots& \underbrace{0}_{n\times n_{frwd}} \end{pmatrix}$and $B = \begin{pmatrix} \underbrace{0}_{n \times n_{static}}&\vdots & \underbrace{0}_{n \times n_{pred}} & \vdots & \underbrace{f_{y^{**}_{+}}}_{n \times n_{sfwrd}} \end{pmatrix}$ and $z=(y_{-}^{\ast}; y; y_{+}^{\ast\ast}; u)$ denotes the dynamic model variables as in `M_.lead_lag_incidence`, $y^\ast$ denote state variables, $y^{\ast\ast}$ denote forward looking variables, $y_+$ denote the variables with a lead, $y_{-}$ denote variables with a lag, $y_0$ denote variables at period t, $f$ the model equations, and $f_z$ the first-order dynamic model derivatives, $f_{zz}$ the second-order dynamic derivatives, and $f_{zzz}$ the third-order dynamic model derivatives. Then: $z_{x} = \begin{pmatrix}I\\g_{x}\\g^{**}_{x} g^{*}_{x}\\0\end{pmatrix}$, $z_{u} =\begin{pmatrix}0\\g_{u}\\g^{**}_{x} \cdot g^{*}_{u}\\I\end{pmatrix}$, $z_{u_{+}} =\begin{pmatrix}0\\0\\g^{**}_{u}\\0\end{pmatrix}$ $z_{xx} = \begin{pmatrix} 0\\g_{xx}\\g^{**}_{x} \left( g^{*}_x \otimes g^{*}_{x} \right) + g^{**}_{x} g^{*}_{x}\\0\end{pmatrix}$, $z_{xu} =\begin{pmatrix}0\\g_{xu}\\g^{**}_{xx} \left( g^{*}_x \otimes g^{*}_{u} \right) + g^{**}_{x} g^{*}_{xu}\\0\end{pmatrix}$, $z_{uu} =\begin{pmatrix}0\\g_{uu}\\g^{**}_{xx} \left( g^{*}_u \otimes g^{*}_{u} \right) + g^{**}_{x} g^{*}_{uu}\\0\end{pmatrix}$, $z_{xu_{+}} =\begin{pmatrix}0\\0\\g^{**}_{xu} \left( g^{*}_x \otimes I \right)\\0\end{pmatrix}$, $z_{uu_{+}} =\begin{pmatrix}0\\0\\g^{**}_{xu} \left( g^{*}_{u} \otimes I \right)\\0\end{pmatrix}$, $z_{u_{+}u_{+}} =\begin{pmatrix}0\\0\\g^{\ast\ast}_{uu}\\0\end{pmatrix}$, $z_{\sigma\sigma} = \begin{pmatrix}0\\ g_{\sigma\sigma}\\ g^{\ast\ast}_{x}g^{\ast}_{\sigma\sigma} + g^{\ast\ast}_{\sigma\sigma}\\0 \end{pmatrix}$ $P$ are permutation matrices that can be computed using Matlab's `ipermute` function. ## Parameter derivatives of perturbation solutions First, we need the parameter derivatives of first, second, third, and fourth derivatives of the dynamic model (i.e. g1,g2,g3,g4 in dynamic files), I make use of the implicit function theorem: Let $f_{z^k}$ denote the kth derivative (wrt all dynamic variables) of the dynamic model, then let $df_{z^k}$ denote the first-derivative (wrt all model parameters) of $f_{z^k}$ evaluated at the steady state. Note that $f_{z^k}$ is a function of both the model parameters $\theta$ and of the steady state of all dynamic variables $\bar{z}$, which also depend on the parameters. Hence, implicitly $f_{z^k}=f_{z^k}(\theta,\bar{z}(\theta))$ and $df_{z^k}$ consists of two parts: 1. direct derivative wrt to all model parameters given by the preprocessor in the `_params_derivs.m` files 2. contribution of derivative of steady state of dynamic variables (wrt all model parameters): $f_{z^{k+1}} \cdot d\bar{z}$ Note that we already have functionality to compute $d\bar{z}$ analytically. Having this, the above perturbation systems are basically equations of the following types $AX +BXC = RHS$ or $AX = RHS$ Now when taking the derivative (wrt to one single parameter $\theta_j$), we get $A\mathrm{d}\{X\} + B\mathrm{d}\{X\}C = \mathrm{d}\{RHS\} - \mathrm{d}\{A\}X - \mathrm{d}\{B\}XC - BX\mathrm{d}\{C\}$ or $A\mathrm{d}\{X\} = \mathrm{d}\{RHS\} - \mathrm{d}\{A\}X$ The first one is a Sylvester type equation, the second one can be solved by taking the inverse of $A$. The only diffculty and tedious work arrises in computing (the highly sparse) derivatives of $RHS$. *** # New functions: ` ## get_perturbation_params_derivs.m`and `get_perturbation_params_derivs_numerical_objective.m` * The parameter derivatives up to third order are computed in the new function`get_perturbation_params_derivs.m` both analytically and numerically. For numerical derivatives `get_perturbation_params_derivs_numerical_objective.m` is the objective for `fjaco.m` or `hessian_sparse.m` or `hessian.m`. * `get_perturbation_params_derivs.m` is basically an extended version of the previous `get_first_order_solution_params_derivs.m` function. * * `get_perturbation_params_derivs_numerical_objective.m`builds upon `identification_numerical_objective.m`. It is used for numerical derivatives, whenever `analytic_derivation_mode=-1|-2`. It takes from `identification_numerical_objective.m` the parts that compute numerical parameter Jacobians of steady state, dynamic model equations, and perturbation solution matrices. Hence, these parts are removed in `identification_numerical_objective.m` and it only computes numerical parameter Jacobian of moments and spectrum which are needed for identification analysis in `get_identification_jacobians.m`, when `analytic_derivation_mode=-1` only. * Detailed changes: * Most important: notation of this function is now in accordance to the k_order_solver, i.e. we do not compute derivatives of Kalman transition matrices A and B, but rather the solution matrices ghx,ghu,ghxx,ghxu,ghuu,ghs2,ghxxx,ghxxu,ghxuu,ghuuu,ghxss,ghuss in accordance with notation used in `oo_.dr`. As a byproduct at first-order, focusing on ghx and ghu instead of Kalman transition matrices A and B makes the computations slightly faster for large models (e.g. for Quest the computations were faster by a couple of seconds, not much, but okay). * Removed use of `kstate`, see also Dynare/dynare#1653 and Dynare/dynare!1656 * Output arguments are stored in a structure `DERIVS`, there is also a flag `d2flag` that computes parameter hessians needed only in `dsge_likelihood.m`. * Removed `kronflag` as input. `options_.analytic_derivation_mode` is now used instead of `kronflag`. * Removed `indvar`, an index that was used to selected specific variables in the derivatives. This is not needed, as we always compute the parameter derivatives for all variables first and then select a subset of variables. The selection now takes place in other functions, like `dsge_likelihood.m`. * Introduced some checks: (i) deterministic exogenous variables are not supported, (ii) Kronecker method only compatible with first-order approximation so reset to sylvester method, (iii) for purely backward or forward models we need to be careful with the rows in `M_.lead_la g_incidence`, (iv) if `_params_derivs.m` files are missing an error is thrown. * For numerical derivatives, if mod file does not contain an `estimated_params_block`, a temporary one with the most important parameter information is created. ## `unfold_g4.m` * When evaluating g3 and g4 one needs to take into account that these do not contain symmetric elements, so one needs to use `unfold_g3.m` and the new function `unfold_g4.m`. This returns an unfolded version of the same matrix (i.e. with symmetric elements). *** # New test models `.gitignore` and `Makefile.am` are changed accordingly. Also now it is possible to run test suite on analytic_derivatives, i.e. run `make check m/analytic_derivatives` ## `analytic_derivatives/BrockMirman_PertParamsDerivs.mod` * This is the Brock Mirman model, where we know the exact policy function $g$ for capital and consumption. As this does not imply a nonzero $g_{\sigma\sigma}$, $g_{x\sigma\sigma}$, $g_{u\sigma\sigma}$ I added some artificial equations to get nonzero solution matrices with respect to $\sigma$. The true perturbation solution matrices $g_x$ , $g_u$, $g_{xx}$, $g_{xu}$, $g_{uu}$, $g_{\sigma\sigma}$, $g_{xxx}$, $g_{xxu}$, $g_{xuu}$, $g_{uuu}$, $g_{x\sigma\sigma}$, $g_{u\sigma\sigma}$ are then computed analytically with Matlab's symbolic toolbox and saved in `nBrockMirmanSYM.mat`. There is a preprocessor flag that recreates these analytical computations if changes are needed (and to check whether I made some errors here ;-) ) * Then solution matrices up to third order and their parameter Jacobians are then compared to the ones computed by Dynare's `k_order_solver` and by `get_perturbation_params_derivs` for all `analytic_derivation_mode`'s. There will be an error if the maximum absolute deviation is too large, i.e. for numerical derivatives (`analytic_derivation_mode=-1|-2`) the tolerance is choosen lower (around 1e-5); for analytical methods we are stricter: around 1e-13 for first-order, 1e-12 for second order, and 1e-11 for third-order. * As a side note, this mod file also checks Dynare's `k_order_solver` algorithm and throws an error if something is wrong. * This test model shows that the new functionality works well. And analytical derivatives perform way better and accurate than numerical ones, even for this small model. ## `analytic_derivatives/burnside_3_order_PertParamsDerivs.mod` * This builds upon `tests/k_order_perturbation/burnside_k_order.mod` and computes the true parameter derivatives analytically by hand. * This test model also shows that the new functionality works well. ## `analytic_derivatives/LindeTrabandt2019.mod` * Shows that the new functionality also works for medium-sized models, i.e. a SW type model solved at third order with 35 variables (11 states). 2 shocks and 20 parameters. * This mod file can be used to tweak the speed of the computations in the future. * Compares numerical versus analytical parameter derivatives (for first, second and third order). Note that this model clearly shows that numerical ones are quite different than analytical ones even at first order! ## `identification/LindeTrabandt2019_xfail.mod` * This model is a check for issue Dynare/dynare#1595, see fjaco.m below, and will fail. * Removed `analytic_derivatives/ls2003.mod` as this mod file is neither in the testsuite nor does it work. *** # Detailed changes in other functions ## `get_first_order_solution_params_derivs.m` * Deleted, or actually, renamed to `get_perturbation_params_derivs.m`, as this function now is able to compute the derivatives up to third order ## `identification_numerical_objective.m` * `get_perturbation_params_derivs_numerical_objective.m`builds upon `identification_numerical_objective.m`. It takes from `identification_numerical_objective.m` the parts that compute numerical parameter Jacobians of steady state, dynamic model equations, and perturbation solution matrices. Hence, these parts are removed in `identification_numerical_objective.m` and it only computes numerical parameter Jacobian of moments and spectrum which are needed for identification analysis in `get_identification_jacobians.m`, when `analytic_derivation_mode=-1` only. ## `dsge_likelihood.m` * As `get_first_order_solution_params_derivs.m`is renamed to `get_perturbation_params_derivs.m`, the call is adapted. That is,`get_perturbation_params_derivs` does not compute the derivatives of the Kalman transition `T`matrix anymore, but instead of the dynare solution matrix `ghx`. So we recreate `T` here as this amounts to adding some zeros and focusing on selected variables only. * Added some checks to make sure the first-order approximation is selected. * Removed `kron_flag` as input, as `get_perturbation_params_derivs` looks into `options_.analytic_derivation_mode` for `kron_flag`. ## `dynare_identification.m` * make sure that setting `analytic_derivation_mode` is set both in `options_ident` and `options_`. Note that at the end of the function we restore the `options_` structure, so all changes are local. In a next merge request, I will remove the global variables to make all variables local. ## `get_identification_jacobians.m` * As `get_first_order_solution_params_derivs.m`is renamed to `get_perturbation_params_derivs.m`, the call is adapted. That is,`get_perturbation_params_derivs` does not compute the derivatives of the Kalman transition `A` and `B` matrix anymore, but instead of the dynare solution matrix `ghx` and `ghu`. So we recreate these matrices here instead of in `get_perturbation_params_derivs.m`. * Added `str2func` for better function handles in `fjaco.m`. ## `fjaco.m` * make `tol`an option, which can be adjusted by changing `options_.dynatol.x`for identification and parameter derivatives purposes. * include a check and an informative error message, if numerical derivatives (two-sided finite difference method) yield errors in `resol.m` for identification and parameter derivatives purposes. This closes issue Dynare/dynare#1595. * Changed year of copyright to 2010-2017,2019 *** # Further suggestions and questions * Ones this is merged, I will merge request an improvement of the identification toolbox, which will work up to third order using the pruned state space. This will also remove some issues and bugs, and also I will remove global variables in this request. * The third-order derivatives can be further improved by taking sparsity into account and use mex versions for kronecker products etc. I leave this for further testing (and if anybody actually uses this ;-) )
2019-12-17 19:17:09 +01:00
D2T = DERIVS.d2KalmanA(indD2T(iv,iv),:);
D2Om = DERIVS.d2Om(dyn_vech(indD2Om(iv,iv)),:);
D2Yss = DERIVS.d2Yss(iv,:,:);
else
DERIVS = get_perturbation_params_derivs(M_, options_, estim_params_, dr, endo_steady_state, exo_steady_state, exo_det_steady_state, indparam, indexo, [], false);
end
DT = zeros(M_.endo_nbr, M_.endo_nbr, size(DERIVS.dghx,3));
DT(:,M_.nstatic+(1:M_.nspred),:) = DERIVS.dghx;
Added parameter derivatives of perturbation solution up to 3 order # Preliminary comments I finished the identification toolbox at orders two and three using the pruned state space system, but before I merge request this, I decided to first merge the new functionality to compute parameter derivatives of perturbation solution matrices at higher orders. So after this is approved, I merge the identification toolbox. I guess @rattoma, @sebastien, and @michel are best choices to review this. I outline the main idea first and then provide some more detailed changes I made to the functions. *** # Main idea This merge request is concerned with the *analytical*computation of the parameter derivatives of first, second and third order perturbation solution matrices, i.e. using _closed-form_ expressions to efficiently compute the derivative of $g_x$ , $g_u$, $g_{xx}$, $g_{xu}$, $g_{uu}$, $g_{\sigma\sigma}$, $g_{xxx}$, $g_{xxu}$, $g_{xuu}$, $g_{uuu}$, $g_{x\sigma\sigma}$, $g_{u\sigma\sigma}$ *with respect to model parameters* $\theta$. Note that $\theta$ contains model parameters, stderr and corr parameters of shocks. stderr and corr parameters of measurement errors are not yet supported, (they can easily be included as exogenous shocks). The availability of such derivatives is beneficial in terms of more reliable analysis of model sensitivity and parameter identifiability as well as more efficient estimation methods, in particular for models solved up to third order, as it is well-known that numerical derivatives are a tricky business, especially for large models. References for my approach are: * Iskrev (2008, 2010) and Schmitt-Grohé and Uribe (2012, Appendix) who were the first to compute the parameter derivatives analytically at first order, however, using inefficient (sparse) Kronecker products. * Mutschler (2015) who provides the expressions for a second-order, but again using inefficient (sparse) Kronecker products. * Ratto and Iskrev (2012) who show how the first-order system can be solved accurately, fast and efficiently using existing numerical algorithms for generalized Sylvester equations by taking the parameter derivative with respect to each parameter separately. * Julliard and Kamenik (2004) who provide the perturbation solution equation system in tensor notation at any order k. * Levintal (2017) who introduces permutation matrices to express the perturbation solution equation system in matrix notation up to fifth order. Note that @rattoma already implemented the parameter derivatives of $g_x$ and $g_u$ analytically (and numerically), and I rely heavily on his work in `get_first_order_solution_params_derivs.m` (previously `getH.m`). My additions are mainly to this function and thus it is renamed to `get_perturbation_params_derivs.m`. The basic idea of this merge request is to take the second- and third-order perturbation solution systems in Julliard and Kamenik (2004), unfold these into an equivalent matrix representation using permutation matrices as in Levintal (2017). Then extending Ratto and Iskrev (2012) one takes the derivative with respect to each parameter separately and gets a computational problem that is linear, albeit large, as it involves either solving generalized Sylvester equations or taking inverses of highly sparse matrices. I will now briefly summarize the perturbation solution system at third order and the system that results when taking the derivative with respect to parameters. ## Perturbation Solution The following systems arise at first, second, and third order: $(ghx): f_{x} z_{x} = f_{y_{-}^*} + f_{y_0} g_{x} + f_{y_{+}^{**}} g^{**}_{x} g^{*}_{x}= A g_{x} + f_{y_{-}^*}=0$ $(ghu): f_{z} z_{u} = f_{y_0} g_{u} + f_{y_{+}^{**}} g^{**}_{x} g^{*}_{u} + f_{u}= A g_u + f_u = 0$ $(ghxx) : A g_{xx} + B g_{xx} \left(g^{*}_{x} \otimes g^{*}_{x}\right) + f_{zz} \left( z_{x} \otimes z_{x} \right) = 0$ $(ghxu) : A g_{xu} + B g_{xx} \left(g^{*}_{x} \otimes g^{*}_{u}\right) + f_{zz} \left( z_{x} \otimes z_{u} \right) = 0$ $(ghuu) : A g_{uu} + B g_{xx} \left(g^{*}_{u} \otimes g^{*}_{u}\right) + f_{zz} \left( z_{u} \otimes z_{u} \right) = 0$ $(ghs2) : (A+B) g_{\sigma\sigma} + \left( f_{y^{**}_{+}y^{**}_{+}} \left(g^{**}_{u} \otimes g^{**}_{u}\right) + f_{y^{**}_{+}} g^{**}_{uu}\right)vec(\Sigma) = 0$ $(ghxxx) : A g_{xxx} + B g_{xxx} \left(g^{*}_{x} \otimes g^{*}_{x} \otimes g^{*}_{x}\right) + f_{y_{+}}g^{**}_{xx} \left(g^{*}_x \otimes g^{*}_{xx}\right)P_{x\_xx} + f_{zz} \left( z_{x} \otimes z_{xx} \right)P_{x\_xx} + f_{zzz} \left( z_{x} \otimes z_{x} \otimes z_{x} \right) = 0$ $(ghxxu) : A g_{xxu} + B g_{xxx} \left(g^{*}_{x} \otimes g^{*}_{x} \otimes g^{*}_{u}\right) + f_{zzz} \left( z_{x} \otimes z_{x} \otimes z_{u} \right) + f_{zz} \left( \left( z_{x} \otimes z_{xu} \right)P_{x\_xu} + \left(z_{xx} \otimes z_{u}\right) \right) + f_{y_{+}}g^{**}_{xx} \left( \left(g^{*}_{x} \otimes g^{*}_{xu}\right)P_{x\_xu} + \left(g^{*}_{xx} \otimes g^{*}_{u}\right) \right) = 0$ $(ghxuu) : A g_{xuu} + B g_{xxx} \left(g^{*}_{x} \otimes g^{*}_{u} \otimes g^{*}_{u}\right) + f_{zzz} \left( z_{x} \otimes z_{u} \otimes z_{u} \right)+ f_{zz} \left( \left( z_{xu} \otimes z_{u} \right)P_{xu\_u} + \left(z_{x} \otimes z_{uu}\right) \right) + f_{y_{+}}g^{**}_{xx} \left( \left(g^{*}_{xu} \otimes g^{*}_{u}\right)P_{xu\_u} + \left(g^{*}_{x} \otimes g^{*}_{uu}\right) \right) = 0$ $(ghuuu) : A g_{uuu} + B g_{xxx} \left(g^{*}_{u} \otimes g^{*}_{u} \otimes g^{*}_{u}\right) + f_{zzz} \left( z_{u} \otimes z_{u} \otimes z_{u} \right)+ f_{zz} \left( z_{u} \otimes z_{uu} \right)P_{u\_uu} + f_{y_{+}}g^{**}_{xx} \left(g^{*}_{u} \otimes g^{*}_{uu}\right)P_{u\_uu} = 0$ $(ghx\sigma\sigma) : A g_{x\sigma\sigma} + B g_{x\sigma\sigma} g^{*}_x + f_{y_{+}} g^{**}_{xx}\left(g^{*}_{x} \otimes g^{*}_{\sigma\sigma}\right) + f_{zz} \left(z_{x} \otimes z_{\sigma\sigma}\right) + F_{xu_{+}u_{+}}\left(I_{n_x} \otimes vec(\Sigma)\right) = 0$ $F_{xu_{+}u_{+}} = f_{y_{+}^{\ast\ast}} g_{xuu}^{\ast\ast} (g_x^{\ast} \otimes I_{n_u^2}) + f_{zz} \left( \left( z_{xu_{+}} \otimes z_{u_{+}} \right)P_{xu\_u} + \left(z_{x} \otimes z_{u_{+}u_{+}}\right) \right) + f_{zzz}\left(z_{x} \otimes z_{u_{+}} \otimes z_{u_{+}}\right)$ $(ghu\sigma\sigma) : A g_{u\sigma\sigma} + B g_{x\sigma\sigma} g^{*}_{u} + f_{y_{+}} g^{**}_{xx}\left(g^{*}_{u} \otimes g^{*}_{\sigma\sigma}\right) + f_{zz} \left(z_{u} \otimes z_{\sigma\sigma}\right) + F_{uu_{+}u_{+}}\left(I_{n_u} \otimes vec(\Sigma_u)\right) = 0$ $F_{uu_{+}u_{+}} = f_{y_{+}^{\ast\ast}} g_{xuu}^{\ast\ast} (g_u^{\ast} \otimes I_{n_u^2}) + f_{zz} \left( \left( z_{uu_{+}} \otimes z_{u_{+}} \right)P_{uu\_u} + \left(z_{u} \otimes z_{u_{+}u_{+}}\right) \right) + f_{zzz}\left(z_{u} \otimes z_{u_{+}} \otimes z_{u_{+}}\right)$ A and B are the common perturbation matrices: $A = f_{y_0} + \begin{pmatrix} \underbrace{0}_{n\times n_{static}} &\vdots& \underbrace{f_{y^{**}_{+}} \cdot g^{**}_{x}}_{n \times n_{spred}} &\vdots& \underbrace{0}_{n\times n_{frwd}} \end{pmatrix}$and $B = \begin{pmatrix} \underbrace{0}_{n \times n_{static}}&\vdots & \underbrace{0}_{n \times n_{pred}} & \vdots & \underbrace{f_{y^{**}_{+}}}_{n \times n_{sfwrd}} \end{pmatrix}$ and $z=(y_{-}^{\ast}; y; y_{+}^{\ast\ast}; u)$ denotes the dynamic model variables as in `M_.lead_lag_incidence`, $y^\ast$ denote state variables, $y^{\ast\ast}$ denote forward looking variables, $y_+$ denote the variables with a lead, $y_{-}$ denote variables with a lag, $y_0$ denote variables at period t, $f$ the model equations, and $f_z$ the first-order dynamic model derivatives, $f_{zz}$ the second-order dynamic derivatives, and $f_{zzz}$ the third-order dynamic model derivatives. Then: $z_{x} = \begin{pmatrix}I\\g_{x}\\g^{**}_{x} g^{*}_{x}\\0\end{pmatrix}$, $z_{u} =\begin{pmatrix}0\\g_{u}\\g^{**}_{x} \cdot g^{*}_{u}\\I\end{pmatrix}$, $z_{u_{+}} =\begin{pmatrix}0\\0\\g^{**}_{u}\\0\end{pmatrix}$ $z_{xx} = \begin{pmatrix} 0\\g_{xx}\\g^{**}_{x} \left( g^{*}_x \otimes g^{*}_{x} \right) + g^{**}_{x} g^{*}_{x}\\0\end{pmatrix}$, $z_{xu} =\begin{pmatrix}0\\g_{xu}\\g^{**}_{xx} \left( g^{*}_x \otimes g^{*}_{u} \right) + g^{**}_{x} g^{*}_{xu}\\0\end{pmatrix}$, $z_{uu} =\begin{pmatrix}0\\g_{uu}\\g^{**}_{xx} \left( g^{*}_u \otimes g^{*}_{u} \right) + g^{**}_{x} g^{*}_{uu}\\0\end{pmatrix}$, $z_{xu_{+}} =\begin{pmatrix}0\\0\\g^{**}_{xu} \left( g^{*}_x \otimes I \right)\\0\end{pmatrix}$, $z_{uu_{+}} =\begin{pmatrix}0\\0\\g^{**}_{xu} \left( g^{*}_{u} \otimes I \right)\\0\end{pmatrix}$, $z_{u_{+}u_{+}} =\begin{pmatrix}0\\0\\g^{\ast\ast}_{uu}\\0\end{pmatrix}$, $z_{\sigma\sigma} = \begin{pmatrix}0\\ g_{\sigma\sigma}\\ g^{\ast\ast}_{x}g^{\ast}_{\sigma\sigma} + g^{\ast\ast}_{\sigma\sigma}\\0 \end{pmatrix}$ $P$ are permutation matrices that can be computed using Matlab's `ipermute` function. ## Parameter derivatives of perturbation solutions First, we need the parameter derivatives of first, second, third, and fourth derivatives of the dynamic model (i.e. g1,g2,g3,g4 in dynamic files), I make use of the implicit function theorem: Let $f_{z^k}$ denote the kth derivative (wrt all dynamic variables) of the dynamic model, then let $df_{z^k}$ denote the first-derivative (wrt all model parameters) of $f_{z^k}$ evaluated at the steady state. Note that $f_{z^k}$ is a function of both the model parameters $\theta$ and of the steady state of all dynamic variables $\bar{z}$, which also depend on the parameters. Hence, implicitly $f_{z^k}=f_{z^k}(\theta,\bar{z}(\theta))$ and $df_{z^k}$ consists of two parts: 1. direct derivative wrt to all model parameters given by the preprocessor in the `_params_derivs.m` files 2. contribution of derivative of steady state of dynamic variables (wrt all model parameters): $f_{z^{k+1}} \cdot d\bar{z}$ Note that we already have functionality to compute $d\bar{z}$ analytically. Having this, the above perturbation systems are basically equations of the following types $AX +BXC = RHS$ or $AX = RHS$ Now when taking the derivative (wrt to one single parameter $\theta_j$), we get $A\mathrm{d}\{X\} + B\mathrm{d}\{X\}C = \mathrm{d}\{RHS\} - \mathrm{d}\{A\}X - \mathrm{d}\{B\}XC - BX\mathrm{d}\{C\}$ or $A\mathrm{d}\{X\} = \mathrm{d}\{RHS\} - \mathrm{d}\{A\}X$ The first one is a Sylvester type equation, the second one can be solved by taking the inverse of $A$. The only diffculty and tedious work arrises in computing (the highly sparse) derivatives of $RHS$. *** # New functions: ` ## get_perturbation_params_derivs.m`and `get_perturbation_params_derivs_numerical_objective.m` * The parameter derivatives up to third order are computed in the new function`get_perturbation_params_derivs.m` both analytically and numerically. For numerical derivatives `get_perturbation_params_derivs_numerical_objective.m` is the objective for `fjaco.m` or `hessian_sparse.m` or `hessian.m`. * `get_perturbation_params_derivs.m` is basically an extended version of the previous `get_first_order_solution_params_derivs.m` function. * * `get_perturbation_params_derivs_numerical_objective.m`builds upon `identification_numerical_objective.m`. It is used for numerical derivatives, whenever `analytic_derivation_mode=-1|-2`. It takes from `identification_numerical_objective.m` the parts that compute numerical parameter Jacobians of steady state, dynamic model equations, and perturbation solution matrices. Hence, these parts are removed in `identification_numerical_objective.m` and it only computes numerical parameter Jacobian of moments and spectrum which are needed for identification analysis in `get_identification_jacobians.m`, when `analytic_derivation_mode=-1` only. * Detailed changes: * Most important: notation of this function is now in accordance to the k_order_solver, i.e. we do not compute derivatives of Kalman transition matrices A and B, but rather the solution matrices ghx,ghu,ghxx,ghxu,ghuu,ghs2,ghxxx,ghxxu,ghxuu,ghuuu,ghxss,ghuss in accordance with notation used in `oo_.dr`. As a byproduct at first-order, focusing on ghx and ghu instead of Kalman transition matrices A and B makes the computations slightly faster for large models (e.g. for Quest the computations were faster by a couple of seconds, not much, but okay). * Removed use of `kstate`, see also Dynare/dynare#1653 and Dynare/dynare!1656 * Output arguments are stored in a structure `DERIVS`, there is also a flag `d2flag` that computes parameter hessians needed only in `dsge_likelihood.m`. * Removed `kronflag` as input. `options_.analytic_derivation_mode` is now used instead of `kronflag`. * Removed `indvar`, an index that was used to selected specific variables in the derivatives. This is not needed, as we always compute the parameter derivatives for all variables first and then select a subset of variables. The selection now takes place in other functions, like `dsge_likelihood.m`. * Introduced some checks: (i) deterministic exogenous variables are not supported, (ii) Kronecker method only compatible with first-order approximation so reset to sylvester method, (iii) for purely backward or forward models we need to be careful with the rows in `M_.lead_la g_incidence`, (iv) if `_params_derivs.m` files are missing an error is thrown. * For numerical derivatives, if mod file does not contain an `estimated_params_block`, a temporary one with the most important parameter information is created. ## `unfold_g4.m` * When evaluating g3 and g4 one needs to take into account that these do not contain symmetric elements, so one needs to use `unfold_g3.m` and the new function `unfold_g4.m`. This returns an unfolded version of the same matrix (i.e. with symmetric elements). *** # New test models `.gitignore` and `Makefile.am` are changed accordingly. Also now it is possible to run test suite on analytic_derivatives, i.e. run `make check m/analytic_derivatives` ## `analytic_derivatives/BrockMirman_PertParamsDerivs.mod` * This is the Brock Mirman model, where we know the exact policy function $g$ for capital and consumption. As this does not imply a nonzero $g_{\sigma\sigma}$, $g_{x\sigma\sigma}$, $g_{u\sigma\sigma}$ I added some artificial equations to get nonzero solution matrices with respect to $\sigma$. The true perturbation solution matrices $g_x$ , $g_u$, $g_{xx}$, $g_{xu}$, $g_{uu}$, $g_{\sigma\sigma}$, $g_{xxx}$, $g_{xxu}$, $g_{xuu}$, $g_{uuu}$, $g_{x\sigma\sigma}$, $g_{u\sigma\sigma}$ are then computed analytically with Matlab's symbolic toolbox and saved in `nBrockMirmanSYM.mat`. There is a preprocessor flag that recreates these analytical computations if changes are needed (and to check whether I made some errors here ;-) ) * Then solution matrices up to third order and their parameter Jacobians are then compared to the ones computed by Dynare's `k_order_solver` and by `get_perturbation_params_derivs` for all `analytic_derivation_mode`'s. There will be an error if the maximum absolute deviation is too large, i.e. for numerical derivatives (`analytic_derivation_mode=-1|-2`) the tolerance is choosen lower (around 1e-5); for analytical methods we are stricter: around 1e-13 for first-order, 1e-12 for second order, and 1e-11 for third-order. * As a side note, this mod file also checks Dynare's `k_order_solver` algorithm and throws an error if something is wrong. * This test model shows that the new functionality works well. And analytical derivatives perform way better and accurate than numerical ones, even for this small model. ## `analytic_derivatives/burnside_3_order_PertParamsDerivs.mod` * This builds upon `tests/k_order_perturbation/burnside_k_order.mod` and computes the true parameter derivatives analytically by hand. * This test model also shows that the new functionality works well. ## `analytic_derivatives/LindeTrabandt2019.mod` * Shows that the new functionality also works for medium-sized models, i.e. a SW type model solved at third order with 35 variables (11 states). 2 shocks and 20 parameters. * This mod file can be used to tweak the speed of the computations in the future. * Compares numerical versus analytical parameter derivatives (for first, second and third order). Note that this model clearly shows that numerical ones are quite different than analytical ones even at first order! ## `identification/LindeTrabandt2019_xfail.mod` * This model is a check for issue Dynare/dynare#1595, see fjaco.m below, and will fail. * Removed `analytic_derivatives/ls2003.mod` as this mod file is neither in the testsuite nor does it work. *** # Detailed changes in other functions ## `get_first_order_solution_params_derivs.m` * Deleted, or actually, renamed to `get_perturbation_params_derivs.m`, as this function now is able to compute the derivatives up to third order ## `identification_numerical_objective.m` * `get_perturbation_params_derivs_numerical_objective.m`builds upon `identification_numerical_objective.m`. It takes from `identification_numerical_objective.m` the parts that compute numerical parameter Jacobians of steady state, dynamic model equations, and perturbation solution matrices. Hence, these parts are removed in `identification_numerical_objective.m` and it only computes numerical parameter Jacobian of moments and spectrum which are needed for identification analysis in `get_identification_jacobians.m`, when `analytic_derivation_mode=-1` only. ## `dsge_likelihood.m` * As `get_first_order_solution_params_derivs.m`is renamed to `get_perturbation_params_derivs.m`, the call is adapted. That is,`get_perturbation_params_derivs` does not compute the derivatives of the Kalman transition `T`matrix anymore, but instead of the dynare solution matrix `ghx`. So we recreate `T` here as this amounts to adding some zeros and focusing on selected variables only. * Added some checks to make sure the first-order approximation is selected. * Removed `kron_flag` as input, as `get_perturbation_params_derivs` looks into `options_.analytic_derivation_mode` for `kron_flag`. ## `dynare_identification.m` * make sure that setting `analytic_derivation_mode` is set both in `options_ident` and `options_`. Note that at the end of the function we restore the `options_` structure, so all changes are local. In a next merge request, I will remove the global variables to make all variables local. ## `get_identification_jacobians.m` * As `get_first_order_solution_params_derivs.m`is renamed to `get_perturbation_params_derivs.m`, the call is adapted. That is,`get_perturbation_params_derivs` does not compute the derivatives of the Kalman transition `A` and `B` matrix anymore, but instead of the dynare solution matrix `ghx` and `ghu`. So we recreate these matrices here instead of in `get_perturbation_params_derivs.m`. * Added `str2func` for better function handles in `fjaco.m`. ## `fjaco.m` * make `tol`an option, which can be adjusted by changing `options_.dynatol.x`for identification and parameter derivatives purposes. * include a check and an informative error message, if numerical derivatives (two-sided finite difference method) yield errors in `resol.m` for identification and parameter derivatives purposes. This closes issue Dynare/dynare#1595. * Changed year of copyright to 2010-2017,2019 *** # Further suggestions and questions * Ones this is merged, I will merge request an improvement of the identification toolbox, which will work up to third order using the pruned state space. This will also remove some issues and bugs, and also I will remove global variables in this request. * The third-order derivatives can be further improved by taking sparsity into account and use mex versions for kronecker products etc. I leave this for further testing (and if anybody actually uses this ;-) )
2019-12-17 19:17:09 +01:00
DT = DT(iv,iv,:);
DOm = DERIVS.dOm(iv,iv,:);
DYss = DERIVS.dYss(iv,:);
options_.order = old_order; %make sure order is reset (not sure if necessary)
options_.analytic_derivation_mode = old_analytic_derivation_mode;%make sure analytic_derivation_mode is reset (not sure if necessary)
else
DT = derivatives_info.DT(iv,iv,:);
DOm = derivatives_info.DOm(iv,iv,:);
DYss = derivatives_info.DYss(iv,:);
2017-05-16 12:42:01 +02:00
if isfield(derivatives_info,'full_Hess')
full_Hess = derivatives_info.full_Hess;
end
2017-05-16 12:42:01 +02:00
if full_Hess
2017-05-16 15:10:20 +02:00
D2T = derivatives_info.D2T;
D2Om = derivatives_info.D2Om;
D2Yss = derivatives_info.D2Yss;
end
2017-05-16 12:42:01 +02:00
if isfield(derivatives_info,'no_DLIK')
no_DLIK = derivatives_info.no_DLIK;
end
clear('derivatives_info');
end
DYss = [zeros(size(DYss,1),offset) DYss];
DH=zeros([length(H),length(H),length(xparam1)]);
DQ=zeros([size(Q),length(xparam1)]);
DP=zeros([size(T),length(xparam1)]);
2017-05-16 12:42:01 +02:00
if full_Hess
for j=1:size(D2Yss,1)
2017-05-16 15:10:20 +02:00
tmp(j,:,:) = blkdiag(zeros(offset,offset), squeeze(D2Yss(j,:,:)));
end
D2Yss = tmp;
D2H=sparse(size(D2Om,1),size(D2Om,2)); %zeros([size(H),length(xparam1),length(xparam1)]);
D2P=sparse(size(D2Om,1),size(D2Om,2)); %zeros([size(T),length(xparam1),length(xparam1)]);
jcount=0;
end
if options_.lik_init==1
for i=1:estim_params_.nvx
k =estim_params_.var_exo(i,1);
2017-05-16 15:10:20 +02:00
DQ(k,k,i) = 2*sqrt(Q(k,k));
dum = lyapunov_symm(T,DOm(:,:,i),options_.lyapunov_fixed_point_tol,options_.qz_criterium,options_.lyapunov_complex_threshold,[],options_.debug);
2017-05-16 15:10:20 +02:00
% kk = find(abs(dum) < 1e-12);
% dum(kk) = 0;
DP(:,:,i)=dum;
if full_Hess
for j=1:i
jcount=jcount+1;
dum = lyapunov_symm(T,dyn_unvech(D2Om(:,jcount)),options_.lyapunov_fixed_point_tol,options_.qz_criterium,options_.lyapunov_complex_threshold,[],options_.debug);
2017-05-16 15:10:20 +02:00
% kk = (abs(dum) < 1e-12);
% dum(kk) = 0;
D2P(:,jcount)=dyn_vech(dum);
% D2P(:,:,j,i)=dum;
end
end
end
end
offset = estim_params_.nvx;
for i=1:estim_params_.nvn
k = estim_params_.var_endo(i,1);
DH(k,k,i+offset) = 2*sqrt(H(k,k));
if full_Hess
2017-05-16 15:10:20 +02:00
D2H(k,k,i+offset,i+offset) = 2;
end
end
offset = offset + estim_params_.nvn;
if options_.lik_init==1
for j=1:estim_params_.np
dum = lyapunov_symm(T,DT(:,:,j+offset)*Pstar*T'+T*Pstar*DT(:,:,j+offset)'+DOm(:,:,j+offset),options_.lyapunov_fixed_point_tol,options_.qz_criterium,options_.lyapunov_complex_threshold,[],options_.debug);
2017-05-16 15:10:20 +02:00
% kk = find(abs(dum) < 1e-12);
% dum(kk) = 0;
DP(:,:,j+offset)=dum;
if full_Hess
DTj = DT(:,:,j+offset);
DPj = dum;
for i=1:j+offset
jcount=jcount+1;
DTi = DT(:,:,i);
DPi = DP(:,:,i);
D2Tij = reshape(D2T(:,jcount),size(T));
D2Omij = dyn_unvech(D2Om(:,jcount));
tmp = D2Tij*Pstar*T' + T*Pstar*D2Tij' + DTi*DPj*T' + DTj*DPi*T' + T*DPj*DTi' + T*DPi*DTj' + DTi*Pstar*DTj' + DTj*Pstar*DTi' + D2Omij;
dum = lyapunov_symm(T,tmp,options_.lyapunov_fixed_point_tol,options_.qz_criterium,options_.lyapunov_complex_threshold,[],options_.debug);
2017-05-16 15:10:20 +02:00
% dum(abs(dum)<1.e-12) = 0;
D2P(:,jcount) = dyn_vech(dum);
% D2P(:,:,j+offset,i) = dum;
end
end
end
end
2017-05-16 12:42:01 +02:00
if analytic_derivation==1
analytic_deriv_info={analytic_derivation,DT,DYss,DOm,DH,DP,asy_Hess};
else
analytic_deriv_info={analytic_derivation,DT,DYss,DOm,DH,DP,D2T,D2Yss,D2Om,D2H,D2P};
2017-05-16 12:42:01 +02:00
clear DT DYss DOm DP D2T D2Yss D2Om D2H D2P
end
else
analytic_deriv_info={0};
end
%------------------------------------------------------------------------------
% 4. Likelihood evaluation
%------------------------------------------------------------------------------
if options_.heteroskedastic_filter
Q=Qvec;
end
singularity_has_been_detected = false;
2016-06-14 17:28:52 +02:00
% First test multivariate filter if specified; potentially abort and use univariate filter instead
if ((kalman_algo==1) || (kalman_algo==3))% Multivariate Kalman Filter
if no_missing_data_flag && ~options_.occbin.likelihood.status
if options_.fast_kalman_filter
if diffuse_periods
%kalman_algo==3 requires no diffuse periods (stationary
%observables) as otherwise FE matrix will not be positive
%definite
fval = Inf;
info(1) = 55;
info(4) = 0.1;
2017-05-16 15:10:20 +02:00
exit_flag = 0;
return
end
[LIK,lik] = kalman_filter_fast(Y,diffuse_periods+1,size(Y,2), ...
a_0_given_tm1,Pstar, ...
kalman_tol, riccati_tol, ...
options_.presample, ...
T,Q,R,H,Z,mm,pp,rr,Zflag,diffuse_periods, ...
analytic_deriv_info{:});
else
if options_.kalman_filter_mex
[LIK,lik] = kalman_filter_mex(Y,a_0_given_tm1,Pstar, ...
kalman_tol, riccati_tol, ...
T,Q,R,Z,Zflag,H,diffuse_periods, ...
options_.presample);
else
[LIK,lik] = kalman_filter(Y,diffuse_periods+1,size(Y,2), ...
a_0_given_tm1,Pstar, ...
kalman_tol, riccati_tol, ...
options_.rescale_prediction_error_covariance, ...
options_.presample, ...
T,Q,R,H,Z,mm,pp,rr,Zflag,diffuse_periods, ...
analytic_deriv_info{:});
end
end
else
[LIK,lik] = missing_observations_kalman_filter(dataset_info.missing.aindex,dataset_info.missing.number_of_observations,dataset_info.missing.no_more_missing_observations,Y,diffuse_periods+1,size(Y,2), ...
a_0_given_tm1, Pstar, ...
kalman_tol, options_.riccati_tol, ...
options_.rescale_prediction_error_covariance, ...
options_.presample, ...
T,Q,R,H,Z,mm,pp,rr,Zflag,diffuse_periods, occbin_);
if occbin_.status && isinf(LIK)
fval = Inf;
info(1) = 320;
info(4) = 0.1;
exit_flag = 0;
return
end
end
2017-05-16 12:42:01 +02:00
if analytic_derivation
LIK1=LIK;
LIK=LIK1{1};
lik1=lik;
lik=lik1{1};
end
if isinf(LIK)
if options_.use_univariate_filters_if_singularity_is_detected
singularity_has_been_detected = true;
if kalman_algo == 1
kalman_algo = 2;
else
kalman_algo = 4;
end
else
2016-06-16 10:40:21 +02:00
fval = Inf;
info(1) = 50;
info(4) = 0.1;
exit_flag = 0;
return
end
else
if options_.lik_init==3
LIK = LIK + dLIK;
2017-05-16 12:42:01 +02:00
if analytic_derivation==0 && nargout>3
if ~singular_diffuse_filter
lik = [dlik; lik];
else
lik = [sum(dlik,2); lik];
end
end
end
end
end
if (kalman_algo==2) || (kalman_algo==4)
% Univariate Kalman Filter
% resetting measurement error covariance matrix when necessary following DK (2012), Section 6.4.3 %
if isequal(H,0)
H1 = zeros(pp,1);
mmm = mm;
2017-05-16 12:42:01 +02:00
if analytic_derivation
DH = zeros(pp,length(xparam1));
end
else
if all(all(abs(H-diag(diag(H)))<1e-14))% ie, the covariance matrix is diagonal...
H1 = diag(H);
mmm = mm;
clear('tmp')
2017-05-16 12:42:01 +02:00
if analytic_derivation
for j=1:pp
tmp(j,:)=DH(j,j,:);
end
DH=tmp;
end
else
if ~expanded_state_vector_for_univariate_filter
Z1=zeros(pp,size(T,2));
for jz=1:length(Z)
Z1(jz,Z(jz))=1;
end
Z = [Z1, eye(pp)];
Zflag=1;
T = blkdiag(T,zeros(pp));
if options_.heteroskedastic_filter
clear Q
for kv=1:size(Qvec,3)
Q(:,:,kv) = blkdiag(Qvec(:,:,kv),H);
end
else
Q = blkdiag(Q,H);
end
R = blkdiag(R,eye(pp));
Pstar = blkdiag(Pstar,H);
Pinf = blkdiag(Pinf,zeros(pp));
H1 = zeros(pp,1);
Zflag=1;
end
mmm = mm+pp;
if singularity_has_been_detected
a_tmp = zeros(mmm,1);
a_tmp(1:length(a_0_given_tm1)) = a_0_given_tm1;
a_0_given_tm1 = a_tmp;
elseif ~expanded_state_vector_for_univariate_filter
a_0_given_tm1 = [a_0_given_tm1; zeros(pp,1)];
end
end
end
2017-05-16 12:42:01 +02:00
if analytic_derivation
analytic_deriv_info{5}=DH;
end
[LIK, lik] = univariate_kalman_filter(dataset_info.missing.aindex,dataset_info.missing.number_of_observations,dataset_info.missing.no_more_missing_observations,Y,diffuse_periods+1,size(Y,2), ...
a_0_given_tm1,Pstar, ...
options_.kalman_tol, ...
options_.riccati_tol, ...
options_.presample, ...
T,Q,R,H1,Z,mmm,pp,rr,Zflag,diffuse_periods,analytic_deriv_info{:});
2017-05-16 12:42:01 +02:00
if analytic_derivation
LIK1=LIK;
LIK=LIK1{1};
lik1=lik;
lik=lik1{1};
end
if options_.lik_init==3
LIK = LIK+dLIK;
2017-05-16 12:42:01 +02:00
if analytic_derivation==0 && nargout>3
lik = [dlik; lik];
end
end
end
if analytic_derivation
if no_DLIK==0
DLIK = LIK1{2};
end
2017-05-16 12:42:01 +02:00
if full_Hess
Hess = -LIK1{3};
end
2017-05-16 12:42:01 +02:00
if asy_Hess
Hess = LIK1{3};
end
end
if isnan(LIK)
2016-06-01 18:22:51 +02:00
fval = Inf;
info(1) = 45;
info(4) = 0.1;
exit_flag = 0;
return
end
if imag(LIK)~=0
2016-06-01 18:22:51 +02:00
fval = Inf;
info(1) = 46;
info(4) = 0.1;
exit_flag = 0;
return
end
2016-06-14 17:28:52 +02:00
if isinf(LIK)~=0
fval = Inf;
info(1) = 50;
info(4) = 0.1;
exit_flag = 0;
return
end
likelihood = LIK;
% ------------------------------------------------------------------------------
% 5. Adds prior if necessary
% ------------------------------------------------------------------------------
if analytic_derivation
2017-05-16 12:42:01 +02:00
if full_Hess
[lnprior, dlnprior, d2lnprior] = priordens(xparam1,bayestopt_.pshape,bayestopt_.p6,bayestopt_.p7,bayestopt_.p3,bayestopt_.p4);
Hess = Hess - d2lnprior;
else
[lnprior, dlnprior] = priordens(xparam1,bayestopt_.pshape,bayestopt_.p6,bayestopt_.p7,bayestopt_.p3,bayestopt_.p4);
end
if no_DLIK==0
DLIK = DLIK - dlnprior';
end
2017-05-16 12:42:01 +02:00
if outer_product_gradient
dlik = lik1{2};
dlik=[- dlnprior; dlik(start:end,:)];
Hess = dlik'*dlik;
end
else
lnprior = priordens(xparam1,bayestopt_.pshape,bayestopt_.p6,bayestopt_.p7,bayestopt_.p3,bayestopt_.p4);
end
if options_.endogenous_prior==1
if options_.lik_init==2 || options_.lik_init==3
2017-05-16 15:10:20 +02:00
error('Endogenous prior not supported with non-stationary models')
else
[lnpriormom] = endogenous_prior(Y,dataset_info,Pstar,bayestopt_,H);
2017-05-16 15:10:20 +02:00
fval = (likelihood-lnprior-lnpriormom);
end
else
2017-05-16 15:10:20 +02:00
fval = (likelihood-lnprior);
end
if options_.prior_restrictions.status
tmp = feval(options_.prior_restrictions.routine, M_, dr, endo_steady_state, exo_steady_state, exo_det_steady_state, options_, dataset_, dataset_info);
fval = fval - tmp;
end
if isnan(fval)
2016-06-01 18:22:51 +02:00
fval = Inf;
info(1) = 47;
info(4) = 0.1;
exit_flag = 0;
return
end
if imag(fval)~=0
2016-06-01 18:22:51 +02:00
fval = Inf;
info(1) = 48;
info(4) = 0.1;
exit_flag = 0;
return
end
if ~options_.kalman.keep_kalman_algo_if_singularity_is_detected
% Update options_.kalman_algo.
options_.kalman_algo = kalman_algo;
end
2017-05-16 12:42:01 +02:00
if analytic_derivation==0 && nargout>3
lik=lik(start:end,:);
DLIK=[-lnprior; lik(:)];
end
function a=set_Kalman_starting_values(a,M_,dr,options_,bayestopt_)
% function a=set_Kalman_starting_values(a,M_,dr,options_,bayestopt_)
% Sets initial states guess for Kalman filter/smoother based on M_.filter_initial_state
%
% INPUTS
% o a [double] (p*1) vector of states
% o M_ [structure] decribing the model
% o dr [structure] storing the decision rules
% o options_ [structure] describing the options
% o bayestopt_ [structure] describing the priors
%
% OUTPUTS
% o a [double] (p*1) vector of set initial states
if isfield(M_,'filter_initial_state') && ~isempty(M_.filter_initial_state)
state_indices=dr.order_var(dr.restrict_var_list(bayestopt_.mf0));
for ii=1:size(state_indices,1)
if ~isempty(M_.filter_initial_state{state_indices(ii),1})
if options_.loglinear && ~options_.logged_steady_state
a(bayestopt_.mf0(ii)) = log(eval(M_.filter_initial_state{state_indices(ii),2})) - log(dr.ys(state_indices(ii)));
elseif ~options_.loglinear && ~options_.logged_steady_state
a(bayestopt_.mf0(ii)) = eval(M_.filter_initial_state{state_indices(ii),2}) - dr.ys(state_indices(ii));
else
error(['The steady state is logged. This should not happen. Please contact the developers'])
end
end
end
2021-02-18 14:26:10 +01:00
end
function occbin_options = set_occbin_options(options_)
% this builds the opts_simul options field needed by occbin.solver
occbin_options.opts_simul = options_.occbin.simul;
occbin_options.opts_simul.curb_retrench = options_.occbin.likelihood.curb_retrench;
occbin_options.opts_simul.maxit = options_.occbin.likelihood.maxit;
occbin_options.opts_simul.periods = options_.occbin.likelihood.periods;
occbin_options.opts_simul.check_ahead_periods = options_.occbin.likelihood.check_ahead_periods;
occbin_options.opts_simul.max_check_ahead_periods = options_.occbin.likelihood.max_check_ahead_periods;
occbin_options.opts_simul.periodic_solution = options_.occbin.likelihood.periodic_solution;
occbin_options.opts_simul.restrict_state_space = options_.occbin.likelihood.restrict_state_space;
occbin_options.opts_simul.full_output = options_.occbin.likelihood.full_output;
occbin_options.opts_simul.piecewise_only = options_.occbin.likelihood.piecewise_only;
if ~isempty(options_.occbin.smoother.init_binding_indicator)
occbin_options.opts_simul.init_binding_indicator = options_.occbin.likelihood.init_binding_indicator;
occbin_options.opts_simul.init_regime_history=options_.occbin.likelihood.init_regime_history;
end