missing_DiffuseKalmanSmootherH1_Z.m 15.8 KB
Newer Older
1
function [alphahat,epsilonhat,etahat,atilde,P,aK,PK,decomp,V] = missing_DiffuseKalmanSmootherH1_Z(a_initial,T,Z,R,Q,H,Pinf1,Pstar1,Y,pp,mm,smpl,data_index,nk,kalman_tol,diffuse_kalman_tol,decomp_flag,state_uncertainty_flag)
2

3
% function [alphahat,epsilonhat,etahat,a,aK,PK,decomp] = DiffuseKalmanSmoother1(T,Z,R,Q,H,Pinf1,Pstar1,Y,pp,mm,smpl,data_index,nk,kalman_tol,diffuse_kalman_tol,decomp_flag,state_uncertainty_flag)
4
% Computes the diffuse kalman smoother without measurement error, in the case of a non-singular var-cov matrix.
5
6
%
% INPUTS
7
%    a_initial:mm*1 vector of initial (predicted) states
8
%    T:        mm*mm matrix
9
%    Z:        pp*mm matrix
10
11
%    R:        mm*rr matrix
%    Q:        rr*rr matrix
12
%    H:        pp*pp matrix variance of measurement errors
13
14
15
16
17
18
19
%    Pinf1:    mm*mm diagonal matrix with with q ones and m-q zeros
%    Pstar1:   mm*mm variance-covariance matrix with stationary variables
%    Y:        pp*1 vector
%    pp:       number of observed variables
%    mm:       number of state variables
%    smpl:     sample size
%    data_index                   [cell]      1*smpl cell of column vectors of indices.
20
21
%    nk        number of forecasting periods
%    kalman_tol   tolerance for reciprocal condition number
22
%    diffuse_kalman_tol   tolerance for reciprocal condition number (for Finf) and the rank of Pinf
23
%    decomp_flag  if true, compute filter decomposition
24
25
%    state_uncertainty_flag     if true, compute uncertainty about smoothed
%                               state estimate
26
%
27
28
29
30
31
32
33
34
35
36
37
38
% OUTPUTS
%    alphahat: smoothed variables (a_{t|T})
%    epsilonhat:smoothed measurement errors
%    etahat:   smoothed shocks
%    atilde:   matrix of updated variables (a_{t|t})
%    aK:       3D array of k step ahead filtered state variables (a_{t+k|t)
%              (meaningless for periods 1:d)
%    P:        3D array of one-step ahead forecast error variance
%              matrices
%    PK:       4D array of k-step ahead forecast error variance
%              matrices (meaningless for periods 1:d)
%    decomp:   decomposition of the effect of shocks on filtered values
39
%    V:        3D array of state uncertainty matrices
40
%
41
42
43
44
% Notes:
%   Outputs are stored in decision-rule order, i.e. to get variables in order of declaration
%   as in M_.endo_names, ones needs code along the lines of:
%   variables_declaration_order(dr.order_var,:) = alphahat
45
%
46
47
% SPECIAL REQUIREMENTS
%   See "Filtering and Smoothing of State Vector for Diffuse State Space
48
49
%   Models", S.J. Koopman and J. Durbin (2003, in Journal of Time Series
%   Analysis, vol. 24(1), pp. 85-98).
50
51
%   Durbin/Koopman (2012): "Time Series Analysis by State Space Methods", Oxford University Press,
%   Second Edition, Ch. 5
52

53
% Copyright (C) 2004-2018 Dynare Team
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
%
% 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 <http://www.gnu.org/licenses/>.

% modified by M. Ratto:
% new output argument aK (1-step to k-step predictions)
% new options_.nk: the max step ahed prediction in aK (default is 4)
% new crit1 value for rank of Pinf
74
% it is assured that P is symmetric
75
76
77
78
79
80
81

