File: write.m

package info (click to toggle)
dynare 4.6.3-4
  • links: PTS, VCS
  • area: main
  • in suites: bullseye
  • size: 74,896 kB
  • sloc: cpp: 98,057; ansic: 28,929; pascal: 13,844; sh: 5,947; objc: 4,236; yacc: 4,215; makefile: 2,583; lex: 1,534; fortran: 877; python: 647; ruby: 291; lisp: 152; xml: 22
file content (331 lines) | stat: -rw-r--r-- 13,684 bytes parent folder | download
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
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
function write(DynareModel)

% Writes the nonlinear problem to be solved for computing the growth
% rates and levels along the Balanced Growth Path. Note that for
% the variables growing along the BGP, the identified levels are
% meaningless, only the relative levels of the growing variables
% sharing the same trend(s) are relevant.
%
% INPUTS
% - DynareModel     [struct]   Dynare generated stucture describing the model (M_).
%
% OUTPUTS
% None
%
% REMARKS
% - The trends are assumed to be multiplicative.

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

if DynareModel.maximum_lag  && ~DynareModel.maximum_lead
    i0 = transpose(DynareModel.lead_lag_incidence(1,:)); % Indices of the lagged variables.
    i1 = transpose(DynareModel.lead_lag_incidence(2,:)); % Indices of the current variables.
    i2 = [];                                             % Indices of the leaded variables.
elseif DynareModel.maximum_lag  && DynareModel.maximum_lead
    i0 = transpose(DynareModel.lead_lag_incidence(1,:)); % Indices of the lagged variables.
    i1 = transpose(DynareModel.lead_lag_incidence(2,:)); % Indices of the current variables.
    i2 = transpose(DynareModel.lead_lag_incidence(3,:)); % Indices of the leaded variables.
elseif ~DynareModel.maximum_lag  && DynareModel.maximum_lead
    i0 = [];                                             % Indices of the lagged variables.
    i1 = transpose(DynareModel.lead_lag_incidence(1,:)); % Indices of the current variables.
    i2 = transpose(DynareModel.lead_lag_incidence(2,:)); % Indices of the leaded variables.
else
    error('The model is static. The BGP is trivial.')
end

n0 = length(find(i0)); % Number of lagged variables.
n1 = length(find(i1)); % Number of current variables.
n2 = length(find(i2)); % Number of leaded variables.

purely_backward_model = logical(n0 && n1 && ~n2);
purely_forward_model = logical(~n0 && n1 && n2);

if purely_backward_model
    I0 = i0;
    I1 = i1;
    I2 = [];
elseif purely_forward_model
    I0 = i1;
    I1 = i2;
    I2 = [];
else
    % Model has both leads and lags.
    I0 = i0;
    I1 = i1;
    I2 = i2;
end

% Create function in mod namespace.
fid = fopen(sprintf('+%s/bgpfun.m', DynareModel.fname), 'w');

% Write header. 
fprintf(fid, 'function [F, JAC] = bgpfun(z)\n\n');
fprintf(fid, '%% This file has been generated by dynare (%s).\n\n', datestr(now));

% The function admits a unique vector as input argument. The first
% half of the elements are for the levels of the endogenous
% variables, the second half for the growth factors. 
fprintf(fid, 'y = z(1:%u);\n\n', DynareModel.endo_nbr);
fprintf(fid, 'g = z(%u:%u);\n', DynareModel.endo_nbr+1, 2*DynareModel.endo_nbr);

% Define the point where the dynamic model is to be evaluated.
fprintf(fid, 'Y = zeros(%u, 1);\n', 2*(n0+n1+n2));
for i=1:length(I0) % period t equations, lagged variables.
    if I0(i)
        fprintf(fid, 'Y(%u) = y(%u);\n', I0(i), i);
    end
end
for i=1:length(I1) % period t equations, current variables.
    if I1(i)
        fprintf(fid, 'Y(%u) = y(%u)*g(%u);\n', I1(i), i, i);
    end
end
for i=1:length(I2) % period t equations, leaded variables.
    if I2(i)
        fprintf(fid, 'Y(%u) = y(%u)*g(%u)*g(%u);\n', I2(i), i, i, i);
    end
