File: posterior_moments.m

package info (click to toggle)
dynare 4.5.7-1
  • links: PTS, VCS
  • area: main
  • in suites: buster
  • size: 49,408 kB
  • sloc: cpp: 84,998; ansic: 29,058; pascal: 13,843; sh: 4,833; objc: 4,236; yacc: 3,622; makefile: 2,278; lex: 1,541; python: 236; lisp: 69; xml: 8
file content (96 lines) | stat: -rw-r--r-- 3,764 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
function [post_mean, post_median, post_var, hpd_interval, post_deciles, density] = posterior_moments(xx,info,mh_conf_sig,kernel_options)
% Computes posterior mean, median, variance, HPD interval, deciles, and density from posterior draws.
%
% INPUTS
%    xx            [double]    Vector of posterior draws (or prior draws ;-)
%    info          [integer]   If equal to one the posterior density is estimated.
%    mh_config_sig [double]    Scalar between 0 and 1 specifying the size of the HPD interval.
%
%
% OUTPUTS
%    post_mean     [double]    Scalar, posterior mean.
%    post_median   [double]    Scalar, posterior median.
%    post_var      [double]    Scalar, posterior variance.
%    hpd_interval  [double]    Vector (1*2), Highest Probability Density interval
%    post_deciles  [double]    Vector (9*1), deciles of the posterior distribution.
%    density       [double]    Matrix (n*2), non parametric estimate of the posterior density. First and second
%                              columns are respectively abscissa and ordinate coordinates.
%
% SPECIAL REQUIREMENTS
%    Other matlab routines distributed with Dynare: mh_optimal_bandwidth.m
%                                                   kernel_density_estimate.m.
%

% Copyright (C) 2005-2017 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 <http://www.gnu.org/licenses/>.

if nargin<4
    number_of_grid_points = 2^9;      % 2^9 = 512 !... Must be a power of two.
    bandwidth = 0;                    % Rule of thumb optimal bandwidth parameter.
    kernel_function = 'gaussian';     % Gaussian kernel for Fast Fourrier Transform approximaton.
else
    number_of_grid_points = kernel_options.gridpoints;
    bandwidth = kernel_options.bandwidth;
    kernel_function = kernel_options.kernel_function;
end

xx = xx(:);
xx = sort(xx);


post_mean = mean(xx);
post_median = median(xx);
post_var = var(xx);

number_of_draws = length(xx);
hpd_draws = round((1-mh_conf_sig)*number_of_draws);

if hpd_draws>2
    kk = zeros(hpd_draws,1);
    jj = number_of_draws-hpd_draws;
    for ii = 1:hpd_draws
        kk(ii) = xx(jj)-xx(ii);
        jj = jj + 1;
    end
    [kmin,idx] = min(kk);
    hpd_interval = [xx(idx) xx(idx)+kmin];
else
    hpd_interval=NaN(1,2);
end
if length(xx)>9
    post_deciles = xx([round(0.1*number_of_draws) ...
                       round(0.2*number_of_draws) ...
                       round(0.3*number_of_draws) ...
                       round(0.4*number_of_draws) ...
                       round(0.5*number_of_draws) ...
                       round(0.6*number_of_draws) ...
                       round(0.7*number_of_draws) ...
                       round(0.8*number_of_draws) ...
                       round(0.9*number_of_draws)]);
else
    post_deciles=NaN(9,1);
end
density = [];
if info
    if post_var > 1e-12
        optimal_bandwidth = mh_optimal_bandwidth(xx,number_of_draws,bandwidth,kernel_function);
        [density(:,1),density(:,2)] = kernel_density_estimate(xx,number_of_grid_points,...
                                                          number_of_draws,optimal_bandwidth,kernel_function);
    else
        density = NaN(number_of_grid_points,2);
    end
end