Home > . > DiffuseLikelihood3.m

DiffuseLikelihood3

PURPOSE ^

M. Ratto added lik in output [October 2005]

SYNOPSIS ^

function [LIK, lik] = DiffuseLikelihood3(T,R,Q,Pinf,Pstar,Y,trend,start)%//Z,T,R,Q,Pinf,Pstar,Y)

DESCRIPTION ^

 M. Ratto added lik in output [October 2005]
 changes by M. Ratto [April 2005]
 introduced new options options_.diffuse_d  for termination of DKF
 new icc counter for Finf steps in DKF
 new termination for DKF
 likelihood terms for Fstar must be cumulated in DKF also when Pinf is non
 zero. 

 [4/5/2005] correctyed bug in the modified verson of Ratto for rank of Pinf 
 introduced a specific crit1 for the DKF termination

 stepane.adjemian@cepremap.cnrs.fr [07-19-2004]
 
   See "Filtering and Smoothing of State Vector for Diffuse State Space
   Models", S.J. Koopman and J. Durbin (2003, in Journal of Time Series 
   Analysis, vol. 24(1), pp. 85-98).  

    Case where F_{\infty,t} is singular ==> Univariate treatment of multivariate
    time series.

   THE PROBLEM:
 
   y_t =   Z_t * \alpha_t + \varepsilon_t
   \alpha_{t+1} = T_t  * \alpha_t + R_t * \eta_t
 
   with:
 
   \alpha_1 = a + A*\delta + R_0*\eta_0

   m*q matrix A and m*(m-q) matrix R_0 are selection matrices (their
   columns constitue all the columns of the m*m identity matrix) so that 

       A'*R_0 = 0 and A'*\alpha_1 = \delta

   We assume that the vector \delta is distributed as a N(0,\kappa*I_q)
   for a given  \kappa > 0. So that the expectation of \alpha_1 is a and
   its variance is P, with

       P = \kappa*P_{\infty} + P_{\star}

           P_{\infty} = A*A'
           P_{\star}  = R_0*Q_0*R_0'

   P_{\infty} is a m*m diagonal matrix with q ones and m-q zeros. 

 
   and where:
 
   y_t             is a pp*1 vector
   \alpha_t        is a mm*1 vector
   \varepsilon_t   is a pp*1 multivariate random variable (iid N(0,H_t))
   \eta_t          is a rr*1 multivariate random variable (iid N(0,Q_t))
   a_1             is a mm*1 vector
 
   Z_t     is a pp*mm matrix
   T_t     is a mm*mm matrix
   H_t     is a pp*pp matrix
   R_t     is a mm*rr matrix
   Q_t     is a rr*rr matrix
   P_1     is a mm*mm matrix
 
 
   FILTERING EQUATIONS:
 
   v_t = y_t - Z_t* a_t
   F_t = Z_t * P_t * Z_t' + H_t
   K_t = T_t * P_t * Z_t' * F_t^{-1}
   L_t = T_t - K_t * Z_t
   a_{t+1} = T_t * a_t + K_t * v_t
   P_{t+1} = T_t * P_t * L_t' + R_t*Q_t*R_t'
 

   DIFFUSE FILTERING EQUATIONS:

   a_{t+1} = T_t*a_t + K_{\ast,t}v_t
   P_{\infty,t+1} = T_t*P_{\infty,t}*T_t'
   P_{\ast,t+1}  = T_t*P_{\ast,t}*L_{\ast,t}' + R_t*Q_t*R_t'
   K_{\ast,t}   = T_t*P_{\ast,t}*Z_t'*F_{\ast,t}^{-1}
   v_t = y_t - Z_t*a_t
   L_{\ast,t} = T_t - K_{\ast,t}*Z_t
   F_{\ast,t}  = Z_t*P_{\ast,t}*Z_t' + H_t

CROSS-REFERENCE INFORMATION ^

This function calls: This function is called by:

SOURCE CODE ^

