Move the location of various generated files on the filesystem

- M and MEX files are now under +${MODELNAME}/
- bytecode, C source and JSON now under ${MODELNAME}/model/
time-shift
Sébastien Villemot 2018-06-27 17:02:13 +02:00
parent f5bf76deb5
commit a1b8bd39b2
58 changed files with 148 additions and 166 deletions

View File

@ -772,30 +772,21 @@ file. By default (unless @code{use_dll} option has been given to
@table @file
@item @var{FILENAME}.m
@item +@var{FILENAME}/driver.m
Contains variable declarations, and computing tasks
@item @var{FILENAME}_dynamic.m
@item +@var{FILENAME}/dynamic.m
@vindex M_.lead_lag_incidence
Contains the dynamic model equations. Note that Dynare might introduce auxiliary equations and variables (@pxref{Auxiliary variables}). Outputs are the residuals of the dynamic model equations in the order the equations were declared and the Jacobian of the dynamic model equations. For higher order approximations also the Hessian and the third-order derivatives are provided. When computing the Jacobian of the dynamic model, the order of the endogenous variables in the columns is stored in @code{M_.lead_lag_incidence}. The rows of this matrix represent time periods: the first row denotes a lagged (time t-1) variable, the second row a contemporaneous (time t) variable, and the third row a leaded (time t+1) variable. The columns of the matrix represent the endogenous variables in their order of declaration. A zero in the matrix means that this endogenous does not appear in the model in this time period. The value in the @code{M_.lead_lag_incidence} matrix corresponds to the column of that variable in the Jacobian of the dynamic model. Example: Let the second declared variable be @code{c} and the @code{(3,2)} entry of @code{M_.lead_lag_incidence} be @code{15}. Then the @code{15}th column of the Jacobian is the derivative with respect to @code{c(+1)}.
@item @var{FILENAME}_static.m
@item +@var{FILENAME}/static.m
Contains the long run static model equations. Note that Dynare might introduce auxiliary equations and variables (@pxref{Auxiliary variables}). Outputs are the residuals of the static model equations in the order the equations were declared and the Jacobian of the static equations. Entry @code{(i,j)} of the Jacobian represents the derivative of the @code{i}th static model equation with respect to the @code{j}th model variable in declaration order.
@end table
@noindent
These files may be looked at to understand errors reported at the simulation stage.
@code{dynare} will then run the computing tasks by executing @file{@var{FILENAME}.m}.
A few words of warning is warranted here: the filename of the
@file{.mod} file should be chosen in such a way that the generated
@file{.m} files described above do not conflict with @file{.m} files
provided by MATLAB/Octave or by Dynare. Not respecting this rule could
cause crashes or unexpected behaviour. In particular, it means that
the @file{.mod} file cannot be given the name of a MATLAB/Octave or
Dynare command. Under Octave, it also means that the @file{.mod} file
cannot be named @file{test.mod}.
@code{dynare} will then run the computing tasks by executing @file{+@var{FILENAME}/driver.m}.
@optionshead
@ -869,16 +860,16 @@ Suppresses all warnings.
@item json=parse|transform|compute
Causes the preprocessor to output a version of the @file{.mod} file in JSON
format. If @code{parse} is passed, the output will be written after the parsing
of the @file{.mod} file to a file called @file{@var{FILENAME}.json}. If
of the @file{.mod} file to a file called @file{@var{FILENAME}/model/json/modfile.json}. If
@code{transform} is passed, the JSON output of the transformed model (maximum
lead of 1, minimum lag of -1, expectation operators substituted, etc.) will be
written to a file called @file{@var{FILENAME}.json} and the original,
untransformed model will be written in @file{@var{FILENAME}_original.json}. And
written to a file called @file{@var{FILENAME}/model/json/modfile.json} and the original,
untransformed model will be written in @file{@var{FILENAME}/model/json/modfile-original.json}. And
if @code{compute} is passed, the output is written after the computing pass. In
this case, the transformed model is written to @file{@var{FILENAME}.json}, the
original model is written to @file{@var{FILENAME}_original.json}, and the
dynamic and static files are written to @file{@var{FILENAME}_dynamic.json} and
@file{@var{FILENAME}_static.json}.
this case, the transformed model is written to @file{@var{FILENAME}/model/json/modfile.json}, the
original model is written to @file{@var{FILENAME}/model/json/model-original.json}, and the
dynamic and static files are written to @file{@var{FILENAME}/model/json/dynamic.json} and
@file{@var{FILENAME}/model/json/static.json}.
@item jsonstdout
Instead of writing output requested by @ref{json} to files, write to standard
@ -889,8 +880,8 @@ Quit processing once the output requested by @ref{json} has been written.
@item jsonderivsimple
Print a simplified version (excluding variable name(s) and lag information) of the
static and dynamic files in @file{@var{FILENAME}_static.json} and
@file{@var{FILENAME}_dynamic.json}.
static and dynamic files in @file{@var{FILENAME}/model/json/static.json} and
@file{@var{FILENAME}/model/json/dynamic.json}.
@item warn_uninit
Display a warning for each variable or parameter which is not
@ -3549,7 +3540,7 @@ is described below in more details. See also @file{fs2000.mod} in the
@file{examples} directory for an example.
The steady state file generated by Dynare will be called
@file{@var{FILENAME}_steadystate2.m}.
@file{+@var{FILENAME}/steadystate.m}.
@item
You can write the corresponding MATLAB function by hand. If your

View File

@ -233,16 +233,16 @@ else
'dataset_info',dataset_info);
% which files have to be copied to run remotely
NamFileInput(1,:) = {'',[M_.fname '_static.m']};
NamFileInput(2,:) = {'',[M_.fname '_dynamic.m']};
NamFileInput(1,:) = {'',[M_.fname '.static.m']};
NamFileInput(2,:) = {'',[M_.fname '.dynamic.m']};
if M_.set_auxiliary_variables
NamFileInput(3,:) = {'',[M_.fname '_set_auxiliary_variables.m']};
NamFileInput(3,:) = {'',[M_.fname '.set_auxiliary_variables.m']};
end
if options_.steadystate_flag
if options_.steadystate_flag == 1
NamFileInput(length(NamFileInput)+1,:)={'',[M_.fname '_steadystate.m']};
else
NamFileInput(length(NamFileInput)+1,:)={'',[M_.fname '_steadystate2.m']};
NamFileInput(length(NamFileInput)+1,:)={'',[M_.fname '.steadystate.m']};
end
end
[fout] = masterParallel(options_.parallel, 1, B,NamFileInput,'PosteriorIRF_core1', localVars, globalVars, options_.parallel_info);

View File

@ -82,7 +82,7 @@ ModelInversion.x_free_id = freeinnovations_id;
ModelInversion.J_id = [ModelInversion.y_free_id ; sum(DynareModel.lead_lag_incidence(:)>0)+ModelInversion.x_free_id];
% Get the name of the dynamic model routines.
model_dynamic = str2func([DynareModel.fname,'_dynamic']);
model_dynamic = str2func([DynareModel.fname,'.dynamic']);
model_dtransf = str2func('dynamic_backward_model_for_inversion');
% Initialization of the returned simulations (endogenous variables).

View File

@ -36,7 +36,7 @@ function [residuals, info] = calibrateresiduals(dbase, info, DynareModel)
% along with Dynare. If not, see <http://www.gnu.org/licenses/>.
% Get function handle for the dynamic model
model_dynamic = str2func([DynareModel.fname,'_dynamic']);
model_dynamic = str2func([DynareModel.fname,'.dynamic']);
% Get data for all the endogenous variables.
ydata = dbase{info.endonames{:}}.data;