d = 0;
decomp = [];
spinf           = size(Pinf1);
spstar          = size(Pstar1);
v               = zeros(pp,smpl);
a               = zeros(mm,smpl+1);
82
a(:,1)          = a_initial;
83
84
85
86
87
atilde          = zeros(mm,smpl);
aK              = zeros(nk,mm,smpl+nk);
PK              = zeros(nk,mm,mm,smpl+nk);
iF              = zeros(pp,pp,smpl);
Fstar           = zeros(pp,pp,smpl);
88
iFstar          = zeros(pp,pp,smpl);
89
90
91
92
iFinf           = zeros(pp,pp,smpl);
K               = zeros(mm,pp,smpl);
L               = zeros(mm,mm,smpl);
Linf            = zeros(mm,mm,smpl);
93
Lstar           = zeros(mm,mm,smpl);
94
Kstar           = zeros(mm,pp,smpl);
95
Kinf            = zeros(mm,pp,smpl);
96
P               = zeros(mm,mm,smpl+1);
97
Pstar           = zeros(spstar(1),spstar(2),smpl+1);
98
Pstar(:,:,1)    = Pstar1;
99
Pinf            = zeros(spinf(1),spinf(2),smpl+1);
100
Pinf(:,:,1)     = Pinf1;
101
102
103
104
105
rr              = size(Q,1);
QQ              = R*Q*transpose(R);
QRt             = Q*transpose(R);
alphahat        = zeros(mm,smpl);
etahat          = zeros(rr,smpl);
106
epsilonhat      = zeros(rr,smpl);
107
r               = zeros(mm,smpl+1);
108
Finf_singular   = zeros(1,smpl);
109
110
111
112
113
114
if state_uncertainty_flag
    V               = zeros(mm,mm,smpl);
    N               = zeros(mm,mm,smpl+1);
else
    V=[];
end
115
116

t = 0;
117
while rank(Pinf(:,:,t+1),diffuse_kalman_tol) && t<smpl
118
119
120
    t = t+1;
    di = data_index{t};
    if isempty(di)
121
122
123
124
        %no observations, propagate forward without updating based on
        %observations
        atilde(:,t)     = a(:,t);
        a(:,t+1)        = T*atilde(:,t);
125
126
127
128
        Linf(:,:,t)     = T;
        Pstar(:,:,t+1)  = T*Pstar(:,:,t)*T' + QQ;
        Pinf(:,:,t+1)   = T*Pinf(:,:,t)*T';
    else
129
130
131
        ZZ = Z(di,:);                                                       %span selector matrix
        v(di,t)= Y(di,t) - ZZ*a(:,t);                                       %get prediction error v^(0) in (5.13) DK (2012)
        Finf = ZZ*Pinf(:,:,t)*ZZ';                                          % (5.7) in DK (2012)
132
        if rcond(Finf) < diffuse_kalman_tol                                 %F_{\infty,t} = 0
133
            if ~all(abs(Finf(:)) < diffuse_kalman_tol)                      %rank-deficient but not rank 0
134
                                                                            % The univariate diffuse kalman filter should be used.
135
                alphahat = Inf;
136
                return
137
            else                                                            %rank of F_{\infty,t} is 0
138
                Finf_singular(1,t) = 1;
139
140
141
                Fstar(di,di,t)  = ZZ*Pstar(:,:,t)*ZZ' + H(di,di);             % (5.7) in DK (2012)
                if rcond(Fstar(di,di,t)) < kalman_tol                         %F_{*} is singular
                    if ~all(all(abs(Fstar(di,di,t))<kalman_tol))
142
                        % The univariate diffuse kalman filter should be used.
143
                        alphahat = Inf;
144
                        return
145
                    else %rank 0
146
                        a(:,t+1) = T*a(:,t);
147
148
149
150
                        Pstar(:,:,t+1) = T*Pstar(:,:,t)*transpose(T)+QQ;
                        Pinf(:,:,t+1)  = T*Pinf(:,:,t)*transpose(T);
                    end
                else
151
152
                    iFstar(di,di,t) = inv(Fstar(di,di,t));
                    Kstar(:,di,t)   = Pstar(:,:,t)*ZZ'*iFstar(di,di,t);       %(5.15) of DK (2012) with Kstar=T^{-1}*K^(0)
153
                    Pinf(:,:,t+1)   = T*Pinf(:,:,t)*transpose(T);           % DK (2012), 5.16
154
155
                    Lstar(:,:,t)    = T - T*Kstar(:,di,t)*ZZ;               %L^(0) in DK (2012), eq. 5.12
                    Pstar(:,:,t+1)  = T*Pstar(:,:,t)*Lstar(:,:,t)'+QQ;      % (5.17) DK (2012)
156
                    a(:,t+1)        = T*(a(:,t)+Kstar(:,di,t)*v(di,t));       % (5.13) DK (2012)
157
158
159
                end
            end
        else
160
            %see notes in kalman_filter_d.m for details of computations