end
for i=1:length(I0) % period t+1 equations lagged variables.
    if I0(i)
        fprintf(fid, 'Y(%u) = y(%u)*g(%u);\n', n0+n1+n2+I0(i), i, i);
    end
end
for i=1:length(I1) % period t+1 equations current variables.
    if I1(i)
        fprintf(fid, 'Y(%u) = y(%u)*g(%u)*g(%u);\n', n0+n1+n2+I1(i), i, i, i);
    end
end
for i=1:length(I2) % period t+1 equations leaded variables.
    if I2(i)
        fprintf(fid, 'Y(%u) = y(%u)*g(%u)*g(%u)*g(%u);\n', n0+n1+n2+I2(i), i, i, i, i);
    end
end
fprintf(fid, '\n');

% Define the vector of parameters.
fprintf(fid, 'p = zeros(%u, 1);\n', DynareModel.param_nbr);
for i = 1:DynareModel.param_nbr
    fprintf(fid, 'p(%u) = %16.12f;\n', i, DynareModel.params(i));
end
fprintf(fid, '\n');

% Initialize the vector holding the residuals over the two periods.
fprintf(fid, 'F = NaN(%u, 1);\n', 2*DynareModel.endo_nbr);

% Set vector of exogenous variables to 0.
fprintf(fid, 'x = zeros(1, %u);\n\n', DynareModel.exo_nbr);

% Evaluate the residuals and jacobian of the dynamic model in periods t and t+1.
fprintf(fid, 'if nargout>1\n');
fprintf(fid, '    J = zeros(%u, %u);\n', 2*DynareModel.endo_nbr, n0+n1+n2+DynareModel.endo_nbr);
fprintf(fid, '    [F(1:%u), tmp] = %s.dynamic(Y(1:%u), x, p, y, 1);\n', DynareModel.endo_nbr, DynareModel.fname, n0+n1+n2);
fprintf(fid, '    J(1:%u,1:%u) = tmp(:,1:%u);\n', DynareModel.endo_nbr, n0+n1+n2, n0+n1+n2);
fprintf(fid, '    [F(%u:%u), tmp] = %s.dynamic(Y(1+%u:%u), x, p, y, 1);\n', DynareModel.endo_nbr+1, 2*DynareModel.endo_nbr, DynareModel.fname, n0+n1+n2, 2*(n0+n1+n2));
fprintf(fid, '    J(%u:%u,1:%u) = tmp(:,1:%u);\n', DynareModel.endo_nbr+1, 2*DynareModel.endo_nbr, n0+n1+n2, n0+n1+n2);
fprintf(fid, 'else\n');
fprintf(fid, '    F(1:%u) = %s.dynamic(Y(1:%u), x, p, y, 1);\n', DynareModel.endo_nbr, DynareModel.fname, n0+n1+n2);
fprintf(fid, '    F(%u:%u) = %s.dynamic(Y(1+%u:%u), x, p, y, 1);\n', DynareModel.endo_nbr+1, 2*DynareModel.endo_nbr, DynareModel.fname, n0+n1+n2, 2*(n0+n1+n2));
fprintf(fid, 'end\n\n');

% Compute the jacobian if required.
fprintf(fid, 'if nargout>1\n');
fprintf(fid, '    JAC = zeros(%u,%u);\n', 2*DynareModel.endo_nbr, 2*DynareModel.endo_nbr);

% Compute the derivatives of the first block of equations (period t)
% with respect to the endogenous variables.
if purely_backward_model || purely_forward_model
    for i=1:DynareModel.eq_nbr
        for j=1:DynareModel.endo_nbr
            if I1(j)
                if I0(j)
                    fprintf(fid, '    JAC(%u,%u) = J(%u,%u)+J(%u,%u)*g(%u);\n', i, j, i, I0(j), i, I1(j), j);
                else
                    fprintf(fid, '    JAC(%u,%u) = J(%u,%u)*g(%u);\n', i, j, i, I1(j), j);
                end
            else
                if I0(j)
                    fprintf(fid, '    JAC(%u,%u) = J(%u,%u);\n', i, j, i, I0(j));
                end
            end
        end
    end