View File

@ -31,7 +31,7 @@ if nargin<3
inversionflag = false;
end
set_auxiliary_series = [DynareModel.fname '_dynamic_set_auxiliary_series'];
set_auxiliary_series = [DynareModel.fname '.dynamic_set_auxiliary_series'];
if exist([set_auxiliary_series '.m'])
dbase = feval(set_auxiliary_series, dbase, DynareModel.params);

View File

@ -203,7 +203,7 @@ if nargout>8
idx = 1:DynareModel.endo_nbr;
jdx = idx+ny1;
% Get the name of the dynamic model routine.
model_dynamic = str2func([DynareModel.fname,'_dynamic']);
model_dynamic = str2func([DynareModel.fname,'.dynamic']);
% initialization of vector y.
y = NaN(length(idx)+ny1,1);
end

View File

@ -21,5 +21,5 @@ function [r, g1] = block_mfs_steadystate(y, b, y_all, exo, params, M)
y_all(M.block_structure_stat.block(b).variable) = y;
eval(['[r,g1] = ' M.fname '_static(b, y_all, exo, params);']);
eval(['[r,g1] = ' M.fname '.static(b, y_all, exo, params);']);
g1 = full(g1);

View File

@ -77,7 +77,7 @@ options_.threads.local_state_space_iteration_2 = 1;
options_.jacobian_flag = 1;
% steady state file
if exist([M_.fname '_steadystate2.m'],'file')
if exist(['+' M_.fname '/steadystate.m'],'file')
options_.steadystate_flag = 2;
elseif exist([M_.fname '_steadystate.m'],'file')
options_.steadystate_flag = 1;

View File

