factorize some codes across options and modify the definition of mixtures.

rm-particles^2
Frédéric Karamé 2015-01-23 14:40:55 +01:00
parent 6e7090f1f5
commit ec4f2da2fd
2 changed files with 57 additions and 87 deletions

View File

@ -3,28 +3,24 @@ function [LIK,lik] = gaussian_filter(ReducedForm, Y, start, ParticleOptions, Thr
% predictive (prior) and filtered (posterior) densities for state variables
% by gaussian distributions.
% Gaussian approximation is done by:
% - a Kronrod-Paterson gaussian quadrature with a limited number of nodes.
% Multidimensional quadrature is obtained by the Smolyak operator (ref: Winschel & Kratzig, 2010).
% - a spherical-radial cubature (ref: Arasaratnam & Haykin, 2008,2009).
% - a scaled unscented transform cubature (ref: )
% - 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.
% => The use of nodes is much faster than Monte-Carlo Gaussian particle and standard particles
% filters since it treats a lesser number of particles and there is no need
% 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.
% However, estimations may reveal biaised if the model is truly non-gaussian
% Cons: estimations may be biaised if the model is truly non-gaussian
% since predictive and filtered densities are unimodal.
%
% INPUTS
% reduced_form_model [structure] Matlab's structure describing the reduced form model.
% reduced_form_model.measurement.H [double] (pp x pp) variance matrix of measurement errors.
% reduced_form_model.state.Q [double] (qq x qq) variance matrix of state errors.
% reduced_form_model.state.dr [structure] output of resol.m.
% Y [double] pp*smpl matrix of (detrended) data, where pp is the maximum number of observed variables.
% start [integer] scalar, likelihood evaluation starts at 'start'.
% smolyak_accuracy [integer] scalar.
% 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
@ -34,7 +30,7 @@ function [LIK,lik] = gaussian_filter(ReducedForm, Y, start, ParticleOptions, Thr
%
% NOTES
% The vector "lik" is used to evaluate the jacobian of the likelihood.
% Copyright (C) 2009-2013 Dynare Team
% Copyright (C) 2009-2015 Dynare Team
%
% This file is part of Dynare.
%
@ -110,12 +106,7 @@ const_lik = (2*pi)^(number_of_observed_variables/2) ;
lik = NaN(sample_size,1);
LIK = NaN;
SampleWeights = 1/number_of_particles ;
ks = 0 ;
%Estimate = zeros(number_of_state_variables,sample_size) ;
%V_Estimate = zeros(number_of_state_variables,number_of_state_variables,sample_size) ;
for t=1:sample_size
% build the proposal
[PredictedStateMean,PredictedStateVarianceSquareRoot,StateVectorMean,StateVectorVarianceSquareRoot] = ...
gaussian_filter_bank(ReducedForm,Y(:,t),StateVectorMean,StateVectorVarianceSquareRoot,Q_lower_triangular_cholesky,H_lower_triangular_cholesky,H,ParticleOptions,ThreadsOptions) ;
if ParticleOptions.distribution_approximation.cubature || ParticleOptions.distribution_approximation.unscented
@ -126,40 +117,21 @@ for t=1:sample_size
PredictedStateVarianceSquareRoot,StateParticles,H,const_lik,...
weights2,weights_c2,ReducedForm,ThreadsOptions) ;
SampleWeights = weights2.*IncrementalWeights ;
SumSampleWeights = sum(SampleWeights) ;
lik(t) = log(SumSampleWeights) ;
SampleWeights = SampleWeights./SumSampleWeights ;
StateVectorMean = StateParticles*SampleWeights ;
temp = bsxfun(@minus,StateParticles,StateVectorMean) ;
StateVectorVarianceSquareRoot = reduced_rank_cholesky( bsxfun(@times,SampleWeights',temp)*temp' )';
else % Monte-Carlo draws
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) ;
%SampleWeights = SampleWeights.*IncrementalWeights ;
SampleWeights = IncrementalWeights/number_of_particles ;
SumSampleWeights = sum(SampleWeights) ;
%VarSampleWeights = IncrementalWeights-SumSampleWeights ;
%VarSampleWeights = VarSampleWeights*VarSampleWeights'/(number_of_particles-1) ;
lik(t) = log(SumSampleWeights) ; %+ .5*VarSampleWeights/(number_of_particles*(SumSampleWeights*SumSampleWeights)) ;
SampleWeights = SampleWeights./SumSampleWeights ;
% Neff = neff(SampleWeights) ;
% if (Neff<0.5*sample_size && ParticleOptions.resampling.status.generic) || ParticleOptions.resampling.status.systematic
% ks = ks + 1 ;
% StateParticles = resample(StateParticles',SampleWeights,ParticleOptions)' ;
% StateVectorMean = mean(StateParticles,2) ;
% StateVectorVarianceSquareRoot = reduced_rank_cholesky( (StateParticles*StateParticles')/(number_of_particles-1) - StateVectorMean*(StateVectorMean') )';
% SampleWeights = 1/number_of_particles ;
% elseif ParticleOptions.resampling.status.none
StateVectorMean = StateParticles*SampleWeights ;
temp = bsxfun(@minus,StateParticles,StateVectorMean) ;
StateVectorVarianceSquareRoot = reduced_rank_cholesky( bsxfun(@times,SampleWeights',temp)*temp' )';
%disp(StateVectorVarianceSquareRoot)
% end
end
SumSampleWeights = sum(SampleWeights) ;
lik(t) = log(SumSampleWeights) ;
SampleWeights = SampleWeights./SumSampleWeights ;
StateVectorMean = StateParticles*SampleWeights ;
temp = bsxfun(@minus,StateParticles,StateVectorMean) ;
StateVectorVarianceSquareRoot = reduced_rank_cholesky( bsxfun(@times,SampleWeights',temp)*temp' )';
end
LIK = -sum(lik(start:end));

View File

@ -52,8 +52,7 @@ function [LIK,lik] = gaussian_mixture_filter(ReducedForm,Y,start,ParticleOptions
% You should have received a copy of the GNU General Public License
% along with Dynare. If not, see <http://www.gnu.org/licenses/>.
persistent init_flag mf0 mf1 Gprime Gsecond
persistent init_flag mf0 mf1
persistent nodes weights weights_c I J G number_of_particles
persistent sample_size number_of_state_variables number_of_observed_variables number_of_structural_innovations
@ -71,16 +70,10 @@ if isempty(init_flag)
number_of_observed_variables = length(mf1);
number_of_structural_innovations = length(ReducedForm.Q);
G = ParticleOptions.mixture_state_variables; % number of GM components in state
I = 1 ; %ParticleOptions.mixture_structural_shocks ; % number of GM components in structural noise
J = 1 ; %ParticleOptions.mixture_measurement_shocks ; % number of GM components in observation noise
Gprime = G*I ;
Gsecond = G*I*J ;
number_of_particles = ParticleOptions.number_of_particles;
init_flag = 1;
end
SampleWeights = ones(Gsecond,1)/Gsecond ;
% compute gaussian quadrature nodes and weights on states and shocks
if isempty(nodes)
if ParticleOptions.distribution_approximation.cubature
@ -111,29 +104,54 @@ else
end
Q_lower_triangular_cholesky = reduced_rank_cholesky(Q)';
% Initialize all matrices
% Initialize mixtures
StateWeights = ones(1,G)/G ;
StateMu = ReducedForm.StateVectorMean*ones(1,G) ;
StateMu = ReducedForm.StateVectorMean ;
StateSqrtP = zeros(number_of_state_variables,number_of_state_variables,G) ;
temp = reduced_rank_cholesky(ReducedForm.StateVectorVariance)' ;
StateMu = bsxfun(@plus,StateMu,bsxfun(@times,diag(temp),(-(G-1)/2:1:(G-1)/2))/10) ;
for g=1:G
StateSqrtP(:,:,g) = temp ;
end
StructuralShocksWeights = ones(1,I)/I ;
StructuralShocksMu = zeros(number_of_structural_innovations,I) ;
if ParticleOptions.proposal_approximation.cubature
[StructuralShocksMu,StructuralShocksWeights] = spherical_radial_sigma_points(number_of_structural_innovations);
StructuralShocksWeights = ones(size(StructuralShocksMu,1),1)*StructuralShocksWeights ;
elseif ParticleOptions.proposal_approximation.unscented
[StructuralShocksMu,StructuralShocksWeights,raf] = unscented_sigma_points(number_of_structural_innovations,ParticleOptions);
else
if ~ParticleOptions.distribution_approximation.montecarlo
error('Estimation: This approximation for the proposal is not implemented or unknown!')
end
end
I = size(StructuralShocksWeights,1) ;
StructuralShocksMu = Q_lower_triangular_cholesky*(StructuralShocksMu') ;
StructuralShocksSqrtP = zeros(number_of_structural_innovations,number_of_structural_innovations,I) ;
for i=1:I
StructuralShocksSqrtP(:,:,i) = Q_lower_triangular_cholesky ;
end
ObservationShocksWeights = ones(1,J)/J ;
ObservationShocksMu = zeros(number_of_observed_variables,J) ;
if ParticleOptions.proposal_approximation.cubature
[ObservationShocksMu,ObservationShocksWeights] = spherical_radial_sigma_points(number_of_observed_variables);
ObservationShocksWeights = ones(size(ObservationShocksMu,1),1)*ObservationShocksWeights;
elseif ParticleOptions.proposal_approximation.unscented
[ObservationShocksMu,ObservationShocksWeights,raf] = unscented_sigma_points(number_of_observed_variables,ParticleOptions);
else
if ~ParticleOptions.distribution_approximation.montecarlo
error('Estimation: This approximation for the proposal is not implemented or unknown!')
end
end
J = size(ObservationShocksWeights,1) ;
ObservationShocksMu = H_lower_triangular_cholesky*(ObservationShocksMu') ;
ObservationShocksSqrtP = zeros(number_of_observed_variables,number_of_observed_variables,J) ;
for j=1:J
ObservationShocksSqrtP(:,:,j) = H_lower_triangular_cholesky ;
end
Gprime = G*I ;
Gsecond = G*I*J ;
SampleWeights = ones(Gsecond,1)/Gsecond ;
StateWeightsPrior = zeros(1,Gprime) ;
StateMuPrior = zeros(number_of_state_variables,Gprime) ;
StateSqrtPPrior = zeros(number_of_state_variables,number_of_state_variables,Gprime) ;
@ -142,10 +160,8 @@ StateWeightsPost = zeros(1,Gsecond) ;
StateMuPost = zeros(number_of_state_variables,Gsecond) ;
StateSqrtPPost = zeros(number_of_state_variables,number_of_state_variables,Gsecond) ;
%estimate = zeros(sample_size,number_of_state_variables,3) ;
const_lik = (2*pi)^(.5*number_of_observed_variables) ;
ks = 0 ;
lik = NaN(sample_size,1);
LIK = NaN;
for t=1:sample_size
@ -154,13 +170,13 @@ for t=1:sample_size
for i=1:I
for j=1:J
for g=1:G ;
a = g + (j-1)*G ;
b = a + (i-1)*Gprime ;
[StateMuPrior(:,a),StateSqrtPPrior(:,:,a),StateWeightsPrior(1,a),...
StateMuPost(:,b),StateSqrtPPost(:,:,b),StateWeightsPost(1,b)] =...
gaussian_mixture_filter_bank(ReducedForm,Y(:,t),StateMu(:,g),StateSqrtP(:,:,g),StateWeights(1,g),...
StructuralShocksMu(:,i),StructuralShocksSqrtP(:,:,i),StructuralShocksWeights(1,i),...
ObservationShocksMu(:,j),ObservationShocksSqrtP(:,:,j),ObservationShocksWeights(1,j),...
gprime = g + (i-1)*G ;
gsecond = gprime + (j-1)*Gprime ;
[StateMuPrior(:,gprime),StateSqrtPPrior(:,:,gprime),StateWeightsPrior(1,gprime),...
StateMuPost(:,gsecond),StateSqrtPPost(:,:,gsecond),StateWeightsPost(1,gsecond)] =...
gaussian_mixture_filter_bank(ReducedForm,Y(:,t),StateMu(:,g),StateSqrtP(:,:,g),StateWeights(g),...
StructuralShocksMu(:,i),StructuralShocksSqrtP(:,:,i),StructuralShocksWeights(i),...
ObservationShocksMu(:,j),ObservationShocksSqrtP(:,:,j),ObservationShocksWeights(j),...
H,H_lower_triangular_cholesky,const_lik,ParticleOptions,ThreadsOptions) ;
end
end
@ -197,29 +213,11 @@ for t=1:sample_size
StateParticles,H,const_lik,1/number_of_particles,...
1/number_of_particles,ReducedForm,ThreadsOptions) ;
% calculate importance weights of particles
% SampleWeights = SampleWeights.*IncrementalWeights ;
SampleWeights = IncrementalWeights/number_of_particles ;
SumSampleWeights = sum(SampleWeights,1) ;
SampleWeights = SampleWeights./SumSampleWeights ;
lik(t) = log(SumSampleWeights) ;
% First possible state point estimates
%estimate(t,:,1) = SampleWeights*StateParticles' ;
% Resampling if needed of required
% Neff = neff(SampleWeights) ;
% if (ParticleOptions.resampling.status.generic && Neff<.5*sample_size) || ParticleOptions.resampling.status.systematic
% ks = ks + 1 ;
% StateParticles = resample(StateParticles',SampleWeights,ParticleOptions)' ;
% StateVectorMean = mean(StateParticles,2) ;
% StateVectorVarianceSquareRoot = reduced_rank_cholesky( (StateParticles*StateParticles')/number_of_particles - StateVectorMean*(StateVectorMean') )';
% SampleWeights = 1/number_of_particles ;
% elseif ParticleOptions.resampling.status.none
% StateVectorMean = StateParticles*sampleWeights ;
% temp = sqrt(SampleWeights').*StateParticles ;
% StateVectorVarianceSquareRoot = reduced_rank_cholesky( temp*temp' - StateVectorMean*(StateVectorMean') )';
% end
% Use the information from particles to update the gaussian mixture on state variables
[StateMu,StateSqrtP,StateWeights] = fit_gaussian_mixture(StateParticles,StateMu,StateSqrtP,StateWeights,0.001,10,1) ;
%estimate(t,:,3) = StateWeights*StateMu' ;
end
end