File: perfect_foresight_solver_core.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 (153 lines) | stat: -rw-r--r-- 6,987 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
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
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
function [y, success, maxerror, iter, per_block_status] = perfect_foresight_solver_core(y, exo_simul, steady_state, exo_steady_state, M_, options_)

% Core function calling solvers for perfect foresight model
%
% INPUTS
% - y                   [matrix] initial path of endogenous (typically oo_.endo_simul)
% - exo_simul           [matrix] path of exogenous
% - steady_state        [vector] steady state of endogenous variables
% - exo_steady_state    [vector] steady state of exogenous variables
% - M_                  [struct] contains a description of the model.
% - options_            [struct] contains various options.
%
% OUTPUTS
% - y                   [double array] path for the endogenous variables (solution)
% - success             [logical] Whether a solution was found
% - maxerror            [double] contains the maximum absolute error
% - iter                [integer] Number of iterations of the underlying nonlinear solver (empty for non-iterative methods)
% - per_block_status    [struct] In the case of block decomposition, provides per-block solver status information (empty if no block decomposition)

% Copyright © 2015-2023 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/>.

if options_.lmmcp.status
    options_.stack_solve_algo=7;
    options_.solve_algo = 10;
elseif options_.stack_solve_algo==7 && options_.solve_algo == 11
    options_.lmmcp.status = 1; %Path solver
end

periods = options_.periods;

if options_.linear_approximation && ~(isequal(options_.stack_solve_algo,0) || isequal(options_.stack_solve_algo,7))
    error('perfect_foresight_solver: Option linear_approximation is only available with option stack_solve_algo equal to 0 or 7.')
end

if options_.endogenous_terminal_period && options_.stack_solve_algo ~= 0
    error('perfect_foresight_solver: option endogenous_terminal_period is only available with option stack_solve_algo equal to 0')
end

if options_.linear && (isequal(options_.stack_solve_algo, 0) || isequal(options_.stack_solve_algo, 7))
    options_.linear_approximation = true;
end

maxerror = [];
iter = [];
per_block_status = [];

if options_.block
    if M_.block_structure.time_recursive
        error('Internal error: can''t perform stacked perfect foresight simulation with time-recursive block decomposition')
    end
    if options_.bytecode && ~ismember(options_.stack_solve_algo, [1 6])
        try
            y = bytecode('dynamic', 'block_decomposed', M_, options_, y, exo_simul, M_.params, steady_state, periods);
            success = true;
        catch ME
            if options_.verbosity >= 1
                disp(ME.message)
            end
            success = false;
        end
    else
        [y, success, maxerror, per_block_status] = solve_block_decomposed_problem(y, exo_simul, steady_state, options_, M_);
    end
else
    if options_.bytecode && ~ismember(options_.stack_solve_algo, [1 6])
        try
            y = bytecode('dynamic', M_, options_, y, exo_simul, M_.params, repmat(steady_state, 1, periods+2), periods);
            success = true;
        catch ME
            if options_.verbosity >= 1
                disp(ME.message)
            end
            success = false;
        end
    else
        if M_.maximum_endo_lead == 0 && M_.maximum_endo_lag>0 && ~options_.lmmcp.status % Purely backward model
            [y, success] = sim1_purely_backward(y, exo_simul, steady_state, M_, options_);
        elseif M_.maximum_endo_lag == 0 && M_.maximum_endo_lead>0 && ~options_.lmmcp.status % Purely forward model
            [y, success] = sim1_purely_forward(y, exo_simul, steady_state, M_, options_);
        elseif M_.maximum_endo_lag == 0 && M_.maximum_endo_lead == 0 && ~options_.lmmcp.status % Purely static model
            [y, success] = sim1_purely_static(y, exo_simul, steady_state, M_, options_);
        else % General case
            switch options_.stack_solve_algo
              case 0
                if options_.linear_approximation
                    [y, success, maxerror] = sim1_linear(y, exo_simul, steady_state, exo_steady_state, M_, options_);
                else
                    [y, success, maxerror, iter] = sim1(y, exo_simul, steady_state, M_, options_);
                end
              case {1 6}
                if options_.linear_approximation
                    error('Invalid value of stack_solve_algo option!')
                end
                [y, success, maxerror, iter] = sim1_lbj(y, exo_simul, steady_state, M_, options_);
              case 7
                if options_.linear_approximation
                    if isequal(options_.solve_algo, 10) 
                        if options_.ramsey_policy && isfield(M_,'ramsey_model_constraints') && ~isempty(M_.ramsey_model_constraints)
                            warning('Due to ramsey_constraints you should not specify your model as model(linear)!')
                        elseif options_.lmmcp.status
                            warning('Due to lmmcp option, you should not specify your model as model(linear)!')
                        else
                            warning('It would be more efficient to set option solve_algo equal to 0!')
                        end
                    end
                    [y, success] = solve_stacked_linear_problem(y, exo_simul, steady_state, exo_steady_state, M_, options_);
                else
                    [y, success, maxerror] = solve_stacked_problem(y, exo_simul, steady_state, M_, options_);
                end
              otherwise
                error('Invalid value of stack_solve_algo option!')
            end
        end
    end
end

% Some solvers do not compute the maximum error, so do it here if needed
if nargout > 2 && isempty(maxerror)
    if options_.bytecode
        residuals = bytecode('dynamic', 'evaluate', M_, options_, y, exo_simul, M_.params, steady_state, periods);
    else
        ny = size(y, 1);
        if M_.maximum_lag > 0
            y0 = y(:, M_.maximum_lag);
        else
            y0 = NaN(ny, 1);
        end
        if M_.maximum_lead > 0
            yT = y(:, M_.maximum_lag+periods+1);
        else
            yT = NaN(ny, 1);
        end
        yy = y(:,M_.maximum_lag+(1:periods));

        residuals = perfect_foresight_problem(yy(:), y0, yT, exo_simul, M_.params, steady_state, periods, M_, options_);
    end
    maxerror = max(max(abs(residuals)));
end