


Computes the likelihood of a stationnary state space model (univariate approach).


0001 function [LIK, lik,a,P] = univariate_kalman_filter(data_index,number_of_observations,no_more_missing_observations,Y,start,last,a,P,kalman_tol,riccati_tol,presample,T,Q,R,H,Z,mm,pp,rr,Zflag,diffuse_periods,analytic_derivation,DT,DYss,DOm,DH,DP,D2T,D2Yss,D2Om,D2H,D2P) 0002 % Computes the likelihood of a stationnary state space model (univariate approach). 0003 0004 %@info: 0005 %! @deftypefn {Function File} {[@var{LIK},@var{likk},@var{a},@var{P} ] =} univariate_kalman_filter (@var{data_index}, @var{number_of_observations},@var{no_more_missing_observations}, @var{Y}, @var{start}, @var{last}, @var{a}, @var{P}, @var{kalman_tol}, @var{riccati_tol},@var{presample},@var{T},@var{Q},@var{R},@var{H},@var{Z},@var{mm},@var{pp},@var{rr},@var{Zflag},@var{diffuse_periods}) 0006 %! @anchor{univariate_kalman_filter} 0007 %! @sp 1 0008 %! Computes the likelihood of a stationary state space model, given initial condition for the states (mean and variance). 0009 %! @sp 2 0010 %! @strong{Inputs} 0011 %! @sp 1 0012 %! @table @ @var 0013 %! @item data_index 0014 %! Matlab's cell, 1*T cell of column vectors of indices (in the vector of observed variables). 0015 %! @item number_of_observations 0016 %! Integer scalar, effective number of observations. 0017 %! @item no_more_missing_observations 0018 %! Integer scalar, date after which there is no more missing observation (it is then possible to switch to the steady state kalman filter). 0019 %! @item Y 0020 %! Matrix (@var{pp}*T) of doubles, data. 0021 %! @item start 0022 %! Integer scalar, first period. 0023 %! @item last 0024 %! Integer scalar, last period (@var{last}-@var{first} has to be inferior to T). 0025 %! @item a 0026 %! Vector (@var{mm}*1) of doubles, initial mean of the state vector. 0027 %! @item P 0028 %! Matrix (@var{mm}*@var{mm}) of doubles, initial covariance matrix of the state vector. 0029 %! @item kalman_tol 0030 %! Double scalar, tolerance parameter (rcond, inversibility of the covariance matrix of the prediction errors). 0031 %! @item riccati_tol 0032 %! Double scalar, tolerance parameter (iteration over the Riccati equation). 0033 %! @item presample 0034 %! Integer scalar, presampling if strictly positive (number of initial iterations to be discarded when evaluating the likelihood). 0035 %! @item T 0036 %! Matrix (@var{mm}*@var{mm}) of doubles, transition matrix of the state equation. 0037 %! @item Q 0038 %! Matrix (@var{rr}*@var{rr}) of doubles, covariance matrix of the structural innovations (noise in the state equation). 0039 %! @item R 0040 %! Matrix (@var{mm}*@var{rr}) of doubles, second matrix of the state equation relating the structural innovations to the state variables. 0041 %! @item H 0042 %! Vector (@var{pp}) of doubles, diagonal of covariance matrix of the measurement errors (corelation among measurement errors is handled by a model transformation). 0043 %! @item Z 0044 %! Matrix (@var{pp}*@var{mm}) of doubles or vector of integers, matrix relating the states to the observed variables or vector of indices (depending on the value of @var{Zflag}). 0045 %! @item mm 0046 %! Integer scalar, number of state variables. 0047 %! @item pp 0048 %! Integer scalar, number of observed variables. 0049 %! @item rr 0050 %! Integer scalar, number of structural innovations. 0051 %! @item Zflag 0052 %! Integer scalar, equal to 0 if Z is a vector of indices targeting the obseved variables in the state vector, equal to 1 if Z is a @var{pp}*@var{mm} matrix. 0053 %! @item diffuse_periods 0054 %! Integer scalar, number of diffuse filter periods in the initialization step. 0055 %! @end table 0056 %! @sp 2 0057 %! @strong{Outputs} 0058 %! @sp 1 0059 %! @table @ @var 0060 %! @item LIK 0061 %! Double scalar, value of (minus) the likelihood. 0062 %! @item likk 0063 %! Column vector of doubles, values of the density of each observation. 0064 %! @item a 0065 %! Vector (@var{mm}*1) of doubles, mean of the state vector at the end of the (sub)sample. 0066 %! @item P 0067 %! Matrix (@var{mm}*@var{mm}) of doubles, covariance of the state vector at the end of the (sub)sample. 0068 %! @end table 0069 %! @sp 2 0070 %! @strong{This function is called by:} 0071 %! @sp 1 0072 %! @ref{dsge_likelihood} 0073 %! @sp 2 0074 %! @strong{This function calls:} 0075 %! @sp 1 0076 %! @ref{univariate_kalman_filter_ss} 0077 %! @end deftypefn 0078 %@eod: 0079 0080 % Copyright (C) 2004-2011 Dynare Team 0081 % 0082 % This file is part of Dynare. 0083 % 0084 % Dynare is free software: you can redistribute it and/or modify 0085 % it under the terms of the GNU General Public License as published by 0086 % the Free Software Foundation, either version 3 of the License, or 0087 % (at your option) any later version. 0088 % 0089 % Dynare is distributed in the hope that it will be useful, 0090 % but WITHOUT ANY WARRANTY; without even the implied warranty of 0091 % MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the 0092 % GNU General Public License for more details. 0093 % 0094 % You should have received a copy of the GNU General Public License 0095 % along with Dynare. If not, see <http://www.gnu.org/licenses/>. 0096 0097 % AUTHOR(S) stephane DOT adjemian AT univ DASH lemans DOT fr 0098 0099 if nargin<20 || isempty(Zflag)% Set default value for Zflag ==> Z is a vector of indices. 0100 Zflag = 0; 0101 diffuse_periods = 0; 0102 end 0103 0104 if nargin<21 0105 diffuse_periods = 0; 0106 end 0107 0108 % Get sample size. 0109 smpl = last-start+1; 0110 0111 % Initialize some variables. 0112 QQ = R*Q*transpose(R); % Variance of R times the vector of structural innovations. 0113 t = start; % Initialization of the time index. 0114 lik = zeros(smpl,1); % Initialization of the vector gathering the densities. 0115 LIK = Inf; % Default value of the log likelihood. 0116 oldP = Inf; 0117 l2pi = log(2*pi); 0118 notsteady = 1; 0119 0120 oldK = Inf; 0121 K = NaN(mm,pp); 0122 0123 if analytic_derivation == 0, 0124 DLIK=[]; 0125 Hess=[]; 0126 else 0127 k = size(DT,3); % number of structural parameters 0128 DLIK = zeros(k,1); % Initialization of the score. 0129 Da = zeros(mm,k); % Derivative State vector. 0130 0131 if Zflag==0, 0132 C = zeros(pp,mm); 0133 for ii=1:pp; C(ii,Z(ii))=1;end % SELECTION MATRIX IN MEASUREMENT EQ. (FOR WHEN IT IS NOT CONSTANT) 0134 else 0135 C=Z; 0136 end 0137 dC = zeros(pp,mm,k); % either selection matrix or schur have zero derivatives 0138 if analytic_derivation==2, 0139 Hess = zeros(k,k); % Initialization of the Hessian 0140 D2a = zeros(mm,k,k); % State vector. 0141 d2C = zeros(pp,mm,k,k); 0142 else 0143 Hess=[]; 0144 D2a=[]; 0145 D2T=[]; 0146 D2Yss=[]; 0147 end 0148 LIK={inf,DLIK,Hess}; 0149 end 0150 0151 while notsteady && t<=last 0152 s = t-start+1; 0153 d_index = data_index{t}; 0154 if Zflag 0155 z = Z(d_index,:); 0156 else 0157 z = Z(d_index); 0158 end 0159 oldP = P(:); 0160 for i=1:rows(z) 0161 if Zflag 0162 prediction_error = Y(d_index(i),t) - z(i,:)*a; 0163 PZ = P*z(i,:)'; 0164 Fi = z(i,:)*PZ + H(d_index(i)); 0165 else 0166 prediction_error = Y(d_index(i),t) - a(z(i)); 0167 PZ = P(:,z(i)); 0168 Fi = PZ(z(i)) + H(d_index(i)); 0169 end 0170 if Fi>kalman_tol 0171 Ki = PZ/Fi; 0172 if t>=no_more_missing_observations 0173 K(:,i) = Ki; 0174 end 0175 lik(s) = lik(s) + log(Fi) + prediction_error*prediction_error/Fi + l2pi; 0176 if analytic_derivation, 0177 if analytic_derivation==2, 0178 [Da,DP,DLIKt,D2a,D2P, Hesst] = univariate_computeDLIK(k,i,z(i,:),Zflag,prediction_error,Ki,PZ,Fi,Da,DYss,DP,DH(d_index(i),:),notsteady,D2a,D2Yss,D2P); 0179 else 0180 [Da,DP,DLIKt] = univariate_computeDLIK(k,i,z(i,:),Zflag,prediction_error,Ki,PZ,Fi,Da,DYss,DP,DH(d_index(i),:),notsteady); 0181 end 0182 if t>presample 0183 DLIK = DLIK + DLIKt; 0184 if analytic_derivation==2, 0185 Hess = Hess + Hesst; 0186 end 0187 end 0188 end 0189 a = a + Ki*prediction_error; 0190 P = P - PZ*Ki'; 0191 end 0192 end 0193 if analytic_derivation, 0194 if analytic_derivation==2, 0195 [Da,DP,D2a,D2P] = univariate_computeDstate(k,a,P,T,Da,DP,DT,DOm,notsteady,D2a,D2P,D2T,D2Om); 0196 else 0197 [Da,DP] = univariate_computeDstate(k,a,P,T,Da,DP,DT,DOm,notsteady); 0198 end 0199 end 0200 a = T*a; 0201 P = T*P*T' + QQ; 0202 if t>=no_more_missing_observations 0203 notsteady = max(abs(K(:)-oldK))>riccati_tol; 0204 oldK = K(:); 0205 end 0206 t = t+1; 0207 end 0208 0209 % Divide by two. 0210 lik(1:s) = .5*lik(1:s); 0211 0212 % Call steady state univariate kalman filter if needed. 0213 if t <= last 0214 if analytic_derivation, 0215 if analytic_derivation==2, 0216 [tmp, lik(s+1:end)] = univariate_kalman_filter_ss(Y,t,last,a,P,kalman_tol,T,H,Z,pp,Zflag, ... 0217 analytic_derivation,Da,DT,DYss,DP,DH,D2a,D2T,D2Yss,D2P); 0218 else 0219 [tmp, lik(s+1:end)] = univariate_kalman_filter_ss(Y,t,last,a,P,kalman_tol,T,H,Z,pp,Zflag, ... 0220 analytic_derivation,Da,DT,DYss,DP,DH); 0221 end 0222 DLIK = DLIK + tmp{2}; 0223 if analytic_derivation==2, 0224 Hess = Hess + tmp{3}; 0225 end 0226 else 0227 [tmp, lik(s+1:end)] = univariate_kalman_filter_ss(Y,t,last,a,P,kalman_tol,T,H,Z,pp,Zflag); 0228 end 0229 end 0230 0231 % Compute minus the log-likelihood. 0232 if presample > diffuse_periods 0233 LIK = sum(lik(1+presample-diffuse_periods:end)); 0234 else 0235 LIK = sum(lik); 0236 end 0237 0238 if analytic_derivation, 0239 DLIK = DLIK/2; 0240 if analytic_derivation==2, 0241 Hess = -Hess/2; 0242 LIK={LIK, DLIK, Hess}; 0243 else 0244 LIK={LIK, DLIK}; 0245 end 0246 end