161
            iFinf(di,di,t)  = inv(Finf);
162
163
164
165
            Kinf(:,di,t)    = Pinf(:,:,t)*ZZ'*iFinf(di,di,t);               %define Kinf=T^{-1}*K_0 with M_{\infty}=Pinf*Z'
            atilde(:,t)     = a(:,t) + Kinf(:,di,t)*v(di,t);
            Linf(:,:,t)     = T - T*Kinf(:,di,t)*ZZ;                        %L^(0) in DK (2012), eq. 5.12
            Fstar(di,di,t)  = ZZ*Pstar(:,:,t)*ZZ' + H(di,di);               %(5.7) DK(2012)
166
            Kstar(:,di,t)   = (Pstar(:,:,t)*ZZ'-Kinf(:,di,t)*Fstar(di,di,t))*iFinf(di,di,t); %(5.12) DK(2012) with Kstar=T^{-1}*K^(1); note that there is a typo in DK (2003) with "+ Kinf" instead of "- Kinf", but it is correct in their appendix
167
168
            Pstar(:,:,t+1)  = T*Pstar(:,:,t)*Linf(:,:,t)'-T*Kinf(:,di,t)*Finf*Kstar(:,di,t)'*T' + QQ; %(5.14) DK(2012)
            Pinf(:,:,t+1)   = T*Pinf(:,:,t)*Linf(:,:,t)';                   %(5.14) DK(2012)
169
        end
170
        a(:,t+1)            = T*atilde(:,t);
171
172
        aK(1,:,t+1)         = a(:,t+1);
        % isn't a meaningless as long as we are in the diffuse part? MJ
173
        for jnk=2:nk
174
175
176
177
178
179
180
            aK(jnk,:,t+jnk) = T*dynare_squeeze(aK(jnk-1,:,t+jnk-1));
        end
    end
end
d = t;
P(:,:,d+1) = Pstar(:,:,d+1);
iFinf = iFinf(:,:,1:d);
181
iFstar= iFstar(:,:,1:d);
182
Linf  = Linf(:,:,1:d);
183
Lstar = Lstar(:,:,1:d);
184
185
186
187
Kstar = Kstar(:,:,1:d);
Pstar = Pstar(:,:,1:d);
Pinf  = Pinf(:,:,1:d);
notsteady = 1;
188
while notsteady && t<smpl
189
    t = t+1;
190
    P(:,:,t)=tril(P(:,:,t))+transpose(tril(P(:,:,t),-1));                   % make sure P is symmetric
191
192
    di = data_index{t};
    if isempty(di)
193
        atilde(:,t)     = a(:,t);
194
        L(:,:,t)        = T;
195
        P(:,:,t+1)      = T*P(:,:,t)*T' + QQ;                               %p. 111, DK(2012)
196
197
198
199
    else
        ZZ = Z(di,:);
        v(di,t)      = Y(di,t) - ZZ*a(:,t);
        F = ZZ*P(:,:,t)*ZZ' + H(di,di);
200
201
202
        sig=sqrt(diag(F));

        if any(diag(F)<kalman_tol) || rcond(F./(sig*sig')) < kalman_tol
203
            alphahat = Inf;
204
            return
205
        end
206
        iF(di,di,t)   = inv(F./(sig*sig'))./(sig*sig');
207
208
209
210
        PZI         = P(:,:,t)*ZZ'*iF(di,di,t);
        atilde(:,t) = a(:,t) + PZI*v(di,t);
        K(:,di,t)    = T*PZI;
        L(:,:,t)    = T-K(:,di,t)*ZZ;
211
        P(:,:,t+1)  = T*P(:,:,t)*L(:,:,t)' + QQ;
212
213
    end
    a(:,t+1)    = T*atilde(:,t);
214
    Pf          = P(:,:,t+1);
215
216
    aK(1,:,t+1) = a(:,t+1);
    for jnk=1:nk
217
218
219
        if jnk>1
            Pf = T*Pf*T' + QQ;
        end
220
221
222
223
224
        PK(jnk,:,:,t+jnk) = Pf;
        if jnk>1
            aK(jnk,:,t+jnk) = T*dynare_squeeze(aK(jnk-1,:,t+jnk-1));
        end
    end
225
    %    notsteady   = ~(max(max(abs(P(:,:,t+1)-P(:,:,t))))<kalman_tol);
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
end
% $$$ if t<smpl
% $$$     PZI_s = PZI;
% $$$     K_s = K(:,:,t);
% $$$     iF_s = iF(:,:,t);
% $$$     P_s = P(:,:,t+1);
% $$$     P  = cat(3,P(:,:,1:t),repmat(P_s,[1 1 smpl-t]));
% $$$     iF = cat(3,iF(:,:,1:t),repmat(iF_s,[1 1 smpl-t]));
% $$$     L  = cat(3,L(:,:,1:t),repmat(T-K_s*Z,[1 1 smpl-t]));
% $$$     K  = cat(3,K(:,:,1:t),repmat(T*P_s*Z'*iF_s,[1 1 smpl-t]));
% $$$ end
% $$$ while t<smpl
% $$$     t=t+1;
% $$$     v(:,t) = Y(:,t) - Z*a(:,t);
% $$$     atilde(:,t) = a(:,t) + PZI*v(:,t);
% $$$     a(:,t+1) = T*atilde(:,t);
% $$$     Pf          = P(:,:,t);
% $$$     for jnk=1:nk,
% $$$   Pf = T*Pf*T' + QQ;
% $$$         aK(jnk,:,t+jnk) = T^jnk*atilde(:,t);
% $$$   PK(jnk,:,:,t+jnk) = Pf;
% $$$     end
% $$$ end
249
%% backward pass; r_T and N_T, stored in entry (smpl+1) were initialized at 0
250
251
252
253
254
t = smpl+1;
while t>d+1
    t = t-1;
    di = data_index{t};
    if isempty(di)
255
        % in this case, L is simply T due to Z=0, so that DK (2012), eq. 4.93 obtains
256
        r(:,t) = L(:,:,t)'*r(:,t+1);                                        %compute r_{t-1}, DK (2012), eq. 4.38 with Z=0
257
258
259
        if state_uncertainty_flag
            N(:,:,t)=L(:,:,t)'*N(:,:,t+1)*L(:,:,t);                         %compute N_{t-1}, DK (2012), eq. 4.42 with Z=0
        end
260
261
    else
        ZZ = Z(di,:);
262
        r(:,t) = ZZ'*iF(di,di,t)*v(di,t) + L(:,:,t)'*r(:,t+1);              %compute r_{t-1}, DK (2012), eq. 4.38
263
264
265
        if state_uncertainty_flag
            N(:,:,t)=ZZ'*iF(di,di,t)*ZZ+L(:,:,t)'*N(:,:,t+1)*L(:,:,t);      %compute N_{t-1}, DK (2012), eq. 4.42
        end
266
    end
267
268
    alphahat(:,t)       = a(:,t) + P(:,:,t)*r(:,t);                         %DK (2012), eq. 4.35
    etahat(:,t) = QRt*r(:,t);                                               %DK (2012), eq. 4.63
269
270
271
    if state_uncertainty_flag
        V(:,:,t)    = P(:,:,t)-P(:,:,t)*N(:,:,t)*P(:,:,t);                      %DK (2012), eq. 4.43
    end
272
end
273
if d %diffuse periods
274
275
     % initialize r_d^(0) and r_d^(1) as below DK (2012), eq. 5.23
    r0 = zeros(mm,d+1);
276
277
    r0(:,d+1) = r(:,d+1);   %set r0_{d}, i.e. shifted by one period
    r1 = zeros(mm,d+1);     %set r1_{d}, i.e. shifted by one period
278
279
280
281
282
    if state_uncertainty_flag
        %N_0 at (d+1) is N(d+1), so we can use N for continuing and storing N_0-recursion
        N_1=zeros(mm,mm,d+1);   %set N_1_{d}=0, i.e. shifted by one period, below  DK (2012), eq. 5.26
        N_2=zeros(mm,mm,d+1);   %set N_2_{d}=0, i.e. shifted by one period, below  DK (2012), eq. 5.26
    end
283
284
285
286
287
    for t = d:-1:1
        di = data_index{t};
        if isempty(di)
            r1(:,t) = Linf(:,:,t)'*r1(:,t+1);
        else
288
289
290
            if ~Finf_singular(1,t)
                r0(:,t) = Linf(:,:,t)'*r0(:,t+1);                                   % DK (2012), eq. 5.21 where L^(0) is named Linf
                r1(:,t) = Z(di,:)'*(iFinf(di,di,t)*v(di,t)-Kstar(:,di,t)'*T'*r0(:,t+1)) ...
291
                          + Linf(:,:,t)'*r1(:,t+1);                                       % DK (2012), eq. 5.21, noting that i) F^(1)=(F^Inf)^(-1)(see 5.10), ii) where L^(0) is named Linf, and iii) Kstar=T^{-1}*K^(1)
292
293
294
295
296
297
298
299
300
301
302
                if state_uncertainty_flag
                    L_1=(-T*Kstar(:,di,t)*Z(di,:));                                     % noting that Kstar=T^{-1}*K^(1)
                    N(:,:,t)=Linf(:,:,t)'*N(:,:,t+1)*Linf(:,:,t);                       % DK (2012), eq. 5.19, noting that L^(0) is named Linf
                    N_1(:,:,t)=Z(di,:)'*iFinf(di,di,t)*Z(di,:)+Linf(:,:,t)'*N_1(:,:,t+1)*Linf(:,:,t)...
                        +L_1'*N(:,:,t+1)*Linf(:,:,t);                                   % DK (2012), eq. 5.29; note that, compared to DK (2003) this drops the term (L_1'*N(:,:,t+1)*Linf(:,:,t))' in the recursion due to it entering premultiplied by Pinf when computing V, and Pinf*Linf'*N=0
                    N_2(:,:,t)=Z(di,:)'*(-iFinf(di,di,t)*Fstar(di,di,t)*iFinf(di,di,t))*Z(di,:) ...
                        + Linf(:,:,t)'*N_2(:,:,t+1)*Linf(:,:,t)...
                        + Linf(:,:,t)'*N_1(:,:,t+1)*L_1...
                        + L_1'*N_1(:,:,t+1)'*Linf(:,:,t)...
                        + L_1'*N(:,:,t+1)*L_1;                            % DK (2012), eq. 5.29
                end
303
            else
304
                r0(:,t) = Z(di,:)'*iFstar(di,di,t)*v(di,t)-Lstar(:,:,t)'*r0(:,t+1); % DK (2003), eq. (14)
305
                r1(:,t) = T'*r1(:,t+1);                                             % DK (2003), eq. (14)
306
307
                if state_uncertainty_flag
                    N(:,:,t)=Z(di,:)'*iFstar(di,di,t)*Z(di,:)...
308
                             +Lstar(:,:,t)'*N(:,:,t+1)*Lstar(:,:,t);                     % DK (2003), eq. (14)
309
310
311
                    N_1(:,:,t)=T'*N_1(:,:,t+1)*Lstar(:,:,t);                            % DK (2003), eq. (14)
                    N_2(:,:,t)=T'*N_2(:,:,t+1)*T';                                      % DK (2003), eq. (14)
                end
312
            end
313
        end
314
315
        alphahat(:,t)   = a(:,t) + Pstar(:,:,t)*r0(:,t) + Pinf(:,:,t)*r1(:,t);      % DK (2012), eq. 5.23
        etahat(:,t)     = QRt*r0(:,t);                                              % DK (2012), p. 135
316
317
        if state_uncertainty_flag
            V(:,:,t)=Pstar(:,:,t)-Pstar(:,:,t)*N(:,:,t)*Pstar(:,:,t)...
318
319
320
                     -(Pinf(:,:,t)*N_1(:,:,t)*Pstar(:,:,t))'...
                     - Pinf(:,:,t)*N_1(:,:,t)*Pstar(:,:,t)...
                     - Pinf(:,:,t)*N_2(:,:,t)*Pinf(:,:,t);                                   % DK (2012), eq. 5.30
321
        end
322
323
324
    end
end

325
if decomp_flag
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
    decomp = zeros(nk,mm,rr,smpl+nk);
    ZRQinv = inv(Z*QQ*Z');
    for t = max(d,1):smpl
        di = data_index{t};
        % calculate eta_tm1t
        eta_tm1t = QRt*Z(di,:)'*iF(di,di,t)*v(di,t);
        AAA = P(:,:,t)*Z(di,:)'*ZRQinv(di,di)*bsxfun(@times,Z(di,:)*R,eta_tm1t');
        % calculate decomposition
        decomp(1,:,:,t+1) = AAA;
        for h = 2:nk
            AAA = T*AAA;
            decomp(h,:,:,t+h) = AAA;
        end
    end
end

epsilonhat = Y-Z*alphahat;