1function [Gamma_y,stationary_vars] = th_autocovariances(dr,ivar,M_,options_,nodecomposition)
2% Computes the theoretical auto-covariances, Gamma_y, for an AR(p) process
3% with coefficients dr.ghx and dr.ghu and shock variances Sigma_e
4% for a subset of variables ivar.
5% Theoretical HP-filtering and band-pass filtering is available as an option
6%
7% INPUTS
8%   dr:               [structure]    Reduced form solution of the DSGE model  (decisions rules)
9%   ivar:             [integer]      Vector of indices for a subset of variables.
10%   M_                [structure]    Global dynare's structure, description of the DSGE model.
11%   options_          [structure]    Global dynare's structure.
12%   nodecomposition   [integer]      Scalar, if different from zero the variance decomposition is not triggered.
13%
14% OUTPUTS
15%   Gamma_y           [cell]         Matlab cell of nar+3 (second order approximation) or nar+2 (first order approximation) arrays,
16%                                    where nar is the order of the autocorrelation function.
17%                                      Gamma_y{1}       [double]  Covariance matrix.
18%                                      Gamma_y{i+1}     [double]  Autocorrelation function (for i=1,...,options_.nar).
19%                                      Gamma_y{nar+2}   [double]  Variance decomposition.
20%                                      Gamma_y{nar+3}   [double]  Expectation of the endogenous variables associated with a second
21%                                                                 order approximation.
22%   stationary_vars   [integer]      Vector of indices of stationary variables (as a subset of 1:length(ivar))
23%
24% SPECIAL REQUIREMENTS
25%
26% Algorithms
27%   The means at order=2 are based on the pruned state space as
28%   in Kim, Kim, Schaumburg, Sims (2008): Calculating and using second-order accurate
29%   solutions of discrete time dynamic equilibrium models.
30%   The solution at second order can be written as:
31%   \[
32%   \hat x_t = g_x \hat x_{t - 1} + g_u u_t + \frac{1}{2}\left( g_{\sigma\sigma} \sigma^2 + g_{xx}\hat x_t^2 + g_{uu} u_t^2 \right)
33%   \]
34%   Taking expectations on both sides requires to compute E(x^2)=Var(x), which
35%   can be obtained up to second order from the first order solution
36%   \[
37%       \hat x_t = g_x \hat x_{t - 1} + g_u u_t
38%   \]
39%   by solving the corresponding Lyapunov equation.
40%   Given Var(x), the above equation can be solved for E(x_t) as
41%   \[
42%   E(x_t) = (I - {g_x}\right)^{- 1} 0.5\left( g_{\sigma\sigma} \sigma^2 + g_{xx} Var(\hat x_t) + g_{uu} Var(u_t) \right)
43%   \]
44%
45% Copyright (C) 2001-2020 Dynare Team
46%
47% This file is part of Dynare.
48%
49% Dynare is free software: you can redistribute it and/or modify
50% it under the terms of the GNU General Public License as published by
51% the Free Software Foundation, either version 3 of the License, or
52% (at your option) any later version.
53%
54% Dynare is distributed in the hope that it will be useful,
55% but WITHOUT ANY WARRANTY; without even the implied warranty of
56% MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
57% GNU General Public License for more details.
58%
59% You should have received a copy of the GNU General Public License
60% along with Dynare.  If not, see <http://www.gnu.org/licenses/>.
61
62if nargin<5
63    nodecomposition = 0;
64end
65
66if options_.order >= 3
67    error('Theoretical moments not implemented above 2nd order')
68end
69
70local_order = options_.order;
71if local_order~=1 && M_.hessian_eq_zero
72    local_order = 1;
73end
74
75endo_nbr = M_.endo_nbr;
76exo_names_orig_ord  = M_.exo_names_orig_ord;
77if isoctave
78    warning('off', 'Octave:divide-by-zero')
79else
80    warning off MATLAB:dividebyzero
81end
82nar = options_.ar;
83Gamma_y = cell(nar+2,1);
84if isempty(ivar)
85    ivar = [1:endo_nbr]';
86end
87nvar = size(ivar,1);
88
89ghx = dr.ghx;
90ghu = dr.ghu;
91nspred = M_.nspred;
92nstatic = M_.nstatic;
93
94nx = size(ghx,2);
95if ~options_.block
96    %order_var = dr.order_var;
97    inv_order_var = dr.inv_order_var;
98    kstate = dr.kstate;
99    ikx = [nstatic+1:nstatic+nspred];
100    k0 = kstate(find(kstate(:,2) <= M_.maximum_lag+1),:);
101    i0 = find(k0(:,2) == M_.maximum_lag+1);
102    i00 = i0;
103    n0 = length(i0);
104    AS = ghx(:,i0);
105    ghu1 = zeros(nx,M_.exo_nbr);
106    ghu1(i0,:) = ghu(ikx,:);
107    for i=M_.maximum_lag:-1:2
108        i1 = find(k0(:,2) == i);
109        n1 = size(i1,1);
110        j1 = zeros(n1,1);
111        for k1 = 1:n1
112            j1(k1) = find(k0(i00,1)==k0(i1(k1),1));
113        end
114        AS(:,j1) = AS(:,j1)+ghx(:,i1);
115        i0 = i1;
116    end
117else
118    ghu1 = zeros(nx,M_.exo_nbr);
119    trend = 1:M_.endo_nbr;
120    inv_order_var = trend(M_.block_structure.variable_reordered);
121    ghu1(1:length(dr.state_var),:) = ghu(dr.state_var,:);
122end
123b = ghu1*M_.Sigma_e*ghu1';
124
125
126if ~options_.block
127    ipred = nstatic+(1:nspred)';
128else
129    ipred = dr.state_var;
130end
131% state space representation for state variables only
132[A,B] = kalman_transition_matrix(dr,ipred,1:nx,M_.exo_nbr);
133% Compute stationary variables (before HP filtering),
134% and compute 2nd order mean correction on stationary variables (in case of
135% HP filtering, this mean correction is computed *before* filtering)
136if local_order == 2 || options_.hp_filter == 0
137    [vx, u] =  lyapunov_symm(A,B*M_.Sigma_e*B',options_.lyapunov_fixed_point_tol,options_.qz_criterium,options_.lyapunov_complex_threshold,[],options_.debug);
138    if ~options_.block
139        iky = inv_order_var(ivar);
140    else
141        iky = ivar;
142    end
143    stationary_vars = (1:length(ivar))';
144    if ~isempty(u)
145        x = abs(ghx*u);
146        iky = iky(find(all(x(iky,:) < options_.Schur_vec_tol,2)));
147        stationary_vars = find(all(x(inv_order_var(ivar(stationary_vars)),:) < options_.Schur_vec_tol,2));
148    end
149    aa = ghx(iky,:);
150    bb = ghu(iky,:);
151    if local_order == 2         % mean correction for 2nd order
152        if ~isempty(ikx)
153            Ex = (dr.ghs2(ikx)+dr.ghxx(ikx,:)*vx(:)+dr.ghuu(ikx,:)*M_.Sigma_e(:))/2;
154            Ex = (eye(n0)-AS(ikx,:))\Ex;
155            Gamma_y{nar+3} = NaN*ones(nvar, 1);
156            Gamma_y{nar+3}(stationary_vars) = AS(iky,:)*Ex+(dr.ghs2(iky)+dr.ghxx(iky,:)*vx(:)+...
157                                                            dr.ghuu(iky,:)*M_.Sigma_e(:))/2;
158        else %no static and no predetermined
159            Gamma_y{nar+3} = NaN*ones(nvar, 1);
160            Gamma_y{nar+3}(stationary_vars) = (dr.ghs2(iky)+ dr.ghuu(iky,:)*M_.Sigma_e(:))/2;
161        end
162    end
163end
164if options_.hp_filter == 0 && ~options_.bandpass.indicator
165    v = NaN*ones(nvar,nvar);
166    v(stationary_vars,stationary_vars) = aa*vx*aa'+ bb*M_.Sigma_e*bb';
167    k = find(abs(v) < 1e-12);
168    v(k) = 0;
169    Gamma_y{1} = v;
170    % autocorrelations
171    if nar > 0
172        vxy = (A*vx*aa'+ghu1*M_.Sigma_e*bb');
173        sy = sqrt(diag(Gamma_y{1}));
174        sy = sy(stationary_vars);
175        sy = sy *sy';
176        Gamma_y{2} = NaN*ones(nvar,nvar);
177        Gamma_y{2}(stationary_vars,stationary_vars) = aa*vxy./sy;
178        for i=2:nar
179            vxy = A*vxy;
180            Gamma_y{i+1} = NaN*ones(nvar,nvar);
181            Gamma_y{i+1}(stationary_vars,stationary_vars) = aa*vxy./sy;
182        end
183    end
184    % variance decomposition
185    if ~nodecomposition && M_.exo_nbr > 0 && size(stationary_vars, 1) > 0
186        if M_.exo_nbr == 1
187            Gamma_y{nar+2} = ones(nvar,1);
188        else
189            Gamma_y{nar+2} = NaN(nvar,M_.exo_nbr);
190            SS(exo_names_orig_ord,exo_names_orig_ord)=M_.Sigma_e+1e-14*eye(M_.exo_nbr);
191            cs = chol(SS)';
192            b1(:,exo_names_orig_ord) = ghu1;
193            b1 = b1*cs;
194            b2(:,exo_names_orig_ord) = ghu(iky,:);
195            b2 = b2*cs;
196            vx  = lyapunov_symm(A,b1*b1',options_.lyapunov_fixed_point_tol,options_.qz_criterium,options_.lyapunov_complex_threshold,1,options_.debug);
197            vv = diag(aa*vx*aa'+b2*b2');
198            vv2 = 0;
199            for i=1:M_.exo_nbr
200                vx1 = lyapunov_symm(A,b1(:,i)*b1(:,i)',options_.lyapunov_fixed_point_tol,options_.qz_criterium,options_.lyapunov_complex_threshold,2,options_.debug);
201                vx2 = abs(diag(aa*vx1*aa'+b2(:,i)*b2(:,i)'));
202                Gamma_y{nar+2}(stationary_vars,i) = vx2;
203                vv2 = vv2 +vx2;
204            end
205            if max(abs(vv2-vv)./vv) > 1e-4
206                warning(['Aggregate variance and sum of variances by shocks ' ...
207                         'differ by more than 0.01 %'])
208            end
209            for i=1:M_.exo_nbr
210                Gamma_y{nar+2}(stationary_vars,i) = Gamma_y{nar+ ...
211                                    2}(stationary_vars,i)./vv2;
212            end
213        end
214    end
215else% ==> Theoretical filters.
216    % By construction, all variables are stationary when HP filtered
217    iky = inv_order_var(ivar);
218    stationary_vars = (1:length(ivar))';
219    aa = ghx(iky,:); %R in Uhlig (2001)
220    bb = ghu(iky,:); %S in Uhlig (2001)
221
222    lambda = options_.hp_filter;
223    ngrid = options_.filtered_theoretical_moments_grid;
224    freqs = 0 : ((2*pi)/ngrid) : (2*pi*(1 - .5/ngrid)); %[0,2*pi)
225    tpos  = exp( sqrt(-1)*freqs); %positive frequencies
226    tneg  =  exp(-sqrt(-1)*freqs); %negative frequencies
227    if options_.bandpass.indicator
228        filter_gain = zeros(1,ngrid);
229        lowest_periodicity=options_.bandpass.passband(2);
230        highest_periodicity=options_.bandpass.passband(1);
231        highest_periodicity=max(2,highest_periodicity); % restrict to upper bound of pi
232        filter_gain(freqs>=2*pi/lowest_periodicity & freqs<=2*pi/highest_periodicity)=1;
233        filter_gain(freqs<=-2*pi/lowest_periodicity+2*pi & freqs>=-2*pi/highest_periodicity+2*pi)=1;
234    else
235        filter_gain = 4*lambda*(1 - cos(freqs)).^2 ./ (1 + 4*lambda*(1 - cos(freqs)).^2);   %HP transfer function
236    end
237    mathp_col = NaN(ngrid,length(ivar)^2);
238    IA = eye(size(A,1));
239    IE = eye(M_.exo_nbr);
240    for ig = 1:ngrid
241        if filter_gain(ig)==0
242            f_hp = zeros(length(ivar),length(ivar));
243        else
244            f_omega  =(1/(2*pi))*([(IA-A*tneg(ig))\ghu1;IE]...
245                                  *M_.Sigma_e*[ghu1'/(IA-A'*tpos(ig)) IE]); % spectral density of state variables; top formula Uhlig (2001), p. 20 with N=0
246            g_omega = [aa*tneg(ig) bb]*f_omega*[aa'*tpos(ig); bb']; % spectral density of selected variables; middle formula Uhlig (2001), p. 20; only middle block, i.e. y_t'
247            f_hp = filter_gain(ig)^2*g_omega; % spectral density of selected filtered series; top formula Uhlig (2001), p. 21;
248        end
249        mathp_col(ig,:) = (f_hp(:))';    % store as matrix row for ifft
250    end
251    % Covariance of filtered series
252    imathp_col = real(ifft(mathp_col))*(2*pi); % Inverse Fast Fourier Transformation; middle formula Uhlig (2001), p. 21;
253    Gamma_y{1} = reshape(imathp_col(1,:),nvar,nvar);
254    % Autocorrelations
255    if nar > 0
256        sy = sqrt(diag(Gamma_y{1}));
257        sy = sy *sy';
258        for i=1:nar
259            Gamma_y{i+1} = reshape(imathp_col(i+1,:),nvar,nvar)./sy;
260        end
261    end
262    % Variance decomposition
263    if ~nodecomposition && M_.exo_nbr > 0
264        if M_.exo_nbr == 1
265            Gamma_y{nar+2} = ones(nvar,1);
266        else
267            Gamma_y{nar+2} = zeros(nvar,M_.exo_nbr);
268            SS(exo_names_orig_ord,exo_names_orig_ord) = M_.Sigma_e+1e-14*eye(M_.exo_nbr); %make sure Covariance matrix is positive definite
269            cs = chol(SS)';
270            SS = cs*cs';
271            b1(:,exo_names_orig_ord) = ghu1;
272            b2(:,exo_names_orig_ord) = ghu(iky,:);
273            mathp_col = NaN(ngrid,length(ivar)^2);
274            IA = eye(size(A,1));
275            IE = eye(M_.exo_nbr);
276            for ig = 1:ngrid
277                if filter_gain(ig)==0
278                    f_hp = zeros(length(ivar),length(ivar));
279                else
280                    f_omega  =(1/(2*pi))*( [(IA-A*tneg(ig))\b1;IE]...
281                                           *SS*[b1'/(IA-A'*tpos(ig)) IE]); % spectral density of state variables; top formula Uhlig (2001), p. 20 with N=0
282                    g_omega = [aa*tneg(ig) b2]*f_omega*[aa'*tpos(ig); b2']; % spectral density of selected variables; middle formula Uhlig (2001), p. 20; only middle block, i.e. y_t'
283                    f_hp = filter_gain(ig)^2*g_omega;  % spectral density of selected filtered series; top formula Uhlig (2001), p. 21;
284                end
285                mathp_col(ig,:) = (f_hp(:))';    % store as matrix row for ifft
286            end
287            imathp_col = real(ifft(mathp_col))*(2*pi);
288            vv = diag(reshape(imathp_col(1,:),nvar,nvar));
289            for i=1:M_.exo_nbr
290                mathp_col = NaN(ngrid,length(ivar)^2);
291                SSi = cs(:,i)*cs(:,i)';
292                for ig = 1:ngrid
293                    if filter_gain(ig)==0
294                        f_hp = zeros(length(ivar),length(ivar));
295                    else
296                        f_omega  =(1/(2*pi))*( [(IA-A*tneg(ig))\b1;IE]...
297                                               *SSi*[b1'/(IA-A'*tpos(ig)) IE]); % spectral density of state variables; top formula Uhlig (2001), p. 20 with N=0
298                        g_omega = [aa*tneg(ig) b2]*f_omega*[aa'*tpos(ig); b2']; % spectral density of selected variables; middle formula Uhlig (2001), p. 20; only middle block, i.e. y_t'
299                        f_hp = filter_gain(ig)^2*g_omega; % spectral density of selected filtered series; top formula Uhlig (2001), p. 21;
300                    end
301                    mathp_col(ig,:) = (f_hp(:))';    % store as matrix row for ifft
302                end
303                imathp_col = real(ifft(mathp_col))*(2*pi);
304                Gamma_y{nar+2}(:,i) = abs(diag(reshape(imathp_col(1,:),nvar,nvar)))./vv;
305            end
306        end
307    end
308end
309if isoctave
310    warning('on', 'Octave:divide-by-zero')
311else
312    warning on MATLAB:dividebyzero
313end
314