0001 function [LIK, lik] = DiffuseLikelihood3(T,R,Q,Pinf,Pstar,Y,trend,start)%//Z,T,R,Q,Pinf,Pstar,Y)
0002 % M. Ratto added lik in output [October 2005]
0003 % changes by M. Ratto [April 2005]
0004 % introduced new options options_.diffuse_d  for termination of DKF
0005 % new icc counter for Finf steps in DKF
0006 % new termination for DKF
0007 % likelihood terms for Fstar must be cumulated in DKF also when Pinf is non
0008 % zero.
0009 %
0010 % [4/5/2005] correctyed bug in the modified verson of Ratto for rank of Pinf
0011 % introduced a specific crit1 for the DKF termination
0012 %
0013 % stepane.adjemian@cepremap.cnrs.fr [07-19-2004]
0014 %
0015 %   See "Filtering and Smoothing of State Vector for Diffuse State Space
0016 %   Models", S.J. Koopman and J. Durbin (2003, in Journal of Time Series
0017 %   Analysis, vol. 24(1), pp. 85-98).
0018 %
0019 %    Case where F_{\infty,t} is singular ==> Univariate treatment of multivariate
0020 %    time series.
0021 %
0022 %   THE PROBLEM:
0023 %
0024 %   y_t =   Z_t * \alpha_t + \varepsilon_t
0025 %   \alpha_{t+1} = T_t  * \alpha_t + R_t * \eta_t
0026 %
0027 %   with:
0028 %
0029 %   \alpha_1 = a + A*\delta + R_0*\eta_0
0030 %
0031 %   m*q matrix A and m*(m-q) matrix R_0 are selection matrices (their
0032 %   columns constitue all the columns of the m*m identity matrix) so that
0033 %
0034 %       A'*R_0 = 0 and A'*\alpha_1 = \delta
0035 %
0036 %   We assume that the vector \delta is distributed as a N(0,\kappa*I_q)
0037 %   for a given  \kappa > 0. So that the expectation of \alpha_1 is a and
0038 %   its variance is P, with
0039 %
0040 %       P = \kappa*P_{\infty} + P_{\star}
0041 %
0042 %           P_{\infty} = A*A'
0043 %           P_{\star}  = R_0*Q_0*R_0'
0044 %
0045 %   P_{\infty} is a m*m diagonal matrix with q ones and m-q zeros.
0046 %
0047 %
0048 %   and where:
0049 %
0050 %   y_t             is a pp*1 vector
0051 %   \alpha_t        is a mm*1 vector
0052 %   \varepsilon_t   is a pp*1 multivariate random variable (iid N(0,H_t))
0053 %   \eta_t          is a rr*1 multivariate random variable (iid N(0,Q_t))
0054 %   a_1             is a mm*1 vector
0055 %
0056 %   Z_t     is a pp*mm matrix
0057 %   T_t     is a mm*mm matrix
0058 %   H_t     is a pp*pp matrix
0059 %   R_t     is a mm*rr matrix
0060 %   Q_t     is a rr*rr matrix
0061 %   P_1     is a mm*mm matrix
0062 %
0063 %
0064 %   FILTERING EQUATIONS:
0065 %
0066 %   v_t = y_t - Z_t* a_t
0067 %   F_t = Z_t * P_t * Z_t' + H_t
0068 %   K_t = T_t * P_t * Z_t' * F_t^{-1}
0069 %   L_t = T_t - K_t * Z_t
0070 %   a_{t+1} = T_t * a_t + K_t * v_t
0071 %   P_{t+1} = T_t * P_t * L_t' + R_t*Q_t*R_t'
0072 %
0073 %
0074 %   DIFFUSE FILTERING EQUATIONS:
0075 %
0076 %   a_{t+1} = T_t*a_t + K_{\ast,t}v_t
0077 %   P_{\infty,t+1} = T_t*P_{\infty,t}*T_t'
0078 %   P_{\ast,t+1}  = T_t*P_{\ast,t}*L_{\ast,t}' + R_t*Q_t*R_t'
0079 %   K_{\ast,t}   = T_t*P_{\ast,t}*Z_t'*F_{\ast,t}^{-1}
0080 %   v_t = y_t - Z_t*a_t
0081 %   L_{\ast,t} = T_t - K_{\ast,t}*Z_t
0082 %   F_{\ast,t}  = Z_t*P_{\ast,t}*Z_t' + H_t
0083 global bayestopt_ options_
0084 
0085 mf = bayestopt_.mf;
0086 pp     = size(Y,1);
0087 mm     = size(T,1);
0088 smpl   = size(Y,2);
0089 a      = zeros(mm,1);
0090 QQ     = R*Q*transpose(R);
0091 t      = 0;
0092 lik       = zeros(smpl+1,1);    
0093 lik(smpl+1) = smpl*pp*log(2*pi);        %% the constant of minus two times the log-likelihood
0094 notsteady     = 1;
0095 crit          = options_.kalman_tol;
0096 crit1          = 1.e-6;
0097 newRank         = rank(Pinf,crit1);
0098 icc=0;
0099 while newRank & t < smpl
0100   t = t+1;
0101   for i=1:pp
0102     v(i)     = Y(i,t)-a(mf(i))-trend(i,t);
0103     Fstar     = Pstar(mf(i),mf(i));
0104     Finf    = Pinf(mf(i),mf(i));
0105     Kstar     = Pstar(:,mf(i));
0106     if Finf > crit & newRank,  %added newRank criterion
0107       icc=icc+1;
0108       Kinf    = Pinf(:,mf(i));
0109       a        = a + Kinf*v(i)/Finf;
0110       Pstar    = Pstar + Kinf*transpose(Kinf)*Fstar/(Finf*Finf) - ...
0111       (Kstar*transpose(Kinf)+Kinf*transpose(Kstar))/Finf;
0112       Pinf    = Pinf - Kinf*transpose(Kinf)/Finf;
0113       lik(t)     = lik(t) + log(Finf);
0114       % start new termination criterion for DKF
0115       if ~isempty(options_.diffuse_d),  
0116     newRank = (icc<options_.diffuse_d);  
0117     %if newRank & any(diag(Pinf(mf,mf))>crit)==0; %  M. Ratto this line is BUGGY
0118     if newRank & (any(diag(Pinf(mf,mf))>crit)==0 & rank(Pinf,crit1)==0); 
0119       options_.diffuse_d = icc;
0120       newRank=0;
0121       disp('WARNING: Change in OPTIONS_.DIFFUSE_D in univariate DKF')
0122       disp(['new OPTIONS_.DIFFUSE_D = ',int2str(icc)])
0123       disp('You may have to reset the optimisation')
0124     end
0125       else
0126     %newRank = any(diag(Pinf(mf,mf))>crit);     % M. Ratto this line is BUGGY
0127     newRank = (any(diag(Pinf(mf,mf))>crit) | rank(Pinf,crit1));                 
0128     if newRank==0, 
0129       P0=    T*Pinf*transpose(T);
0130       %newRank = any(diag(P0(mf,mf))>crit);   % M. Ratto this line is BUGGY
0131       newRank = (any(diag(P0(mf,mf))>crit) | rank(P0,crit1));   % M. Ratto 11/10/2005
0132       if newRank==0, 
0133         options_.diffuse_d = icc;
0134       end
0135     end                    
0136       end,
0137       % end new termination and checks for DKF
0138     elseif Fstar > crit 
0139       %% Note that : (1) rank(Pinf)=0 implies that Finf = 0, (2) outside this loop (when for some i and t the condition
0140       %% rank(Pinf)=0 is satisfied we have P = Pstar and F = Fstar and (3) Finf = 0 does not imply that
0141       %% rank(Pinf)=0. [stphane,11-03-2004].
0142       %if rank(Pinf,crit) == 0
0143       % the likelihood terms should alwasy be cumulated, not only
0144       % when Pinf=0, otherwise the lik would depend on the ordering
0145       % of observed variables
0146       % presample options can be used to ignore initial time points
0147       lik(t) = lik(t) + log(Fstar) + v(i)*v(i)/Fstar;
0148       %end
0149       a    = a + Kstar*v(i)/Fstar;
0150       Pstar = Pstar - Kstar*transpose(Kstar)/Fstar;
0151     else
0152       %disp(['zero F term in DKF for observed ',int2str(i),' ',num2str(Fstar)])
0153     end
0154   end 
0155 %     if all(abs(Pinf(:))<crit),
0156 %         oldRank = 0;
0157 %     else
0158 %         oldRank = rank(Pinf,crit);
0159 %     end
0160     if newRank,
0161         oldRank = rank(Pinf,crit1);
0162     else
0163         oldRank = 0;
0164     end
0165     a         = T*a;
0166     Pstar     = T*Pstar*transpose(T)+QQ;
0167     Pinf    = T*Pinf*transpose(T);
0168 %     if all(abs(Pinf(:))<crit),
0169 %         newRank = 0;
0170 %     else
0171 %         newRank = rank(Pinf,crit);
0172 %     end
0173     if newRank,
0174         newRank = rank(Pinf,crit1);  % new crit1 is used
0175     end
0176     if oldRank ~= newRank
0177         disp('DiffuseLiklihood3 :: T does influence the rank of Pinf!')    
0178     end  
0179 end
0180 if t == smpl                                                           
0181   error(['There isn''t enough information to estimate the initial' ... 
0182      ' conditions of the nonstationary variables']);                   
0183 end   
0184 while notsteady & t < smpl
0185   t = t+1;
0186   oldP = Pstar;
0187   for i=1:pp
0188     v(i) = Y(i,t) - a(mf(i)) - trend(i,t);
0189     Fi   = Pstar(mf(i),mf(i));
0190     if Fi > crit
0191       Ki        = Pstar(:,mf(i));
0192       a        = a + Ki*v(i)/Fi;
0193       Pstar     = Pstar - Ki*transpose(Ki)/Fi;
0194       lik(t)     = lik(t) + log(Fi) + v(i)*v(i)/Fi;
0195     else
0196       %disp(['zero F term for observed ',int2str(i),' ',num2str(Fi)])
0197     end
0198   end    
0199   a             = T*a;
0200   Pstar         = T*Pstar*transpose(T) + QQ;
0201   notsteady     = ~(max(max(abs(Pstar-oldP)))<crit);
0202 end
0203 while t < smpl
0204   t = t+1;
0205   Pstar = oldP;
0206   for i=1:pp
0207     v(i) = Y(i,t) - a(mf(i)) - trend(i,t);
0208     Fi   = Pstar(mf(i),mf(i));
0209     if Fi > crit
0210       Ki         = Pstar(:,mf(i));
0211       a         = a + Ki*v(i)/Fi;
0212       Pstar     = Pstar - Ki*transpose(Ki)/Fi;
0213       lik(t)        = lik(t) + log(Fi) + v(i)*v(i)/Fi;
0214     else
0215       %disp(['zero F term for observed ',int2str(i),' ',num2str(Fi)])
0216     end
0217   end    
0218   a = T*a;
0219 end
0220 
0221 LIK = .5*(sum(lik(start:end))-(start-1)*lik(smpl+1)/smpl);
0222

Generated on Fri 16-Jun-2006 09:09:06 by m2html © 2003