else
    for i=1:DynareModel.eq_nbr
        for j=1:DynareModel.endo_nbr
            if I2(j)
                if I1(j)
                    if I0(j)
                        fprintf(fid, '    JAC(%u,%u) = J(%u,%u)+J(%u,%u)*g(%u)+J(%u,%u)*g(%u)*g(%u);\n', i, j, i, I0(j), i, I1(j), j, i, I2(j), j, j);
                    else
                        fprintf(fid, '    JAC(%u,%u) = J(%u,%u)*g(%u)+J(%u,%u)*g(%u)*g(%u);\n', i, j, i, I1(j), j, i, I2(j), j, j);
                    end
                else
                    if I0(j)
                        fprintf(fid, '    JAC(%u,%u) = J(%u,%u)+J(%u,%u)*g(%u)*g(%u);\n', i, j, i, I0(j), i, I2(j), j, j);
                    else
                        fprintf(fid, '    JAC(%u,%u) = J(%u,%u)*g(%u)*g(%u);\n', i, j, i, I2(j), j, j);
                    end
                end
            else
                if I1(j)
                    if I0(j)
                        fprintf(fid, '    JAC(%u,%u) = J(%u,%u)+J(%u,%u)*g(%u);\n', i, j, i, I0(j), i, I1(j), j);
                    else
                        fprintf(fid, '    JAC(%u,%u) = J(%u,%u)*g(%u);\n', i, j, i, I1(j), j);
                    end
                else
                    if I0(j)
                        fprintf(fid, '    JAC(%u,%u) = J(%u,%u);\n', i, j, i, I0(j));
                    end
                end
            end
        end
    end
end

% Compute the derivatives of the second block of equations (period t+1)
% with respect to the endogenous variables.
if purely_backward_model || purely_forward_model
    for i=DynareModel.eq_nbr+1:2*DynareModel.eq_nbr
        for j=1:DynareModel.endo_nbr
            if I1(j)
                if I0(j)
                    fprintf(fid, '    JAC(%u,%u) = J(%u,%u)*g(%u)+J(%u,%u)*g(%u)*g(%u);\n', i, j, i, I0(j), j, i, I1(j), j, j);
                else
                    fprintf(fid, '    JAC(%u,%u) = J(%u,%u)*g(%u)*g(%u);\n', i, j, i, I1(j), j, j);
                end
            else
                if I0(j)
                    fprintf(fid, '    JAC(%u, %u) = J(%u, %u)*g(%u);\n', i, j, i, I0(j), j);
                end
            end
        end
    end
else
    for i=DynareModel.eq_nbr+1:2*DynareModel.eq_nbr
        for j=1:DynareModel.endo_nbr
            if I2(j)
               if I1(j)
                   if I0(j)
                       fprintf(fid, '    JAC(%u,%u) = J(%u,%u)*g(%u)+J(%u,%u)*g(%u)*g(%u)+J(%u,%u)*g(%u)*g(%u)*g(%u);\n', i, j, i, I0(j), j, i, I1(j), j, j, i, I2(j), j, j, j);
                   else
                       fprintf(fid, '    JAC(%u,%u) = J(%u,%u)*g(%u)*g(%u)+J(%u,%u)*g(%u)*g(%u)*g(%u);\n', i, j, i, I1(j), j, j, i, I2(j), j, j, j);
                   end
               else
                   if I0(j)
                       fprintf(fid, '    JAC(%u,%u) = J(%u,%u)*g(%u)+J(%u,%u)*g(%u)*g(%u)*g(%u);\n', i, j, i, I0(j), j, i, I2(j), j, j, j);
                   else
                       fprintf(fid, '    JAC(%u,%u) = J(%u,%u)*g(%u)*g(%u)*g(%u);\n', i, j, i, I2(j), j, j, j);
                   end
               end
            else
                if I1(j)
                    if I0(j)
                        fprintf(fid, '    JAC(%u,%u) = J(%u,%u)*g(%u)+J(%u,%u)*g(%u)*g(%u);\n', i, j, i, I0(j), j, i, I1(j), j, j);
                    else
                        fprintf(fid, '    JAC(%u,%u) = J(%u,%u)*g(%u)*g(%u);\n', i, j, i, I1(j), j, j);
                    end
                else
                    if I0(j)
                        fprintf(fid, '    JAC(%u, %u) = J(%u, %u)*g(%u);\n', i, j, i, I0(j), j);
                    end
                end
            end
        end
    end
end

