sequential_importance_particle_filter.m 6.96 KB
Newer Older
1
2
3
4
5
6
7
function [LIK,lik] = sequential_importance_particle_filter(ReducedForm,Y,start,DynareOptions)
% Evaluates the likelihood of a nonlinear model with a particle filter (optionally with resampling).

%@info:
%! @deftypefn {Function File} {@var{y}, @var{y_} =} sequential_importance_particle_filter (@var{ReducedForm},@var{Y}, @var{start}, @var{DynareOptions})
%! @anchor{particle/sequential_importance_particle_filter}
%! @sp 1
8
%! Evaluates the likelihood of a nonlinear model with a particle filter (optionally with resampling).
9
10
11
12
13
14
15
16
%!
%! @sp 2
%! @strong{Inputs}
%! @sp 1
%! @table @ @var
%! @item ReducedForm
%! Structure describing the state space model (built in @ref{non_linear_dsge_likelihood}).
%! @item Y
17
%! p*smpl matrix of doubles (p is the number of observed variables), the (detrended) data.
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
%! @item start
%! Integer scalar, likelihood evaluation starts at observation 'start'.
%! @item DynareOptions
%! Structure specifying Dynare's options.
%! @end table
%! @sp 2
%! @strong{Outputs}
%! @sp 1
%! @table @ @var
%! @item LIK
%! double scalar, value of (minus) the logged likelihood.
%! @item lik
%! smpl*1 vector of doubles, density of the observations at each period.
%! @end table
%! @sp 2
%! @strong{This function is called by:}
%! @ref{non_linear_dsge_likelihood}
%! @sp 2
%! @strong{This function calls:}
%!
%! @end deftypefn
%@eod:

41
% Copyright (C) 2011, 2012 Dynare Team
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
%
% 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/>.

58
% AUTHOR(S) frederic DOT karame AT univ DASH lemans DOT fr
59
%           stephane DOT adjemian AT univ DASH lemans DOT fr
60

61
62
persistent init_flag
persistent mf0 mf1
63
64
65
66
67
68
69
70
persistent number_of_particles
persistent sample_size number_of_observed_variables number_of_structural_innovations

% Set default value for start
if isempty(start)
    start = 1;
end

71
72
73
% Set flag for prunning
pruning = DynareOptions.particle.pruning;

74
% Get steady state and mean.
Frédéric Karamé's avatar
Frédéric Karamé committed
75
steadystate = ReducedForm.steadystate;
76
77
78
constant = ReducedForm.constant;
state_variables_steady_state = ReducedForm.state_variables_steady_state;

79
% Set persistent variables (if needed).
80
81
82
83
84
if isempty(init_flag)
    mf0 = ReducedForm.mf0;
    mf1 = ReducedForm.mf1;
    sample_size = size(Y,2);
    number_of_observed_variables = length(mf1);
85
    number_of_structural_innovations = length(ReducedForm.Q);
86
87
88
89
90
91
92
93
94
95
96
97
98
    number_of_particles = DynareOptions.particle.number_of_particles;
    init_flag = 1;
end

% Set local state space model (first order approximation).
ghx  = ReducedForm.ghx;
ghu  = ReducedForm.ghu;

% Set local state space model (second order approximation).
ghxx = ReducedForm.ghxx;
ghuu = ReducedForm.ghuu;
ghxu = ReducedForm.ghxu;

99
% Get covariance matrices.
100
101
102
103
104
105
Q = ReducedForm.Q;
H = ReducedForm.H;
if isempty(H)
    H = 0;
end

106
107
108
109
% Initialization of the likelihood.
const_lik = log(2*pi)*number_of_observed_variables;
lik  = NaN(sample_size,1);

110
111
112
% Get initial condition for the state vector.
StateVectorMean = ReducedForm.StateVectorMean;
StateVectorVarianceSquareRoot = reduced_rank_cholesky(ReducedForm.StateVectorVariance)';
113
114
115
116
117
118
if pruning
    StateVectorMean_ = StateVectorMean;
    StateVectorVarianceSquareRoot_ = StateVectorVarianceSquareRoot;
end

% Get the rank of StateVectorVarianceSquareRoot
119
state_variance_rank = size(StateVectorVarianceSquareRoot,2);
120
121

