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 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135
|
function [LIK,lik] = gaussian_filter(ReducedForm, Y, start, ParticleOptions, ThreadsOptions, options_, M_)
% [LIK,lik] = gaussian_filter(ReducedForm, Y, start, ParticleOptions, ThreadsOptions, options_, M_)
% Evaluates the likelihood of a non-linear model approximating the
% predictive (prior) and filtered (posterior) densities for state variables
% by gaussian distributions.
% Gaussian approximation is done by:
% - a spherical-radial cubature (ref: Arasaratnam & Haykin, 2009).
% - a scaled unscented transform cubature (ref: Julier & Uhlmann 1995)
% - Monte-Carlo draws from a multivariate gaussian distribution.
% First and second moments of prior and posterior state densities are computed
% from the resulting nodes/particles and allows to generate new distributions at the
% following observation.
% Pros: The use of nodes is much faster than Monte-Carlo Gaussian particle and standard particles
% filters since it treats a lesser number of particles. Furthermore, in all cases, there is no need
% of resampling.
% Cons: estimations may be biaised if the model is truly non-gaussian
% since predictive and filtered densities are unimodal.
%
% INPUTS
% Reduced_Form [structure] Matlab's structure describing the reduced form model.
% Y [double] matrix of original observed variables.
% start [double] structural parameters.
% ParticleOptions [structure] Matlab's structure describing options concerning particle filtering.
% ThreadsOptions [structure] Matlab's structure.
%
% OUTPUTS
% LIK [double] scalar, likelihood
% lik [double] vector, density of observations in each period.
%
% REFERENCES
%
% NOTES
% The vector "lik" is used to evaluate the jacobian of the likelihood.
% Copyright © 2009-2019 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
if isempty(start)
start = 1;
end
mf0 = ReducedForm.mf0;
mf1 = ReducedForm.mf1;
sample_size = size(Y,2);
number_of_state_variables = length(mf0);
number_of_observed_variables = length(mf1);
number_of_particles = ParticleOptions.number_of_particles;
% compute gaussian quadrature nodes and weights on states and shocks
if ParticleOptions.distribution_approximation.cubature
[nodes2, weights2] = spherical_radial_sigma_points(number_of_state_variables);
weights_c2 = weights2;
elseif ParticleOptions.distribution_approximation.unscented
[nodes2, weights2, weights_c2] = unscented_sigma_points(number_of_state_variables,ParticleOptions);
else
if ~ParticleOptions.distribution_approximation.montecarlo
error('This approximation for the proposal is unknown!')
end
end
if ParticleOptions.distribution_approximation.montecarlo
options_=set_dynare_seed_local_options(options_,'default');
end
% Get covariance matrices
Q = ReducedForm.Q;
H = ReducedForm.H;
if isempty(H)
H = 0;
H_lower_triangular_cholesky = 0;
else
H_lower_triangular_cholesky = reduced_rank_cholesky(H)';
end
% Get initial condition for the state vector.
StateVectorMean = ReducedForm.StateVectorMean;
StateVectorVarianceSquareRoot = reduced_rank_cholesky(ReducedForm.StateVectorVariance)';
state_variance_rank = size(StateVectorVarianceSquareRoot,2);
Q_lower_triangular_cholesky = reduced_rank_cholesky(Q)';
% Initialization of the likelihood.
const_lik = (2*pi)^(number_of_observed_variables/2) ;
lik = NaN(sample_size,1);
LIK = NaN;
for t=1:sample_size
[PredictedStateMean, PredictedStateVarianceSquareRoot, StateVectorMean, StateVectorVarianceSquareRoot] = ...
gaussian_filter_bank(ReducedForm, Y(:,t), StateVectorMean, StateVectorVarianceSquareRoot, Q_lower_triangular_cholesky, H_lower_triangular_cholesky, ...
H, ParticleOptions, ThreadsOptions, options_, M_);
if ParticleOptions.distribution_approximation.cubature || ParticleOptions.distribution_approximation.unscented
StateParticles = bsxfun(@plus, StateVectorMean, StateVectorVarianceSquareRoot*nodes2');
IncrementalWeights = gaussian_densities(Y(:,t), StateVectorMean, StateVectorVarianceSquareRoot, PredictedStateMean, ...
PredictedStateVarianceSquareRoot, StateParticles, H, const_lik, ...
weights2, weights_c2, ReducedForm, ThreadsOptions, ...
options_, M_);
SampleWeights = weights2.*IncrementalWeights;
else
StateParticles = bsxfun(@plus, StateVectorVarianceSquareRoot*randn(state_variance_rank, number_of_particles), StateVectorMean) ;
IncrementalWeights = gaussian_densities(Y(:,t), StateVectorMean, StateVectorVarianceSquareRoot, PredictedStateMean, ...
PredictedStateVarianceSquareRoot,StateParticles,H,const_lik, ...
1/number_of_particles,1/number_of_particles,ReducedForm,ThreadsOptions, ...
options_, M_);
SampleWeights = IncrementalWeights/number_of_particles;
end
SampleWeights = SampleWeights + 1e-6*ones(size(SampleWeights, 1), 1);
SumSampleWeights = sum(SampleWeights);
lik(t) = log(SumSampleWeights);
SampleWeights = SampleWeights./SumSampleWeights;
if not(ParticleOptions.distribution_approximation.cubature || ParticleOptions.distribution_approximation.unscented)
if (ParticleOptions.resampling.status.generic && neff(SampleWeights)<ParticleOptions.resampling.threshold*sample_size) || ParticleOptions.resampling.status.systematic
StateParticles = resample(StateParticles', SampleWeights, ParticleOptions)';
SampleWeights = ones(number_of_particles, 1)/number_of_particles;
end
end
StateVectorMean = StateParticles*SampleWeights;
temp = bsxfun(@minus, StateParticles, StateVectorMean);
StateVectorVarianceSquareRoot = reduced_rank_cholesky(bsxfun(@times,SampleWeights',temp)*temp')';
end
LIK = -sum(lik(start:end));
|