@ -66,7 +66,7 @@ if options_.steadystate_flag
[junk,M_.params,info] = evaluate_steady_state_file(oo_.steady_state,[oo_.exo_steady_state; oo_.exo_det_steady_state],M_, ...
options_,0);
end
[U,Uy,W] = feval([M_.fname,'_objective_static'],zeros(endo_nbr,1),[], M_.params);
[U,Uy,W] = feval([M_.fname,'.objective.static'],zeros(endo_nbr,1),[], M_.params);
if any(any(Uy~=0))
non_zero_derivs=find(any(Uy~=0));
for ii=1:length(non_zero_derivs)
@ -94,7 +94,7 @@ if exo_nbr == 0
oo_.exo_steady_state = [] ;
end
[junk,jacobia_] = feval([M_.fname '_dynamic'],z, [zeros(size(oo_.exo_simul)) ...
[junk,jacobia_] = feval([M_.fname '.dynamic'],z, [zeros(size(oo_.exo_simul)) ...
oo_.exo_det_simul], M_.params, zeros(endo_nbr,1), it_);
if any(junk~=0)
error(['discretionary_policy: the model must be written in deviation ' ...

View File

@ -73,7 +73,7 @@ end
if (options_.bytecode)
[chck, zz, data]= bytecode('dynamic','evaluate', z, zx, M_.params, dr.ys, 1, data);
else
[r, data] = feval([M_.fname '_dynamic'], options_, M_, oo_, z', zx, M_.params, dr.ys, M_.maximum_lag+1, data);
[r, data] = feval([M_.fname '.dynamic'], options_, M_, oo_, z', zx, M_.params, dr.ys, M_.maximum_lag+1, data);
chck = 0;
end
mexErrCheck('bytecode', chck);

View File

@ -135,7 +135,7 @@ xx(1:M.orig_endo_nbr) = x(1:M.orig_endo_nbr); %set values of original endogenous
if any([M.aux_vars.type] ~= 6) %auxiliary variables other than multipliers
needs_set_auxiliary_variables = 1;
if M.set_auxiliary_variables
fh = str2func([M.fname '_set_auxiliary_variables']);
fh = str2func([M.fname '.set_auxiliary_variables']);
s_a_v_func = @(z) fh(z,...
[oo.exo_steady_state,...
oo.exo_det_steady_state],...
@ -150,7 +150,7 @@ end
% value and Jacobian of objective function
ex = zeros(1,M.exo_nbr);
[U,Uy,Uyy] = feval([fname '_objective_static'],x,ex, params);
[U,Uy,Uyy] = feval([fname '.objective.static'],x,ex, params);
Uyy = reshape(Uyy,endo_nbr,endo_nbr);
% set multipliers and auxiliary variables that
@ -160,7 +160,7 @@ if (options_.bytecode)
params, 'evaluate');
fJ = junk.g1;
else
[res,fJ] = feval([fname '_static'],xx,[oo.exo_steady_state oo.exo_det_steady_state], ...
[res,fJ] = feval([fname '.static'],xx,[oo.exo_steady_state oo.exo_det_steady_state], ...
params);
end
% index of multipliers and corresponding equations
@ -197,7 +197,7 @@ if (options_.bytecode)
[chck, res, junk] = bytecode('static',ys,[oo.exo_steady_state oo.exo_det_steady_state], ...
M.params, 'evaluate');
else
res = feval([M.fname '_static'],ys,[oo.exo_steady_state oo.exo_det_steady_state], ...
res = feval([M.fname '.static'],ys,[oo.exo_steady_state oo.exo_det_steady_state], ...
M.params);
end
if norm(res) < options_.solve_tolf

View File

@ -160,7 +160,7 @@ end
z = repmat(ys,1,3);
z = z(iyr0) ;
[resid1,d1,d2] = feval([M.fname '_dynamic'],z,...
[resid1,d1,d2] = feval([M.fname '.dynamic'],z,...
[oo.exo_simul ...
oo.exo_det_simul], M.params, dr.ys, 2);
if ~isreal(d1) || ~isreal(d2)
@ -204,7 +204,7 @@ v_np = pm.v_np;
% assumed portfolio
dr.ys(v_p) = x;
ys0 = dr.ys(v_np);
f_h =str2func([M.fname '_static']);
f_h =str2func([M.fname '.static']);
[dr.ys(v_np),info] = csolve(@ds_static_model,ys0,[],1e-10,100,f_h,x,pm.eq_np,v_np,v_p, ...
M.endo_nbr,M.exo_nbr,M.params);
if info
@ -220,7 +220,7 @@ iyr0 = find(iyv) ;
z = repmat(dr.ys,1,3);
z = z(iyr0) ;
[resid1,d1,d2] = feval([M.fname '_dynamic'],z,...
[resid1,d1,d2] = feval([M.fname '.dynamic'],z,...
[oo.exo_simul ...
oo.exo_det_simul], M.params, dr.ys, 2);
if ~isreal(d1) || ~isreal(d2)
@ -253,7 +253,7 @@ ys(pm.v_p) = x;
z = repmat(ys,1,3);
z = z(iyr0) ;
[resid1,d1,d2] = feval([M.fname '_dynamic'],z,...
[resid1,d1,d2] = feval([M.fname '.dynamic'],z,...
[oo.exo_simul ...
oo.exo_det_simul], M.params, dr.ys, 2);
if ~isreal(d1) || ~isreal(d2)
@ -294,7 +294,7 @@ end
z = repmat(ys,1,3);
z = z(iyr0) ;
[resid1,d1,d2] = feval([M.fname '_dynamic'],z,...
[resid1,d1,d2] = feval([M.fname '.dynamic'],z,...
[oo.exo_simul ...
oo.exo_det_simul], M.params, dr.ys, 2);
@ -329,7 +329,7 @@ resid = resid1+0.5*(d1(:,i_fwrd_f1)*ghuu(i_fwrd_g,:)+d2(:,i_fwrd_f2)* ...
kron(gu1,gu1))*vec(M.Sigma_e);
if nargout > 1
[resid1,d1,d2,d3] = feval([M.fname '_dynamic'],z,...
[resid1,d1,d2,d3] = feval([M.fname '.dynamic'],z,...
[oo.exo_simul ...
oo.exo_det_simul], M.params, dr.ys, 2);

View File

@ -137,7 +137,7 @@ else
fnamelength = length(fname) - 4;
end
if fnamelength + length('_set_auxiliary_variables') > namelengthmax()
if fnamelength + length('.set_auxiliary_variables') > namelengthmax()
error('The name of your MOD file is too long, please shorten it')
end
@ -283,4 +283,4 @@ end
if ~ isempty(find(abs(fname) == 46))
fname = fname(:,1:find(abs(fname) == 46)-1) ;
end
evalin('base',fname) ;
evalin('base',[fname '.driver']) ;

View File

@ -36,7 +36,7 @@ if options.block && ~options.bytecode
ss(M.block_structure_stat.block(b).variable) = y;
else
n = length(M.block_structure_stat.block(b).variable);
[ss, check] = solve_one_boundary([M.fname '_static_' int2str(b)], ss, exo, ...
[ss, check] = solve_one_boundary([M.fname '.block.static_' int2str(b)], ss, exo, ...
params, [], M.block_structure_stat.block(b).variable, n, 1, 0, b, 0, options.simul.maxit, ...
options.solve_tolf, ...
options.slowc, 0, options.solve_algo, 1, 0, 0,M,options);
@ -46,7 +46,7 @@ if options.block && ~options.bytecode
end
end
end
[r, g1, x] = feval([M.fname '_static'], b, ss, ...
[r, g1, x] = feval([M.fname '.static'], b, ss, ...
exo, params);
end
elseif options.bytecode

View File

@ -17,7 +17,7 @@ function e = euler_equation_error(y0,x,innovations,M,options,oo,pfm,nodes,weight
% You should have received a copy of the GNU General Public License
% along with Dynare. If not, see <http://www.gnu.org/licenses/>.
dynamic_model = str2func([M.fname '_dynamic']);
dynamic_model = str2func([M.fname '.dynamic']);
ep = options.ep;
[y1, info_convergence, endogenousvariablespaths] = extended_path_core(ep.periods, ...
M.endo_nbr, M.exo_nbr, ...

View File

@ -66,7 +66,7 @@ else
end
pfm.i_cols_j = 1:pfm.nd;
pfm.i_upd = pfm.ny+(1:pfm.periods*pfm.ny);
pfm.dynamic_model = str2func([DynareModel.fname,'_dynamic']);
pfm.dynamic_model = str2func([DynareModel.fname,'.dynamic']);
pfm.verbose = DynareOptions.ep.verbosity;
pfm.maxit_ = DynareOptions.simul.maxit;
pfm.tolerance = DynareOptions.dynatol.f;

View File

@ -49,7 +49,7 @@ gu(dr.order_var,:) = dr.ghu;
ys = oo.dr.ys;
[U,Uy,Uyy] = feval([M.fname '_objective_static'],ys,zeros(1,exo_nbr), ...
[U,Uy,Uyy] = feval([M.fname '.objective.static'],ys,zeros(1,exo_nbr), ...
M.params);
%second order terms
Uyy = full(Uyy);

View File

@ -43,7 +43,7 @@ if options.bytecode
exo_ss, params, ys, 1);
mexErrCheck('bytecode', check1);
else
fh_static = str2func([M.fname '_static']);
fh_static = str2func([M.fname '.static']);
if options.block
residuals = zeros(M.endo_nbr,1);
for b = 1:length(M.block_structure_stat.block)

View File

@ -47,7 +47,7 @@ params = M.params;
exo_ss = [oo.exo_steady_state; oo.exo_det_steady_state];
if length(M.aux_vars) > 0 && ~steadystate_flag
h_set_auxiliary_variables = str2func([M.fname '_set_auxiliary_variables']);
h_set_auxiliary_variables = str2func([M.fname '.set_auxiliary_variables']);
if M.set_auxiliary_variables
ys_init = h_set_auxiliary_variables(ys_init,exo_ss,M.params);
end
@ -219,7 +219,7 @@ elseif steadystate_flag
elseif (options.bytecode == 0 && options.block == 0)
if options.linear == 0
% non linear model
static_model = str2func([M.fname '_static']);
static_model = str2func([M.fname '.static']);
[ys,check] = dynare_solve(@static_problem,...
ys_init,...
options, exo_ss, params,...
@ -247,7 +247,7 @@ elseif (options.bytecode == 0 && options.block == 0)
end
else
% linear model
fh_static = str2func([M.fname '_static']);
fh_static = str2func([M.fname '.static']);
[fvec,jacob] = fh_static(ys_init,exo_ss, ...
params);
@ -312,12 +312,12 @@ if M.static_and_dynamic_models_differ
[chck, r, junk]= bytecode('dynamic','evaluate', z, zx, M.params, ys, 1);
mexErrCheck('bytecode', chck);
elseif options.block
[r, oo.dr] = feval([M.fname '_dynamic'], z', zx, M.params, ys, M.maximum_lag+1, oo.dr);
[r, oo.dr] = feval([M.fname '.dynamic'], z', zx, M.params, ys, M.maximum_lag+1, oo.dr);
else
iyv = M.lead_lag_incidence';
iyr0 = find(iyv(:));
xys = z(iyr0);
r = feval([M.fname '_dynamic'], z(iyr0), zx, M.params, ys, M.maximum_lag + 1);
r = feval([M.fname '.dynamic'], z(iyr0), zx, M.params, ys, M.maximum_lag + 1);
end
% Fail if residual greater than tolerance
if max(abs(r)) > options.solve_tolf

View File

@ -53,7 +53,7 @@ if options.steadystate_flag == 1
params1 = evalin('base','M_.params');
else % steadystate_flag == 2
% new format
h_steadystate = str2func([fname '_steadystate2']);
h_steadystate = str2func([fname '.steadystate']);
[ys,params1,check] = h_steadystate(ys_init, exo_ss, params);
end
@ -75,7 +75,7 @@ if M.set_auxiliary_variables
% variables only if the model has auxiliary variables. Otherwise
% Octave may crash (see https://savannah.gnu.org/bugs/?52568) because
% the function does not exist in the path.
h_set_auxiliary_variables = str2func([M.fname '_set_auxiliary_variables']);
h_set_auxiliary_variables = str2func([M.fname '.set_auxiliary_variables']);
end
if isnan(updated_params_flag) || (updated_params_flag && any(isnan(params(~isnan(params))-params1(~isnan(params))))) %checks if new NaNs were added

View File

@ -33,5 +33,5 @@ global it_ M_ oo_
n1 = size(x,1) - M_.exo_nbr;
oo_.exo_simul(it_+M_.maximum_lag-M_.maximum_lag,:) = x(n1+1:end)';
fh = str2func([M_.fname '_static']);
fh = str2func([M_.fname '.static']);
y=feval(fh,x(1:n1),oo_.exo_simul, M_.params);

View File

@ -124,7 +124,7 @@ end
if kronflag==-2
if nargout>5
[residual, g1, g2 ] = feval([M_.fname,'_dynamic'],yy0, oo_.exo_steady_state', ...
[residual, g1, g2 ] = feval([M_.fname,'.dynamic'],yy0, oo_.exo_steady_state', ...
M_.params, oo_.dr.ys, 1);
g22 = hessian_sparse('thet2tau',[M_.params(indx)],options_.gstep,estim_params_,M_, oo_, indx,[],-1);
H2ss=full(g22(1:M_.endo_nbr,:));
@ -144,7 +144,7 @@ if kronflag==-2
g22 = gx22;
clear gx22;
else
[residual, g1 ] = feval([M_.fname,'_dynamic'],yy0, oo_.exo_steady_state', ...
[residual, g1 ] = feval([M_.fname,'.dynamic'],yy0, oo_.exo_steady_state', ...
M_.params, oo_.dr.ys, 1);
end
gp = fjaco('thet2tau',[M_.params(indx)],estim_params_,M_, oo_, indx,[],-1);
@ -154,17 +154,17 @@ if kronflag==-2
else
dyssdtheta=zeros(length(oo_.dr.ys),M_.param_nbr);
d2yssdtheta=zeros(length(oo_.dr.ys),M_.param_nbr,M_.param_nbr);
[residual, gg1] = feval([M_.fname,'_static'],oo_.dr.ys, oo_.exo_steady_state', M_.params);
df = feval([M_.fname,'_static_params_derivs'],oo_.dr.ys, repmat(oo_.exo_steady_state',[M_.maximum_exo_lag+M_.maximum_exo_lead+1]), ...
[residual, gg1] = feval([M_.fname,'.static'],oo_.dr.ys, oo_.exo_steady_state', M_.params);
df = feval([M_.fname,'.static_params_derivs'],oo_.dr.ys, repmat(oo_.exo_steady_state',[M_.maximum_exo_lag+M_.maximum_exo_lead+1]), ...
M_.params);
dyssdtheta = -gg1\df;
if nargout>5
[residual, gg1, gg2] = feval([M_.fname,'_static'],oo_.dr.ys, oo_.exo_steady_state', M_.params);
[residual, g1, g2, g3] = feval([M_.fname,'_dynamic'],yy0, oo_.exo_steady_state', ...
[residual, gg1, gg2] = feval([M_.fname,'.static'],oo_.dr.ys, oo_.exo_steady_state', M_.params);
[residual, g1, g2, g3] = feval([M_.fname,'.dynamic'],yy0, oo_.exo_steady_state', ...
M_.params, oo_.dr.ys, 1);
[nr, nc]=size(gg2);
[df, gpx, d2f] = feval([M_.fname,'_static_params_derivs'],oo_.dr.ys, oo_.exo_steady_state', ...
[df, gpx, d2f] = feval([M_.fname,'.static_params_derivs'],oo_.dr.ys, oo_.exo_steady_state', ...
M_.params);%, oo_.dr.ys, 1, dyssdtheta*0, d2yssdtheta);
d2f = get_all_resid_2nd_derivs(d2f,length(oo_.dr.ys),M_.param_nbr);
if isempty(find(gg2))
@ -207,7 +207,7 @@ else
else
[df, gp] = feval([M_.fname,'_params_derivs'],yy0, repmat(oo_.exo_steady_state',[M_.maximum_exo_lag+M_.maximum_exo_lead+1,1]), ...
M_.params, oo_.dr.ys, 1, dyssdtheta,d2yssdtheta);
[residual, g1, g2 ] = feval([M_.fname,'_dynamic'],yy0, repmat(oo_.exo_steady_state',[M_.maximum_exo_lag+M_.maximum_exo_lead+1,1]), ...
[residual, g1, g2 ] = feval([M_.fname,'.dynamic'],yy0, repmat(oo_.exo_steady_state',[M_.maximum_exo_lag+M_.maximum_exo_lead+1,1]), ...
M_.params, oo_.dr.ys, 1);
end

View File

@ -77,7 +77,7 @@ if info(1)==0
oo0=oo_;
tau=[oo_.dr.ys(oo_.dr.order_var); vec(A); dyn_vech(B*M_.Sigma_e*B')];
yy0=oo_.dr.ys(I);
[residual, g1 ] = feval([M_.fname,'_dynamic'],yy0, ...
[residual, g1 ] = feval([M_.fname,'.dynamic'],yy0, ...
repmat(oo_.exo_steady_state',[M_.maximum_exo_lag+M_.maximum_exo_lead+1]), M_.params, ...
oo_.dr.ys, 1);
vg1 = [oo_.dr.ys(oo_.dr.order_var); vec(g1)];

View File

@ -54,7 +54,7 @@ i_upd = ny+(1:periods*ny);
x = endo_simul(:);
model_dynamic = str2func([M.fname,'_dynamic']);
model_dynamic = str2func([M.fname,'.dynamic']);
z = x(find(lead_lag_incidence'));
[res,A] = model_dynamic(z, exo_simul, params, steady_state,2);
nnzA = nnz(A);

View File

@ -115,7 +115,7 @@ for b=1:nb
int2str(b)]);
end
else
[res,jacob]=feval([M.fname '_static'],dr.ys,exo,M.params);
[res,jacob]=feval([M.fname '.static'],dr.ys,exo,M.params);
end
if any(any(isinf(jacob) | isnan(jacob)))
problem_dummy=1;
@ -208,7 +208,7 @@ if ~options.block
M.params, dr.ys, 1);
jacobia_ = [loc_dr.g1 loc_dr.g1_x loc_dr.g1_xd];
else
[junk,jacobia_] = feval([M.fname '_dynamic'],z(iyr0),exo_simul, ...
[junk,jacobia_] = feval([M.fname '.dynamic'],z(iyr0),exo_simul, ...
M.params, dr.ys, it_);
end
elseif options.order >= 2
@ -217,7 +217,7 @@ if ~options.block
M.params, dr.ys, 1);
jacobia_ = [loc_dr.g1 loc_dr.g1_x];
else
[junk,jacobia_,hessian1] = feval([M.fname '_dynamic'],z(iyr0),...
[junk,jacobia_,hessian1] = feval([M.fname '.dynamic'],z(iyr0),...
exo_simul, ...
M.params, dr.ys, it_);
end

View File

@ -5,15 +5,15 @@ y = z(find(ll(:)));
switch nargout
case 1
r = feval([M.fname '_dynamic'],y,x, ...
r = feval([M.fname '.dynamic'],y,x, ...
M.params, ss, 1);
case 2
[r,g1] = feval([M.fname '_dynamic'],y,x, ...
[r,g1] = feval([M.fname '.dynamic'],y,x, ...
M.params, ss, 1);
case 3
[r,g1,g2] = feval([M.fname '_dynamic'],y,x, ...
[r,g1,g2] = feval([M.fname '.dynamic'],y,x, ...
M.params, ss, 1);
case 4
[r,g1,g2,g3] = feval([M.fname '_dynamic'],y,x, ...
[r,g1,g2,g3] = feval([M.fname '.dynamic'],y,x, ...
M.params, ss, 1);
end

View File

@ -44,11 +44,11 @@ y = [y;ys_(lea_cols)];
if ismac
eval(['[resid,g1]=',M_.fname,'_dynamic(y,x, M_.params, ys_, it_);']);
eval(['[resid,g1]=',M_.fname,'.dynamic(y,x, M_.params, ys_, it_);']);
% Older versions of DYNARE for Mac did not include ys_ in the call structure
%eval(['[resid,g1]=',M_.fname,'_dynamic(y,x, M_.params, it_);']);
%eval(['[resid,g1]=',M_.fname,'.dynamic(y,x, M_.params, it_);']);
else
eval(['[resid,g1]=',M_.fname,'_dynamic(y,x, M_.params, ys_, it_);']);
eval(['[resid,g1]=',M_.fname,'.dynamic(y,x, M_.params, ys_, it_);']);
end

View File

@ -65,9 +65,9 @@ end
dr.ys = ys;
check1 = 0;
% testing for steadystate file
fh = str2func([M_.fname '_static']);
fh = str2func([M_.fname '.static']);
if options_.steadystate_flag
[dr.ys,check1] = feval([M_.fname '_steadystate'],dr.ys,...
[dr.ys,check1] = feval([M_.fname '.steadystate'],dr.ys,...
[oo_.exo_steady_state; ...
oo_.exo_det_steady_state]);
if size(dr.ys,1) < M_.endo_nbr
@ -78,7 +78,7 @@ if options_.steadystate_flag
oo_.exo_det_steady_state,...
M_.params);
else
error([M_.fname '_steadystate.m doesn''t match the model']);
error([M_.fname '.steadystate doesn''t match the model']);
end
end

View File

@ -27,7 +27,7 @@ for i=1:n+1
[exo_steady_state; ...
exo_det_steady_state],params);
else
res = feval([fname '_static'],ys1,...
res = feval([fname '.static'],ys1,...
[exo_steady_state; ...
exo_det_steady_state],params);
end

View File

@ -130,7 +130,7 @@ else
z = repmat(dr.ys,1,klen);
end
z = z(iyr0) ;
[junk,jacobia_] = feval([M_.fname '_dynamic'],z,[oo_.exo_simul ...
[junk,jacobia_] = feval([M_.fname '.dynamic'],z,[oo_.exo_simul ...
oo_.exo_det_simul], M_.params, dr.ys, it_);
if options_.ACES_solver==1 && (length(sim_ruleids)>0 || length(tct_ruleids)>0 )

View File

@ -475,7 +475,7 @@ if pf && ~surprise
if (options_.bytecode)
[chck, zz, data1]= bytecode('dynamic','evaluate', z, zx, M_.params, oo_.steady_state, k, data1);
else
[zz, g1b] = feval([M_.fname '_dynamic'], z', zx, M_.params, oo_.steady_state, k);
[zz, g1b] = feval([M_.fname '.dynamic'], z', zx, M_.params, oo_.steady_state, k);
data1.g1_x = g1b(:,end - M_.exo_nbr + 1:end);
data1.g1 = g1b(:,1 : end - M_.exo_nbr);
chck = 0;
@ -747,7 +747,7 @@ else
if (options_.bytecode)
[chck, zz, data1]= bytecode('dynamic','evaluate', z, zx, M_.params, oo_.steady_state, k, data1);
else
[zz, g1b] = feval([M_.fname '_dynamic'], z', zx, M_.params, oo_.steady_state, k);
[zz, g1b] = feval([M_.fname '.dynamic'], z', zx, M_.params, oo_.steady_state, k);
data1.g1_x = g1b(:,end - M_.exo_nbr + 1:end);
data1.g1 = g1b(:,1 : end - M_.exo_nbr);
chck = 0;

View File

@ -39,7 +39,7 @@ params = M_.params;
endo_simul = oo_.endo_simul;
exo_simul = oo_.exo_simul;
model_dynamic = str2func([M_.fname,'_dynamic']);
model_dynamic = str2func([M_.fname,'.dynamic']);
residuals = zeros(ny,periods);

View File

@ -41,7 +41,7 @@ persistent lead_lag_incidence dynamic_model ny nyp nyf nrs nrc iyf iyp isp is is
if ~nargin && isempty(iflag)% Initialization of the persistent variables.
lead_lag_incidence = M_.lead_lag_incidence;
dynamic_model = [M_.fname '_dynamic'];
dynamic_model = [M_.fname '.dynamic'];
ny = size(oo_.endo_simul,1);
nyp = nnz(lead_lag_incidence(1,:));% number of lagged variables.
nyf = nnz(lead_lag_incidence(3,:));% number of leaded variables.

View File

@ -40,7 +40,7 @@ end
options_.scalv= 1;
if options_.debug
model_static = str2func([M_.fname,'_static']);
model_static = str2func([M_.fname,'.static']);
for ii=1:size(oo_.exo_simul,1)
[residual(:,ii)] = model_static(oo_.steady_state, oo_.exo_simul(ii,:),M_.params);
end
@ -183,7 +183,7 @@ if ~isreal(oo_.endo_simul(:)) %can only happen without bytecode
illi = illi(:,2:3);
[i_cols_J1,junk,i_cols_1] = find(illi(:));
i_cols_T = nonzeros(M_.lead_lag_incidence(1:2,:)');
residuals = perfect_foresight_problem(yy(:),str2func([M_.fname '_dynamic']), y0, yT, ...
residuals = perfect_foresight_problem(yy(:),str2func([M_.fname '.dynamic']), y0, yT, ...
oo_.exo_simul,M_.params,oo_.steady_state, ...
M_.maximum_lag,options_.periods,M_.endo_nbr,i_cols, ...
i_cols_J1, i_cols_1, i_cols_T, i_cols_j, ...

View File

@ -58,7 +58,7 @@ if options_.block
mexErrCheck('bytecode', info);
end
else
oo_ = feval([M_.fname '_dynamic'], options_, M_, oo_);
oo_ = feval([M_.fname '.dynamic'], options_, M_, oo_);
end
else
if options_.bytecode
@ -134,7 +134,7 @@ if nargout>1
if options_.bytecode
[chck, residuals, junk]= bytecode('dynamic','evaluate', oo_.endo_simul, oo_.exo_simul, M_.params, oo_.steady_state, 1);
else
residuals = perfect_foresight_problem(yy(:),str2func([M_.fname '_dynamic']), y0, yT, ...
residuals = perfect_foresight_problem(yy(:),str2func([M_.fname '.dynamic']), y0, yT, ...
oo_.exo_simul,M_.params,oo_.steady_state, ...
M_.maximum_lag,options_.periods,M_.endo_nbr,i_cols, ...
i_cols_J1, i_cols_1, i_cols_T, i_cols_j, ...

View File

@ -75,4 +75,4 @@ illi = M.lead_lag_incidence';
illi = illi(:,2:3);
[i_cols_J1, junk,i_cols_1] = find(illi(:));
i_cols_T = nonzeros(M.lead_lag_incidence(1:2,:)');
dynamicmodel = str2func([M.fname,'_dynamic']);
dynamicmodel = str2func([M.fname,'.dynamic']);

View File

@ -43,7 +43,7 @@ if options_.block
mexErrCheck('bytecode', info);
end
else
oo_ = feval([M_.fname '_dynamic'], options_, M_, oo_);
oo_ = feval([M_.fname '.dynamic'], options_, M_, oo_);
end
else
if options_.bytecode
@ -110,7 +110,7 @@ if nargout>1
if options_.bytecode
[chck, residuals, junk]= bytecode('dynamic','evaluate', oo_.endo_simul, oo_.exo_simul, M_.params, oo_.steady_state, 1);
else
residuals = perfect_foresight_problem(yy(:),str2func([M_.fname '_dynamic']), y0, yT, ...
residuals = perfect_foresight_problem(yy(:),str2func([M_.fname '.dynamic']), y0, yT, ...
oo_.exo_simul,M_.params,oo_.steady_state, ...
M_.maximum_lag,options_.periods,M_.endo_nbr,i_cols, ...
i_cols_J1, i_cols_1, i_cols_T, i_cols_j, ...

View File

@ -76,7 +76,7 @@ if verbose
skipline()
end
model_dynamic = str2func([M.fname,'_dynamic']);
model_dynamic = str2func([M.fname,'.dynamic']);
z = Y(find(lead_lag_incidence'));
[d1,jacobian] = model_dynamic(z, exogenousvariables, params, steadystate,maximum_lag+1);

View File

@ -48,7 +48,7 @@ isf1 = [nyp+ny+1:nyf+nyp+ny+1];
stop = false;
iz = [1:ny+nyp+nyf];
dynamicmodel = sprintf('%s_dynamic', M.fname);
dynamicmodel = sprintf('%s.dynamic', M.fname);
verbose = options.verbosity;

View File

@ -107,7 +107,7 @@ if verbose
skipline()
end
dynamicmodel = str2func([M.fname,'_dynamic']);
dynamicmodel = str2func([M.fname,'.dynamic']);
z = steadystate_y([ip; ic; in]);
x = repmat(transpose(steadystate_x), 1+M.maximum_exo_lag+M.maximum_exo_lead, 1);

View File

@ -34,7 +34,7 @@ if ny0 ~= M.endo_nbr
error('All endogenous variables must appear at the current period!')
end
dynamicmodel = str2func([M.fname,'_dynamic']);
dynamicmodel = str2func([M.fname,'.dynamic']);
info.status = 1;

View File

@ -25,7 +25,7 @@ if ny0 ~= M.endo_nbr
error('All endogenous variables must appear at the current period!')
end
dynamicmodel = str2func([M.fname,'_dynamic']);
dynamicmodel = str2func([M.fname,'.dynamic']);
info.status = 1;

View File

@ -125,16 +125,16 @@ else
% Global variables for parallel routines.
globalVars = struct();
% which files have to be copied to run remotely
NamFileInput(1,:) = {'',[ModelName '_static.m']};
NamFileInput(2,:) = {'',[ModelName '_dynamic.m']};
NamFileInput(1,:) = {'',[ModelName '.static.m']};
NamFileInput(2,:) = {'',[ModelName '.dynamic.m']};
if M_.set_auxiliary_variables
NamFileInput(3,:) = {'',[M_.fname '_set_auxiliary_variables.m']};
NamFileInput(3,:) = {'',[M_.fname '.set_auxiliary_variables.m']};
end
if options_.steadystate_flag
if options_.steadystate_flag == 1
NamFileInput(length(NamFileInput)+1,:)={'',[M_.fname '_steadystate.m']};
else
NamFileInput(length(NamFileInput)+1,:)={'',[M_.fname '_steadystate2.m']};
NamFileInput(length(NamFileInput)+1,:)={'',[M_.fname '.steadystate.m']};
end
end
if (options_.load_mh_file~=0) && any(fline>1)

View File

@ -292,16 +292,16 @@ else
'estim_params_', estim_params_, ...
'oo_', oo_);
% which files have to be copied to run remotely
NamFileInput(1,:) = {'',[M_.fname '_static.m']};
NamFileInput(2,:) = {'',[M_.fname '_dynamic.m']};
NamFileInput(1,:) = {'',[M_.fname '.static.m']};
NamFileInput(2,:) = {'',[M_.fname '.dynamic.m']};
if M.set_auxiliary_variables
NamFileInput(3,:) = {'',[M_.fname '_set_auxiliary_variables.m']};
NamFileInput(3,:) = {'',[M_.fname '.set_auxiliary_variables.m']};
end
if options_.steadystate_flag
if options_.steadystate_flag == 1
NamFileInput(length(NamFileInput)+1,:)={'',[M_.fname '_steadystate.m']};
NamFileInput(length(NamFileInput)+1,:)={'',[M_.fname '.steadystate.m']};
else
NamFileInput(length(NamFileInput)+1,:)={'',[M_.fname '_steadystate2.m']};
NamFileInput(length(NamFileInput)+1,:)={'',[M_.fname '.steadystate.m']};
end
end
[fout] = masterParallel(options_.parallel, 1, B,NamFileInput,'prior_posterior_statistics_core', localVars,globalVars, options_.parallel_info);

View File

@ -70,7 +70,7 @@ end
if options_.block && ~options_.bytecode
z = zeros(M_.endo_nbr,1);
for i = 1:length(M_.block_structure_stat.block)
[r, g, yy, var_indx] = feval([M_.fname '_static'],...
[r, g, yy, var_indx] = feval([M_.fname '.static'],...
i,...
oo_.steady_state,...
[oo_.exo_steady_state; ...
@ -82,7 +82,7 @@ elseif options_.bytecode
[check, z] = bytecode('evaluate','static');
mexErrCheck('bytecode', check);
else
z = feval([M_.fname '_static'],...
z = feval([M_.fname '.static'],...
oo_.steady_state,...
[oo_.exo_steady_state; ...
oo_.exo_det_steady_state], M_.params);

View File

@ -25,7 +25,7 @@ ss = oo_.steady_state;
ss(indx) = y;
eval(['[R,G] = ' M_.fname '_static(ss, x, M_.params);']);
eval(['[R,G] = ' M_.fname '.static(ss, x, M_.params);']);
sR = R(inde);
sG = G(inde,indx);

View File

@ -115,7 +115,7 @@ if local_order == 1
M_.params, dr.ys, 1);
jacobia_ = [loc_dr.g1 loc_dr.g1_x loc_dr.g1_xd];
else
[junk,jacobia_] = feval([M_.fname '_dynamic'],z(iyr0),exo_simul, ...
[junk,jacobia_] = feval([M_.fname '.dynamic'],z(iyr0),exo_simul, ...
M_.params, dr.ys, it_);
end
elseif local_order == 2
@ -124,7 +124,7 @@ elseif local_order == 2
M_.params, dr.ys, 1);
jacobia_ = [loc_dr.g1 loc_dr.g1_x];
else
[junk,jacobia_,hessian1] = feval([M_.fname '_dynamic'],z(iyr0),...
[junk,jacobia_,hessian1] = feval([M_.fname '.dynamic'],z(iyr0),...
exo_simul, ...
M_.params, dr.ys, it_);
end

View File

@ -50,7 +50,7 @@ if flagmoments==0
elseif flagmoments==-1
[I,J]=find(M_.lead_lag_incidence');
yy0=oo_.dr.ys(I);
[residual, g1] = feval([M_.fname,'_dynamic'],yy0, oo_.exo_steady_state', ...
[residual, g1] = feval([M_.fname,'.dynamic'],yy0, oo_.exo_steady_state', ...
M_.params, oo_.dr.ys, 1);
tau=[oo_.dr.ys(oo_.dr.order_var); g1(:)];

View File

@ -32,8 +32,8 @@ function dyn_mex(win_compiler,basename,force)
% You should have received a copy of the GNU General Public License
% along with Dynare. If not, see <http://www.gnu.org/licenses/>.
Dc = dir([basename '_dynamic.c']);
Dmex = dir([basename '_dynamic.' mexext]);
Dc = dir([basename '/mode/src/dynamic.c']);
Dmex = dir(['+' basename '/model/dynamic.' mexext]);
% compile only if date of C file is greater than date of mex file
% and force is not True
@ -50,19 +50,19 @@ if ~exist('OCTAVE_VERSION')
if strcmp(win_compiler,'msvc')
% MATLAB/Windows + Microsoft Visual C++
% Add /TP flag as fix for #1227
eval(['mex -O LINKFLAGS="$LINKFLAGS /export:Dynamic" COMPFLAGS="/TP" ' basename '_dynamic.c ' basename '_dynamic_mex.c'])
eval(['mex -O LINKFLAGS="$LINKFLAGS /export:Static" COMPFLAGS="/TP" ' basename '_static.c ' basename '_static_mex.c'])
eval(['mex -O LINKFLAGS="$LINKFLAGS /export:Dynamic" COMPFLAGS="/TP" ' basename '/model/src/dynamic.c ' basename '/model/src/dynamic_mex.c -output +' basename '/dynamic'])
eval(['mex -O LINKFLAGS="$LINKFLAGS /export:Static" COMPFLAGS="/TP" ' basename '/model/src/static.c ' basename '/model/src/static_mex.c -output +' basename '/static'])
elseif strcmp(win_compiler,'mingw')
eval(['mex -O LINKFLAGS="$LINKFLAGS /export:Dynamic" ' basename '_dynamic.c ' basename '_dynamic_mex.c'])
eval(['mex -O LINKFLAGS="$LINKFLAGS /export:Static" ' basename '_static.c ' basename '_static_mex.c'])
eval(['mex -O LINKFLAGS="$LINKFLAGS /export:Dynamic" ' basename '/model/src/dynamic.c ' basename '/model/src/dynamic_mex.c -output +' basename '/dynamic'])
eval(['mex -O LINKFLAGS="$LINKFLAGS /export:Static" ' basename '/model/src/static.c ' basename '/model/src/static_mex.c -output +' basename '/static'])
elseif strcmp(win_compiler,'cygwin') %legacy support for Cygwin with mexopts.bat
% MATLAB/Windows + Cygwin g++
eval(['mex -O PRELINK_CMDS1="echo EXPORTS > mex.def & echo ' ...
'mexFunction >> mex.def & echo Dynamic >> mex.def" ' ...
basename '_dynamic.c ' basename '_dynamic_mex.c'])
basename '/model/src/dynamic.c ' basename '/model/src/dynamic_mex.c -output +' basename '/dynamic'])
eval(['mex -O PRELINK_CMDS1="echo EXPORTS > mex.def & echo ' ...
'mexFunction >> mex.def & echo Dynamic >> mex.def" ' ...
basename '_static.c ' basename '_static_mex.c'])
basename '/model/src/static.c ' basename '/model/src/static_mex.c -output +' basename '/static'])
else
error(['When using the USE_DLL option, you must give either ' ...
'''cygwin'', ''mingw'' or ''msvc'' option to the ''dynare'' command'])
@ -71,15 +71,15 @@ if ~exist('OCTAVE_VERSION')
% MATLAB/Linux
if matlab_ver_less_than('8.3')
eval(['mex -O LDFLAGS=''-pthread -shared -Wl,--no-undefined'' ' ...
basename '_dynamic.c ' basename '_dynamic_mex.c'])
basename '/model/src/dynamic.c ' basename '/model/src/dynamic_mex.c -output +' basename '/dynamic'])
eval(['mex -O LDFLAGS=''-pthread -shared -Wl,--no-undefined'' ' ...
basename '_static.c ' basename '_static_mex.c'])
basename '/model/src/static.c ' basename '/model/src/static_mex.c -output +' basename '/static'])
elseif matlab_ver_less_than('9.1')
eval(['mex -O LINKEXPORT='''' ' basename '_dynamic.c ' basename '_dynamic_mex.c'])
eval(['mex -O LINKEXPORT='''' ' basename '_static.c ' basename '_static_mex.c'])
eval(['mex -O LINKEXPORT='''' ' basename '/model/src/dynamic.c ' basename '/model/src/dynamic_mex.c -output +' basename '/dynamic'])
eval(['mex -O LINKEXPORT='''' ' basename '/model/src/static.c ' basename '/model/src/static_mex.c -output +' basename '/static'])
else
eval(['mex -O LINKEXPORTVER='''' ' basename '_dynamic.c ' basename '_dynamic_mex.c'])
eval(['mex -O LINKEXPORTVER='''' ' basename '_static.c ' basename '_static_mex.c'])
eval(['mex -O LINKEXPORTVER='''' ' basename '/model/src/dynamic.c ' basename '/model/src/dynamic_mex.c -output +' basename '/dynamic'])
eval(['mex -O LINKEXPORTVER='''' ' basename '/model/src/static.c ' basename '/model/src/static_mex.c -output +' basename '/static'])
end
elseif ismac
% MATLAB/MacOS
@ -87,30 +87,30 @@ if ~exist('OCTAVE_VERSION')
eval(['mex -O LDFLAGS=''-Wl,-twolevel_namespace -undefined ' ...
'error -arch $ARCHS -Wl,-syslibroot,$SDKROOT ' ...
'-mmacosx-version-min=$MACOSX_DEPLOYMENT_TARGET -bundle'' ' ...
basename '_dynamic.c ' basename '_dynamic_mex.c'])
basename '/model/src/dynamic.c ' basename '/model/src/dynamic_mex.c -output +' basename '/dynamic'])
eval(['mex -O LDFLAGS=''-Wl,-twolevel_namespace -undefined ' ...
'error -arch $ARCHS -Wl,-syslibroot,$SDKROOT ' ...
'-mmacosx-version-min=$MACOSX_DEPLOYMENT_TARGET -bundle'' ' ...
basename '_static.c ' basename '_static_mex.c'])
basename '/model/src/static.c ' basename '/model/src/static_mex.c -output +' basename '/static'])
elseif matlab_ver_less_than('8.3')
eval(['mex -O LDFLAGS=''-Wl,-twolevel_namespace -undefined ' ...
'error -arch $ARCHS -Wl,-syslibroot,$MW_SDKROOT ' ...
'-mmacosx-version-min=$MACOSX_DEPLOYMENT_TARGET -bundle'' ' ...
basename '_dynamic.c ' basename '_dynamic_mex.c'])
basename '/model/src/dynamic.c ' basename '/model/src/dynamic_mex.c -output +' basename '/dynamic'])
eval(['mex -O LDFLAGS=''-Wl,-twolevel_namespace -undefined ' ...
'error -arch $ARCHS -Wl,-syslibroot,$MW_SDKROOT ' ...
'-mmacosx-version-min=$MACOSX_DEPLOYMENT_TARGET -bundle'' ' ...
basename '_static.c ' basename '_static_mex.c'])
basename '/model/src/static.c ' basename '/model/src/static_mex.c -output +' basename '/static'])
elseif matlab_ver_less_than('9.1')
eval(['mex -O LINKEXPORT='''' ' basename '_dynamic.c ' basename '_dynamic_mex.c'])
eval(['mex -O LINKEXPORT='''' ' basename '_static.c ' basename '_static_mex.c'])
eval(['mex -O LINKEXPORT='''' ' basename '/model/src/dynamic.c ' basename '/model/src/dynamic_mex.c -output +' basename '/dynamic'])
eval(['mex -O LINKEXPORT='''' ' basename '/model/src/static.c ' basename '/model/src/static_mex.c -output +' basename '/static'])
else
eval(['mex -O LINKEXPORT='''' LINKEXPORTVER='''' ' basename '_dynamic.c ' basename '_dynamic_mex.c'])
eval(['mex -O LINKEXPORT='''' LINKEXPORTVER='''' ' basename '_static.c ' basename '_static_mex.c'])
eval(['mex -O LINKEXPORT='''' LINKEXPORTVER='''' ' basename '/model/src/dynamic.c ' basename '/model/src/dynamic_mex.c -output +' basename '/dynamic'])
eval(['mex -O LINKEXPORT='''' LINKEXPORTVER='''' ' basename '/model/src/static.c ' basename '/model/src/static_mex.c -output +' basename '/static'])
end
end
else
% Octave
eval(['mex ' basename '_dynamic.c ' basename '_dynamic_mex.c'])
eval(['mex ' basename '_static.c ' basename '_static_mex.c'])
eval(['mex ' basename '/model/src/dynamic.c ' basename '/model/src/dynamic_mex.c -o +' basename '/dynamic'])
eval(['mex ' basename '/model/src/static.c ' basename '/model/src/static_mex.c -o +' basename '/static'])
end

View File

@ -577,9 +577,9 @@ void
Interpreter::ReadCodeFile(string file_name, CodeLoad &code)
{
if (steady_state)
file_name += "_static";
file_name += "/model/bytecode/static";
else
file_name += "_dynamic";
file_name += "/model/bytecode/dynamic";
//First read and store in memory the code
code_liste = code.get_op_code(file_name);

View File

@ -524,23 +524,22 @@ dynSparseMatrix::Read_SparseMatrix(string file_name, const int Size, int periods
{
unsigned int eq, var;
int lag;
filename = file_name;
mem_mngr.fixe_file_name(file_name);
/*mexPrintf("steady_state=%d, size=%d, solve_algo=%d, stack_solve_algo=%d, two_boundaries=%d\n",steady_state, Size, solve_algo, stack_solve_algo, two_boundaries);
mexEvalString("drawnow;");*/
if (!SaveCode.is_open())
{
if (steady_state)
SaveCode.open((file_name + "_static.bin").c_str(), ios::in | ios::binary);
SaveCode.open(file_name + "/model/bytecode/static.bin", ios::in | ios::binary);
else
SaveCode.open((file_name + "_dynamic.bin").c_str(), ios::in | ios::binary);
SaveCode.open(file_name + "/model/bytecode/dynamic.bin", ios::in | ios::binary);
if (!SaveCode.is_open())
{
ostringstream tmp;
if (steady_state)
tmp << " in Read_SparseMatrix, " << file_name << "_static.bin cannot be opened\n";
tmp << " in Read_SparseMatrix, " << file_name << "/model/bytecode/static.bin cannot be opened\n";
else
tmp << " in Read_SparseMatrix, " << file_name << "_dynamic.bin cannot be opened\n";
tmp << " in Read_SparseMatrix, " << file_name << "/model/bytecode/dynamic.bin cannot be opened\n";
throw FatalExceptionHandling(tmp.str());
}
}

View File

@ -27,7 +27,7 @@ DynamicModelDLL::DynamicModelDLL(const string &modName) throw (DynareException)
#if !defined(__CYGWIN32__) && !defined(_WIN32)
fName = "./";
#endif
fName += modName + "_dynamic" + MEXEXT;
fName += "+" + modName + "/dynamic" + MEXEXT;
try
{

View File

@ -20,7 +20,7 @@
#include "dynamic_m.hh"
DynamicModelMFile::DynamicModelMFile(const string &modName) throw (DynareException) :
DynamicMFilename(modName + "_dynamic")
DynamicMFilename(modName + ".dynamic")
{
}

@ -1 +1 @@
Subproject commit a2b19cbffff9612aea4d255fd7c6df0591087da8
Subproject commit e376267a2867aff79989503f77a2fd2204ed9c8b

View File

@ -935,15 +935,6 @@ clean-local:
$(patsubst %.trs, %.json, $(O_TRS_FILES)) \
$(patsubst %.trs, %.json, $(O_XFAIL_TRS_FILES))
rm -f $(patsubst %.mod, %.m, $(MODFILES)) \
$(patsubst %.mod, %_static.*, $(MODFILES)) \
$(patsubst %.mod, %_static_*.m, $(MODFILES)) \
$(patsubst %.mod, %_objective_static.m, $(MODFILES)) \
$(patsubst %.mod, %_set_auxiliary_variables.m, $(MODFILES)) \
$(patsubst %.mod, %_steadystate2.m, $(MODFILES)) \
$(patsubst %.mod, %_dynamic.*, $(MODFILES)) \
$(patsubst %.mod, %_dynamic_*.m, $(MODFILES))
rm -f $(patsubst %.mod, %_results.mat, $(MODFILES)) \
$(patsubst %.mod, %_mode.mat, $(MODFILES)) \
$(patsubst %.mod, %_mh_mode.mat, $(MODFILES)) \
@ -955,6 +946,8 @@ clean-local:
rm -rf $(patsubst %.mod, %, $(MODFILES))
rm -rf $(foreach mod, $(MODFILES), $(dir $(mod))+$(basename $(notdir $(mod))))
rm -f $(patsubst %.mod, %*.pdf, $(MODFILES)) \
$(patsubst %.mod, %*.eps, $(MODFILES)) \
$(patsubst %.mod, %*.fig, $(MODFILES))
@ -1013,4 +1006,3 @@ clean-local:
find . -name "*.aux" -type f -delete
find . -name "*.log" -type f -delete
find . -name "*.eps" -type f -delete
find . -name "*.json" -type f -delete

View File

@ -38,6 +38,6 @@ perfect_foresight_solver;
rplot c;
rplot k;
if ~exist('./ramst.json', 'file') || exist('./ramst.log', 'file')
if ~exist('./ramst/model/json/modfile.json', 'file') || exist('./ramst.log', 'file')
error('The dynare command did not honor the options provided in the mod file!')
end

View File

@ -56,7 +56,7 @@ end;
// Write analytical steady state file (without globals)
options_.steadystate_flag = 2;
copyfile('rbcii_steady_state.m','rbcii_steadystate2.m');
copyfile('rbcii_steady_state.m','+rbcii/steadystate.m');
@#if extended_path_version

View File

@ -69,11 +69,11 @@ end;
simul(periods=1000);
newton_solution_is_wrong = abs(evaluate_max_dynamic_residual(str2func('sw_newton_dynamic'), oo_.endo_simul, oo_.exo_simul, M_.params, oo_.steady_state, 1000, size(oo_.endo_simul, 1), 1, M_.lead_lag_incidence))>options_.dynatol.f;
newton_solution_is_wrong = abs(evaluate_max_dynamic_residual(str2func('sw_newton.dynamic'), oo_.endo_simul, oo_.exo_simul, M_.params, oo_.steady_state, 1000, size(oo_.endo_simul, 1), 1, M_.lead_lag_incidence))>options_.dynatol.f;
lmmcp = load('sw_lmmcp_results');
lmmcp_solution_is_wrong = abs(evaluate_max_dynamic_residual(str2func('sw_newton_dynamic'), lmmcp.oo_.endo_simul, lmmcp.oo_.exo_simul, M_.params, oo_.steady_state, 1000, size(oo_.endo_simul, 1), 1, M_.lead_lag_incidence))>options_.dynatol.f;
lmmcp_solution_is_wrong = abs(evaluate_max_dynamic_residual(str2func('sw_newton.dynamic'), lmmcp.oo_.endo_simul, lmmcp.oo_.exo_simul, M_.params, oo_.steady_state, 1000, size(oo_.endo_simul, 1), 1, M_.lead_lag_incidence))>options_.dynatol.f;
solutions_are_different = max(max(abs(lmmcp.oo_.endo_simul-oo_.endo_simul)))>options_.dynatol.x;