% Factorize the covariance matrix of the structural innovations
122
Q_lower_triangular_cholesky = chol(Q)';
123
124
125
[PredictedStateMean,PredictedStateVarianceSquareRoot,StateVectorMean,StateVectorVarianceSquareRoot] = ...
        gaussian_filter_bank(ReducedForm,Y(:,1),StateVectorMean,StateVectorVarianceSquareRoot,Q_lower_triangular_cholesky,Q_lower_triangular_cholesky,H,DynareOptions) ;
StateVectors = bsxfun(@plus,StateVectorVarianceSquareRoot*randn(state_variance_rank,number_of_particles),StateVectorMean) ;
126
127

% Set seed for randn().
128
set_dynare_seed('default');
129
130

% Initialization of the weights across particles.
131
weights = ones(1,number_of_particles)/number_of_particles ;
132
StateVectors = bsxfun(@plus,StateVectorVarianceSquareRoot*randn(state_variance_rank,number_of_particles),StateVectorMean);
133
134
135
136
137
if pruning
    StateVectors_ = StateVectors;
end

% Loop over observations
138
139
140
for t=1:sample_size
    yhat = bsxfun(@minus,StateVectors,state_variables_steady_state);
    epsilon = Q_lower_triangular_cholesky*randn(number_of_structural_innovations,number_of_particles);
141
142
143
144
145
146
    if pruning
        yhat_ = bsxfun(@minus,StateVectors_,state_variables_steady_state);
        [tmp, tmp_] = local_state_space_iteration_2(yhat,epsilon,ghx,ghu,constant,ghxx,ghuu,ghxu,yhat_,steadystate,DynareOptions.threads.local_state_space_iteration_2);
    else
        tmp = local_state_space_iteration_2(yhat,epsilon,ghx,ghu,constant,ghxx,ghuu,ghxu,DynareOptions.threads.local_state_space_iteration_2);
    end
147
    PredictedObservedMean = tmp(mf1,:)*transpose(weights);
148
149
    PredictionError = bsxfun(@minus,Y(:,t),tmp(mf1,:));
    dPredictedObservedMean = bsxfun(@minus,tmp(mf1,:),PredictedObservedMean);
150
    PredictedObservedVariance = bsxfun(@times,dPredictedObservedMean,weights)*dPredictedObservedMean' + H;
151
152
153
    lnw = -.5*(const_lik+log(det(PredictedObservedVariance))+sum(PredictionError.*(PredictedObservedVariance\PredictionError),1));
    dfac = max(lnw);
    wtilde = weights.*exp(lnw-dfac);
154
    lik(t) = log(sum(wtilde))+dfac;
155
    weights = wtilde/sum(wtilde);
156
157
158
159
160
161
    if (strcmp(DynareOptions.particle.resampling.status,'generic') && neff(weights)<DynareOptions.particle.resampling.neff_threshold*sample_size ) || strcmp(DynareOptions.particle.resampling.status,'systematic')
        idx = resample(weights,DynareOptions.particle.resampling.method1,DynareOptions.particle.resampling.method2);
        StateVectors = tmp(mf0,idx);
        if pruning
            StateVectors_ = tmp_(mf0,idx);
        end
Stéphane Adjemian's avatar
Stéphane Adjemian committed
162
        weights = ones(1,number_of_particles)/number_of_particles;
Stéphane Adjemian's avatar
Stéphane Adjemian committed
163
    elseif strcmp(DynareOptions.particle.resampling.status,'smoothed')
Stéphane Adjemian's avatar
Stéphane Adjemian committed
164
        StateVectors = multivariate_smooth_resampling(weights',tmp(mf0,:)',number_of_particles,DynareOptions.particle.resampling.number_of_partitions)';
165
166
167
        if pruning
            StateVectors_ = multivariate_smooth_resampling(weights',tmp_(mf0,:)',number_of_particles,DynareOptions.particle.resampling.number_of_partitions)';
        end
Stéphane Adjemian's avatar
Stéphane Adjemian committed
168
        weights = ones(1,number_of_particles)/number_of_particles;
Stéphane Adjemian's avatar
Stéphane Adjemian committed
169
    elseif strcmp(DynareOptions.particle.resampling.status,'none')
170
        StateVectors = tmp(mf0,:);
171
        if pruning
Frédéric Karamé's avatar
Frédéric Karamé committed
172
            StateVectors_ = tmp_(mf0,:);
173
        end
174
175
176
    end
end

177
178
179
180
181
182
LIK = -sum(lik(start:end));



function n = neff(w)
    n = dot(w,w);