File: conditional_particle_filter.m

package info (click to toggle)
dynare 6.4-2
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid
  • size: 67,648 kB
  • sloc: cpp: 79,109; ansic: 28,917; objc: 12,430; yacc: 4,528; pascal: 1,993; lex: 1,441; sh: 1,129; python: 634; makefile: 626; lisp: 163; xml: 18
file content (110 lines) | stat: -rw-r--r-- 4,319 bytes parent folder | download | duplicates (2)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
function [LIK,lik] = conditional_particle_filter(ReducedForm, Y, s, ParticleOptions, ThreadsOptions, options_, M_)

% Evaluates the likelihood of a non-linear model with a particle filter
%
% INPUTS
% - ReducedForm        [structure]    Matlab's structure describing the reduced form model.
% - Y                  [double]       p×T matrix of (detrended) data, where p is the number of observed variables.
% - s                  [integer]      scalar, likelihood evaluation starts at s (has to be smaller than T, the sample length provided in Y).
% - ParticlesOptions   [struct]
% - ThreadsOptions     [struct]
% - options_           [struct]
% - M_                 [struct]
%
% OUTPUTS
% - LIK                [double]    scalar, likelihood
% - lik                [double]    (T-s+1)×1 vector, density of observations in each period.
%
% REMARKS
% - The proposal is built using the Kalman updating step for each particle.
% - we need draws in the errors distributions
% Whether we use Monte-Carlo draws from a multivariate gaussian distribution
% as in Amisano & Tristani (JEDC 2010).
% Whether we use multidimensional Gaussian sparse grids approximations:
% - a univariate Kronrod-Paterson Gaussian quadrature combined by the Smolyak
% operator (ref: Winschel & Kratzig, 2010).
% - a spherical-radial cubature (ref: Arasaratnam & Haykin, 2009a,2009b).
% - a scaled unscented transform cubature (ref: Julier & Uhlmann 1997, van der
% Merwe & Wan 2003).
%
% Pros:
% - Allows using current observable information in the proposal
% - The use of sparse grids Gaussian approximation is much faster than the Monte-Carlo approach
% Cons:
% - The use of the Kalman updating step may biais the proposal distribution since
% it has been derived in a linear context and is implemented in a nonlinear
% context. That is why particle resampling is performed.

% Copyright © 2009-2020 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/>.

% Set default for third input argument.
if isempty(s)
    s = 1;
end

T = size(Y,2);
p = length(ReducedForm.mf1);
n = ParticleOptions.number_of_particles;

% Get covariance matrices
Q = ReducedForm.Q;
H = ReducedForm.H;
if isempty(H)
    H = 0;
    H_lower_triangular_cholesky = 0;
else
    H_lower_triangular_cholesky = chol(H)';
end

% Get initial condition for the state vector.
StateVectorMean = ReducedForm.StateVectorMean;
StateVectorVarianceSquareRoot = chol(ReducedForm.StateVectorVariance)';
state_variance_rank = size(StateVectorVarianceSquareRoot, 2);
Q_lower_triangular_cholesky = chol(Q)';

% Set seed for randn().
options_=set_dynare_seed_local_options(options_,'default');

% Initialization of the likelihood.
lik = NaN(T, 1);
ks = 0;
StateParticles = bsxfun(@plus, StateVectorVarianceSquareRoot*randn(state_variance_rank, n), StateVectorMean);
SampleWeights = ones(1, n)/n;

for t=1:T
    flags = false(n, 1);
    for i=1:n
        [StateParticles(:,i), SampleWeights(i), flags(i)] = ...
            conditional_filter_proposal(ReducedForm, Y(:,t), StateParticles(:,i), SampleWeights(i), Q_lower_triangular_cholesky, H_lower_triangular_cholesky, H, ParticleOptions, ThreadsOptions, options_, M_);
    end
    if any(flags)
        LIK = -Inf;
        lik(t) = -Inf;
        return 
    end
    SumSampleWeights = sum(SampleWeights);
    lik(t) = log(SumSampleWeights);
    SampleWeights = SampleWeights./SumSampleWeights;
    if (ParticleOptions.resampling.status.generic && neff(SampleWeights)<ParticleOptions.resampling.threshold*T) || ParticleOptions.resampling.status.systematic
        ks = ks + 1;
        StateParticles = resample(StateParticles', SampleWeights', ParticleOptions)';
        SampleWeights = ones(1, n)/n;
    end
end

LIK = -sum(lik(s:end));