% Compute the derivatives of the first block of equations (period t)
% with respect to the growth factors.
if purely_backward_model || purely_forward_model
    for i=1:DynareModel.eq_nbr
        for j=1:DynareModel.endo_nbr
            if I1(j)
                fprintf(fid, '    J(%u,%u) = J(%u,%u)*y(%u);\n', i, n0+n1+n2+j, i, I1(j), j);
            end
        end
    end
else
    for i=1:DynareModel.eq_nbr
        for j=1:DynareModel.endo_nbr
            if I2(j)
                if I1(j)
                    fprintf(fid, '    J(%u,%u) = J(%u,%u)*y(%u)+J(%u,%u)*2*g(%u)*y(%u);\n', i, n0+n1+n2+j, i, I1(j), j, i, I2(j), j, j);
                else
                    fprintf(fid, '    J(%u,%u) = J(%u,%u)*2*g(%u)*y(%u);\n', i, n0+n1+n2+j, i, I2(j), j, j);
                end
            else
                if I1(j)
                    fprintf(fid, '    J(%u,%u) = J(%u,%u)*y(%u);\n', i, n0+n1+n2+j, i, I1(j), j);
                end
            end
        end
    end
end

% Compute the derivatives of the second block of equations (period t+1)
% with respect to the endogenous variables.
if purely_backward_model || purely_forward_model
    for i=DynareModel.eq_nbr+1:2*DynareModel.eq_nbr
        for j=1:DynareModel.endo_nbr
            if I0(j)
                fprintf(fid, '    J(%u,%u) = J(%u,%u)+J(%u,%u)*y(%u);\n', i, n0+n1+n2+j, i, n0+n1+n2+j, i, I0(j), j);
            end
            if I1(j)
                fprintf(fid, '    J(%u,%u) = J(%u,%u)+2*J(%u,%u)*y(%u)*g(%u);\n', i, n0+n1+n2+j, i, n0+n1+n2+j, i, I1(j), j, j);
            end
        end
    end
else
    for i=DynareModel.eq_nbr+1:2*DynareModel.eq_nbr
        for j=1:DynareModel.endo_nbr
            if I2(j)
                if I1(j)
                    if I0(j)
                        fprintf(fid, '    J(%u,%u) = J(%u,%u)*y(%u)+2*J(%u,%u)*y(%u)*g(%u)+3*J(%u,%u)*y(%u)*g(%u)*g(%u);\n', i, n0+n1+n2+j, i, I0(j), j, i, I1(j), j, j, i, I2(j), j, j, j);
                    else
                        fprintf(fid, '    J(%u,%u) = 2*J(%u,%u)*y(%u)*g(%u)+3*J(%u,%u)*y(%u)*g(%u)*g(%u);\n', i, n0+n1+n2+j, i, I1(j), j, j, i, I2(j), j, j, j);
                    end
                else
                    if I0(j)
                        fprintf(fid, '    J(%u,%u) = J(%u,%u)*y(%u)+3*J(%u,%u)*y(%u)*g(%u)*g(%u);\n', i, n0+n1+n2+j, i, I0(j), j, i, I2(j), j, j, j);
                    else
                        fprintf(fid, '    J(%u,%u) = 3*J(%u,%u)*y(%u)*g(%u)*g(%u);\n', i, n0+n1+n2+j, i, I2(j), j, j, j);
                    end
                end
            else
                if I1(j)
                    if I0(j)
                        fprintf(fid, '    J(%u,%u) = J(%u,%u)*y(%u)+2*J(%u,%u)*y(%u)*g(%u);\n', i, n0+n1+n2+j, i, I0(j), j, i, I1(j), j, j);
                    else
                        fprintf(fid, '    J(%u,%u) = 2*J(%u,%u)*y(%u)*g(%u);\n', i, n0+n1+n2+j, i, I1(j), j, j);
                    end
                else
                    if I0(j)
                        fprintf(fid, '    J(%u,%u) = J(%u,%u)*y(%u);\n', i, n0+n1+n2+j, i, I0(j), j);
                    end
                end
            end
        end
    end
end

fprintf(fid, '    JAC(:,%u:%u) = J(:,%u:%u);\n', DynareModel.endo_nbr+1, 2*DynareModel.endo_nbr, n0+n1+n2+1, n0+n1+n2+DynareModel.endo_nbr);

fprintf(fid,'end\n');
fclose(fid);