Forecasting

This page describes forecasting without regime-switching. Click here to learn about forecasting with regime-switching.

Procedure

In the forecast step, we compute smoothed histories, forecast, compute shock decompositions, and compute impulse response functions (IRFs) for states, observables, shocks, and pseudo-observables. To run a forecast on one combination of input parameter type (e.g. modal parameters or full-distribution) and conditional type, call forecast_one. The forecast output is written to the saveroot specified by the model object and can be loaded with the function read_forecast_output.

Main Steps:

  • Prepare forecast inputs: Add required output types, load data, and load draws of parameter vectors saved from the estimation step.

  • Compute forecast outputs: Carry out desired combination of smoothing, forecasting, computing shock decompositions, and computing IRFs. See Forecast Outputs for a list of possible forecast outputs.

  • Save forecast outputs: Save each forecast output as an array to its own file, along with some metadata.

DSGE.forecast_oneFunction
forecast_one(m, input_type, cond_type, output_vars; df = DataFrame(),
    subset_inds = 1:0, forecast_string = "", verbose = :low, ...)

Compute and save output_vars for input draws given by input_type and conditional data case given by cond_type.

Inputs

  • m::AbstractDSGEModel: model object

  • input_type::Symbol: one of:

  - `:mode`: forecast using the modal parameters only
  - `:mean`: forecast using the mean parameters only
  - `:init`: forecast using the initial parameter values only
  - `:full`: forecast using all parameters (full distribution)
  - `:subset`: forecast using a well-defined user-specified subset of draws
  • cond_type::Symbol: one of:
  - `:none`: no conditional data
  - `:semi`: use "semiconditional data" - average of quarter-to-date
    observations for high frequency series
  - `:full`: use "conditional data" - semiconditional plus nowcasts for
    desired observables
  • output_vars::Vector{Symbol}: vector of desired output variables. See ?forecast_one_draw.

Keyword Arguments

  • df::DataFrame: Historical data. If cond_type in [:semi, :full], then the final row of df should be the period containing conditional data. If not provided, will be loaded using load_data with the appropriate cond_type

  • subset_inds::AbstractRange{Int64}: indices specifying the draws we want to use. If a more sophisticated selection criterion is desired, the user is responsible for determining the indices corresponding to that criterion. If input_type is not subset, subset_inds will be ignored

  • forecast_string::String: short string identifying the subset to be appended to the output filenames. If input_type = :subset and forecast_string is empty, an error is thrown.

  • only_filter::Bool: do not run the smoother and only run the filter. This limits the number of output variables which can be calculated.

  • verbose::Symbol: desired frequency of function progress messages printed to standard out. One of :none, :low, or :high.

  • check_empty_columns::Bool = true: check empty columns or not when loading data (if df is empty)

  • bdd_fcast::Bool = true: are we computing the bounded forecasts or not?

  • params::AbstractArray{Float64} = Vector{Float64}(undef, 0): parameter draws for the forecast. If empty, then we load draws from estimation files implied by the settings in m.

  • zlb_method::Symbol: method for enforcing the zero lower bound. Defaults to :shock, meaning we use a monetary policy shock to enforce the ZLB.

    Other available methods:

    1. :temporary_altpolicy -> use a temporary alternative policy to enforce the ZLB.
  • rerun_smoother::Bool = false: if true, rerun the conditional forecast when automatically enforcing the ZLB as a temporary alternative policy.

  • nan_endozlb_failures::Bool = false: if true, failures of an endogenous ZLB (when implemented as a temporary policy) will be handled by throwing NaNs instead of using unanticipated monetary policy shocks.

  • set_regime_vals_altpolicy::Function: Function that adds new regimes to parameters when using temporary alternative policies (if needed). Defaults to identity (which does nothing) This function should take as inputs the model object m and the total number of regimes (after adding the required temporary regimes), i.e. set_regime_vals_altpolicy(m, n). It should then set up regime-switching parameters for these new additional regimes.

  • set_info_sets_altpolicy::Function = auto_temp_altpolicy_info_set: Function that automatically updates the tvis_information_set, e.g. when zlb_method = :temporary_altpolicy.

  • update_regime_eqcond_info!::Function = (x1, x2, x3, x4) -> default_update_regime_eqcond_info(x1, x2, x3, x4, alternative_policy(m): Function that automatically updates the setting :regime_eqcond_info. The arguments of update_regime_eqcond_info! should be (in order) m::AbstractDSGEModel, eqcond_dict::AbstractDict{Int64, EqcondEntry}, zlb_start_regime::Int64, and liftoff_regime::Int64. The last two arguments are the regime numbers of the first regime for which the ZLB applies and the regime after the ZLB ends, respectively. The eqcond_dict argument should specify the EqcondEntry during the historical/conditional horizon regime (if it is desired) but can otherwise be empty. This function should then update eqcond_dict in place to implement a temporary ZLB and any other permanent alternative policies/regime-switching in the forecast horizon (after the conditional horizon). The user should also be careful and make sure update_regime_eqcond_info! handles imperfect awareness properly if they want to implement an imperfectly credible ZLB. For an example, see ?default_update_regime_eqcond_info.

  • show_failed_percent::Bool = false: prints out the number of failed forecasts, which are returned as NaNs. These may occur when the ZLB is not enforced, for example.

  • pegFFR::Bool = false: peg the nominal FFR at the value specified by FFRpeg

  • FFRpeg::Float64 = -0.25/4: value of the FFR peg

  • H::Int = 4: number of horizons for which the FFR is pegged

  • testing_carer_kohn::Bool = false: whether to create a file storing some property of Σ in Carter Kohn

Outputs

None. Output is saved to files returned by get_forecast_output_files(m, input_type, cond_type, output_vars).

source

For example, to do an unconditional forecast of states and observables using the modal parameters, call:

m = AnSchorfheide()
forecast_one(m, :mode, :none, [:forecaststates, forecastobs])

Full-Distribution Forecasts:

Full-distribution forecasts are computed in blocks. The size of each block defaults to 5000 draws (before thinning by get_setting(m, :forecast_jstep)), but can be set using the :forecast_block_size Setting. For each block, draws are read in on the originator process, then computation proceeds in parallel using pmap. When all draws in the block are finished, the forecast outputs are reassembled on the originator process and appended to the HDF5 dataset in their respective output files.

To fully take advantage of the parallelization, the user is responsible for adding processes before calling forecast_one, either by calling addprocs or using one of the functions defined in ClusterManagers.jl. For example, to run a full-distribution unconditional forecast using 10 processes:

my_procs = addprocs(10)
@everywhere using DSGE

m = AnSchorfheide()
forecast_one(m, :full, :none, [:forecaststates, forecastobs])

rmprocs(my_procs)

Notice that it is necessary to load DSGE on all processes using @everywhere using DSGE before calling forecast_one. It is also sometimes necessary to load OrderedCollections on all processes using @everywhere using OrderedCollections.

By default, full-distribution forecasts start from the first block. However, if you want to start the forecast from a later block, you can also do so. For example:

m <= Setting(:forecast_start_block, 2,
    "Block at which to resume forecasting (possibly null)")

Forecast Outputs

A forecast output (i.e. an output_var) is a combination of what we call a "product" and a "class". The possible classes are states (:states), observables (:obs), pseudo-observables (:pseudo), and standardized (:stdshocks) and unstandardized shocks (:shocks). The possible forecast products are:

  • Smoothed histories (:hist): use the smoother specified by forecast_smoother(m) to get smoothed histories of each class.

  • Forecasts (:forecast): iterate the state space forward from the last filtered state, either using a specified set of shock innovations or by drawing these from a distribution. Forecasts in which we enforce the zero lower bound are denoted as :bddforecast.

  • Shock decompositions (:shockdec): decompose a forecast into the contributions from each individual shock (e.g. TFP or monetary policy) or from groups of shocks by exploiting the linearity of the state space system. Starting from an initial state of zero, iterate the state space forward from the first historical period up through the last forecast horizon. To infer the contributions of shocks, use the smoothed historical shocks for one shock at a time during the historical periods. If a modal shock decomposition is computed, then we assume no shocks occur during the forecast periods. If a full-distribution shock decomposition is computed, then we allow shocks to occur during the forecast periods.

  • Deterministic trends (:dettrend): iterate the state space forward from the first historical state up through the last forecast horizon without any shocks.

  • Trends (:trend): for each class, just the constant term in that class's equation, i.e. the CCC vector from the transition equation for states, the DD vector from the measurement equation for observables, and the DD_pseudo vector from the pseuodo-measurement equation for pseudo-observables.

  • IRFs (:irf): see Impulse response. Our IRFs are in response to a shock of size -1 standard deviation.

An output_var is then just a Symbol with a product and class concatenated, e.g. :histstates for smoothed historical states.

It is not necessary to compute all forecast outputs in one call to forecast_one. Which steps are run depends on which output_vars are passed in. Most users will just want to calculate :forecastobs, :histobs, :forecastpseudo, and :histpseudo.

Calculating Shock Decompositions

The user should note that shock decompositions are not the same thing as forecast decompositions in DSGE.jl. The former decomposes a forecast into the contributions of individual or groups of shocks. The latter decomposes the change between two different forecasts into "data revision", "news" and "re-estimation" components.

If the user wants to plot shock decompositions to decompose a forecast's deviations from trend into the contributions of individual or groups of shocks, then the user needs to calculate the shockdec, detttrend, and trend output variables. Otherwise, the plotting script plot_shock_decompositions will fail.

Finally, note that full-distribution shock decompositions are very memory intensive. For Model1002, we typically need 30GB-40GB using the default settings. To avoid memory problems, the user should decrease the setting :forecast_block_size to a smaller number than the default 5000. A smaller block size means that fewer draws from the posterior are loaded into memory at any one time, hence the size of the matrices associated with shock decompositions will also be smaller.

Preparing Forecast Inputs

Adding Required output_vars:

The user specifies what output to compute by passing a Vector{Symbol} of correctly named output variables into forecast_one. Underneath the hood, we will add any additional output variables that are required to complete a forecast using add_requisite_output_vars:

  • If :forecast<class> is in output_vars, then :bddforecast<class> is also added. Hence we always forecast both with and without enforcing the ZLB.
  • If :shockdec<class> is in output_vars, then :dettrend<class> and :trend<class> are also added. This is because to plot shock decompositions, we also need the trend and the deterministic trend.

Loading Data:

This is done the usual way, using load_data with the appropriate cond_type.

Note that if you are running a conditional forecast, then you should check that the appropriate observables are found in cond_full_names(m) and/or cond_semi_names(m). If the observables on which you are conditioning do not match the ones found in these settings exactly, then conditional forecasting will not generate the results you want. To change these observables to, for example, condition on custom_obs1 and custom_obs2, add the line

m <= Setting(:cond_full_names, [custom_obs1, custom_obs2])

to your script for full conditional forecasts. The syntax is similar for semi-conditional forecasts. If you have already created a conditional dataset successfully, then you should be sure to re-create the data set by running

load_data(m; try_disk = false, cond_type = cond_type)

If you do not have try_disk = false, then the code may load a saved data set that does not have the conditional observables you want to use.

Loading Draws:

By default, the draws are loaded from the file whose path is given by get_forecast_input_file. However, you can directly pass a matrix of posterior parameter draws or a vector of the modal parameters via the keyword argument params for forecast_one, or you can override the default input file for a given input type by adding entries to the Dict{Symbol, ASCIIString} returned from forecast_input_file_overrides(m). For example:

overrides = forecast_input_file_overrides(m)
overrides[:mode] = "path/to/input/file.h5"

Note that load_draws expects an HDF5 dataset called either params (for input_type in [:mode, :mean]) or mhparams (for input_type in [:full, :subset]) if the sampling method is :MH. If the sampling method is :SMC, then load_draws expects a jld2 file which has a variable named cloud.

Computing Forecast Outputs

For each draw of parameters, the forecast calculations are run by the lower-level function forecast_one_draw.

This function is also useful on its own when the user wants to run a single-draw forecast without writing the output to data. A common use case is running experiments with different forecast specifications and/or alternative policy rules. Given the appropriate inputs, forecast_one_draw will return a Dict whose keys are the names of output variables (e.g. :forecastobs) and values are the corresponding matrices. This function does not perform transformations, so the units of the output are all model units, which are typically at the quarterly frequency, such as quarterly per-capita GDP growth.

The computations that forecast_one_draw can run are:

Smoothing:

Smoothing is necessary if either:

  • You explicitly want the smoothed histories, or
  • You want to compute shock decompositions or deterministic trends, which use the smoothed historical shocks

It is not necessary to keep track of these cases, however - forecast_one will deduce from the specified output_vars whether or not it is necessary to filter and smooth in order to produce your output_vars.

Forecasting:

Forecasting begins from the last filtered historical state, which is obtained from the Kalman filter. forecast accepts a keyword argument enforce_zlb, which indicates whether to enforce the zero lower bound. If enforce_zlb = true, then if in a given period, the forecasted interest rate goes below forecast_zlb_value(m), we solve for the interest rate shock necessary to push it up to the ZLB. A forecast in which the ZLB is enforced corresponds to the product :bddforecast.

Shock Decompositions, Deterministic Trends, and Trends:

Since shock decompositions have an additional dimension (e.g. nstates x nperiods x nshocks for a single draw of state shock decompositions, compared to nstates x nperiods for a single draw of forecasted states), we usually wish to truncate some periods before returning. This behavior is governed by the Settings :shockdec_starttdate and :shockdec_enddate, which are of type Nullable{Date}.

Deterministic trends are also saved only for date_shockdec_start(m) and date_shockdec_end(m). Trends are not time-dependent.

Note that shock decompositions are memory intensive. For a typical forecast of Model1002 without shock decompositions, only 1-2GB of memory are needed. For a forecast with shock decompositions, roughtly 14-16GB will be needed.

Impulse Response Functions:

Like shock decompositions, IRFs have three dimensions (e.g. nstates x nperiods x nshocks) for each draw.

Saving Forecast Outputs

Forecast outputs are saved in the location specified by get_forecast_output_files(m), which is typically a subdirectory of saveroot(m). Each output_var is saved in its own JLD file, which contains the following datasets:

  • arr::Array: actual array of forecast outputs. For trends, this array is of size ndraws x nvars. For histories, forecasts, and deterministic trends, it is ndraws x nvars x nperiods. For shock decompositions and IRFs, it is ndraws x nvars x nperiods x nshocks. (In all of these, nvars refers to the number of variables of the output class.)

  • date_indices::Dict{Date, Int}: maps Dates to their indices along the nperiods dimension of arr. Not saved for IRFs.

  • <class>_names::Dict{Symbol, Int}: maps names of variables of the output class (e.g. :OutputGap) into their indices along the nvars dimension of arr.

  • <class>_revtransforms::Dict{Symbol, Symbol}: maps names of variables to the names of the reverse transforms (from model units into plotting units) associated with those variables. For example, pseudoobservable_revtransforms[:π_t] = :quartertoannual.

  • shock_names::Dict{Symbol, Int}: for shock decompositions and IRFs only, maps names of shocks into their indices along the nshocks dimension of arr.

Some helpful functions for getting file names, as well as reading and writing forecast outputs, include:

  • get_forecast_input_file
  • get_forecast_filename
  • get_forecast_output_files
  • write_forecast_outputs
  • write_forecast_block
  • write_forecast_metadata
  • read_forecast_metadata
  • read_forecast_output

Forecasting Functions

DSGE.decompose_forecastMethod
decompose_forecast(m_new, m_old, df_new, df_old, input_type, cond_new, cond_old,
    classes; verbose = :low, kwargs...)

decompose_forecast(m_new, m_old, df_new, df_old, params_new, params_old,
    cond_new, cond_old, classes; check = false)

explains the differences between an old forecast and a new forecast by decomposing the differences into three sources:

(1) Data revisions, (2) News (e.g. new data that has become available since the old forecast), (3) Re-estimation (i.e. changes in model parameters).

This function does not compute which shocks explain a forecast. For example, if you want to know whether TFP or financial shocks drive a given forecast, then you want to compute the shock decomposition output variable (see ?shock_decompositions, forecast_one, and compute_meansbands).

Note that this function currently does not work for a model in which there are changes in the degree of "regime-switching" in the TTT, RRR, CCC, ZZ, and DD matrices, e.g. decomposing the changes in the forecast when the monetary policy rule changes or if a temporary policy is implemented that did not occur in the old forecast.

Inputs

  • m_new::M and m_old::M where M<:AbstractDSGEModel
  • df_new::DataFrame and df_old::DataFrame
  • cond_new::Symbol and cond_old::Symbol
  • classes::Vector{Symbol}: some subset of [:states, :obs, :pseudo]

Method 1 only:

  • input_type::Symbol: estimation type to use. Parameters will be loaded using load_draws(m_new, input_type) and load_draws(m_old, input_type) in this method

Method 2 only:

  • params_new::Vector{Float64} and params_old::Vector{Float64}: single parameter draws to use

Keyword Arguments

  • check::Bool: whether to check that the individual components add up to the correct total difference in forecasts. This roughly doubles the runtime

Method 1 only:

  • verbose::Symbol

Outputs

The first method returns nothing. The second method returns decomp::Dict{Symbol, Matrix{Float64}}, which has keys of the form :decomp<component><class> and values of size Ny x Nh, where

  • Ny is the number of variables in the given class
  • Nh is the number of common forecast periods, i.e. periods between date_forecast_start(m_new) and date_forecast_end(m_old)
source
DSGE.forecast_oneMethod
forecast_one(m, input_type, cond_type, output_vars; df = DataFrame(),
    subset_inds = 1:0, forecast_string = "", verbose = :low, ...)

Compute and save output_vars for input draws given by input_type and conditional data case given by cond_type.

Inputs

  • m::AbstractDSGEModel: model object

  • input_type::Symbol: one of:

  - `:mode`: forecast using the modal parameters only
  - `:mean`: forecast using the mean parameters only
  - `:init`: forecast using the initial parameter values only
  - `:full`: forecast using all parameters (full distribution)
  - `:subset`: forecast using a well-defined user-specified subset of draws
  • cond_type::Symbol: one of:
  - `:none`: no conditional data
  - `:semi`: use "semiconditional data" - average of quarter-to-date
    observations for high frequency series
  - `:full`: use "conditional data" - semiconditional plus nowcasts for
    desired observables
  • output_vars::Vector{Symbol}: vector of desired output variables. See ?forecast_one_draw.

Keyword Arguments

  • df::DataFrame: Historical data. If cond_type in [:semi, :full], then the final row of df should be the period containing conditional data. If not provided, will be loaded using load_data with the appropriate cond_type

  • subset_inds::AbstractRange{Int64}: indices specifying the draws we want to use. If a more sophisticated selection criterion is desired, the user is responsible for determining the indices corresponding to that criterion. If input_type is not subset, subset_inds will be ignored

  • forecast_string::String: short string identifying the subset to be appended to the output filenames. If input_type = :subset and forecast_string is empty, an error is thrown.

  • only_filter::Bool: do not run the smoother and only run the filter. This limits the number of output variables which can be calculated.

  • verbose::Symbol: desired frequency of function progress messages printed to standard out. One of :none, :low, or :high.

  • check_empty_columns::Bool = true: check empty columns or not when loading data (if df is empty)

  • bdd_fcast::Bool = true: are we computing the bounded forecasts or not?

  • params::AbstractArray{Float64} = Vector{Float64}(undef, 0): parameter draws for the forecast. If empty, then we load draws from estimation files implied by the settings in m.

  • zlb_method::Symbol: method for enforcing the zero lower bound. Defaults to :shock, meaning we use a monetary policy shock to enforce the ZLB.

    Other available methods:

    1. :temporary_altpolicy -> use a temporary alternative policy to enforce the ZLB.
  • rerun_smoother::Bool = false: if true, rerun the conditional forecast when automatically enforcing the ZLB as a temporary alternative policy.

  • nan_endozlb_failures::Bool = false: if true, failures of an endogenous ZLB (when implemented as a temporary policy) will be handled by throwing NaNs instead of using unanticipated monetary policy shocks.

  • set_regime_vals_altpolicy::Function: Function that adds new regimes to parameters when using temporary alternative policies (if needed). Defaults to identity (which does nothing) This function should take as inputs the model object m and the total number of regimes (after adding the required temporary regimes), i.e. set_regime_vals_altpolicy(m, n). It should then set up regime-switching parameters for these new additional regimes.

  • set_info_sets_altpolicy::Function = auto_temp_altpolicy_info_set: Function that automatically updates the tvis_information_set, e.g. when zlb_method = :temporary_altpolicy.

  • update_regime_eqcond_info!::Function = (x1, x2, x3, x4) -> default_update_regime_eqcond_info(x1, x2, x3, x4, alternative_policy(m): Function that automatically updates the setting :regime_eqcond_info. The arguments of update_regime_eqcond_info! should be (in order) m::AbstractDSGEModel, eqcond_dict::AbstractDict{Int64, EqcondEntry}, zlb_start_regime::Int64, and liftoff_regime::Int64. The last two arguments are the regime numbers of the first regime for which the ZLB applies and the regime after the ZLB ends, respectively. The eqcond_dict argument should specify the EqcondEntry during the historical/conditional horizon regime (if it is desired) but can otherwise be empty. This function should then update eqcond_dict in place to implement a temporary ZLB and any other permanent alternative policies/regime-switching in the forecast horizon (after the conditional horizon). The user should also be careful and make sure update_regime_eqcond_info! handles imperfect awareness properly if they want to implement an imperfectly credible ZLB. For an example, see ?default_update_regime_eqcond_info.

  • show_failed_percent::Bool = false: prints out the number of failed forecasts, which are returned as NaNs. These may occur when the ZLB is not enforced, for example.

  • pegFFR::Bool = false: peg the nominal FFR at the value specified by FFRpeg

  • FFRpeg::Float64 = -0.25/4: value of the FFR peg

  • H::Int = 4: number of horizons for which the FFR is pegged

  • testing_carer_kohn::Bool = false: whether to create a file storing some property of Σ in Carter Kohn

Outputs

None. Output is saved to files returned by get_forecast_output_files(m, input_type, cond_type, output_vars).

source
DSGE.load_drawsMethod
load_draws(m, input_type; subset_inds = 1:0, verbose = :low)

load_draws(m, input_type, block_inds; verbose = :low)

Load and return parameter draws from Metropolis-Hastings or SMC.

Inputs

  • m::AbstractDSGEModel: model object
  • input_type::Symbol: one of the options for input_type described in the documentation for forecast_one
  • block_inds::AbstractRange{Int64}: indices of the current block (already indexed by jstep) to be read in. Only used in second method

Keyword Arguments

  • subset_inds::AbstractRange{Int64}: indices specifying the subset of draws to be read in. Only used in first method
  • verbose::Symbol: desired frequency of function progress messages printed to standard out. One of :none, :low, or :high. If :low or greater, prints location of input file.

Outputs

  • params: first method returns a single parameter draw of type Vector{Float64}. Second method returns a Vector{Vector{Float64}} of parameter draws for this block.
source
DSGE.usual_model_forecastFunction
usual_model_forecast(m, input_type, cond_type,
    output_vars = [:histobs, :histpseudo, :forecastobs, :forecastpseudo];
    est_override = "", forecast_string = "",
    density_bands = [0.5, 0.6, 0.7, 0.8, 0.9],
    mb_matrix = false, check_empty_columns = true, params = [])

Forecast, compute means and bands, and optionally (if mb_matrix) convert MeansBands to matrices. If the path est_override is provided, it will be added to forecast_input_file_overrides(m).

See ?forecast_one for descriptions of the keywords.

source
DSGE.usual_model_settings!Method
usual_model_settings!(m, vint; cdvt = vint, dsid = data_id(m), cdid = cond_id(m),
    fcast_date = Dates.lastdayofquarter(Dates.today()),
    altpolicy = AltPolicy(:historical, eqcond, solve))

Apply usual defaults for the following settings:

  • data_vintage and cond_vintage: given by input argument vint
  • date_forecast_start and date_conditional_end: given by kwarg fcast_date
  • use_population_forecast: true
  • alternative_policy: given by input argument altpolicy. If this argument is specified, then altpolicy_settings! and altpolicy.setup are also called.
source
DSGE.compute_scenario_systemMethod
compute_scenario_system(m, scen::Scenario; apply_altpolicy = false)

Given the current model parameters, compute the state-space system corresponding to model m and alternative scenario scen. This function differs from compute_system in that the CCC, DD, and DD_pseudo vectors are set to zero (since we forecast in deviations from baseline) and shocks that are not in scen.instrument_names are zeroed out in the QQ matrix.

source
DSGE.filter_shocks!Function
filter_shocks!(m, scen, system::Scenario)

Given a scenario draw scen, back out the shocks necessary to hit scen.targets and put them into scen.instruments. This function returns forecastshocks, an nshocks x horizon matrix of filtered and smoothed shocks.

This function checks forecast_uncertainty_override(m) for whether to smooth shocks using the simulation smoother.

source
DSGE.forecastMethod
forecast(m, system, z0; cond_type = :none,
    enforce_zlb = false, shocks = Matrix{S}(undef, 0,0))

forecast(m, altpolicy, z0, states, obs, pseudo, shocks)

forecast(system, z0, shocks; enforce_zlb = false)

forecast(m, system, z0, shocks; cond_type = :none, enforce_zlb = false)

The first method produces a forecast, given a state space system, initial state, and shocks, using information about the desired forecast contained in m. It enforces the ZLB by using monetary policy shocks.

The second method is similar but differs in two ways. First, it produces forecasts specifically when an alternative policy is used. Second, it enforces the ZLB by treating it as a temporary alternative policy.

The third and fourth methods are internal functions used by the first two methods.

Inputs

  • system::System{S}: state-space system matrices
  • z0::Vector{S}: state vector in the final historical period
  • shocks::Matrix{S}: nshocks x nperiods matrix of shocks to use when forecasting. Note that in the first method, nperiods doesn't necessarily have to equal forecast_horizons(m); it will be truncated or padded with zeros appropriately

Method 1 and 2 only:

  • m::AbstractDSGEModel

Method 2 only:

  • altpolicy::Symbol: Which alternative policy is being used
  • obs::Matrix{S <: Real}: matrix of forecasted observables

where S<:AbstractFloat.

Keyword Arguments

  • cond_type::Symbol: one of :none, :semi, or :full, used to determine how many periods to forecast ahead. If cond_type in [:semi, :full], the forecast horizon is reduced by the number of periods of conditional data. Defaults to :none.

  • enforce_zlb::Bool: whether to enforce the zero lower bound. Defaults to false.

  • shocks::Matrix{S}: matrix of size nshocks x shock_horizon of shock innovations under which to forecast. If shock_horizon > horizon, the extra periods of shocks will be ignored; if shock_horizon < horizon, zeros will be filled in for the shocks hitting the remaining forecasted periods.

  • draw_shocks::Bool: if isempty(shocks), indicates whether to draw shocks according to:

    1. If forecast_tdist_shocks(m), draw horizons many shocks from a Distributions.TDist(forecast_tdist_df_val(m))
    2. Otherwise, draw horizons many shocks from a DegenerateMvNormal(zeros(nshocks), sqrt(system[:QQ]))

    or to set shocks to a nshocks x horizon matrix of zeros. Defaults to false. If shocks is provided as a keyword argument, this flag has no effect.

Method 2 only:

  • set_zlb_regime_vals::Function: user-provided function that adds additional regimes to regime-switching parameters if not enough regimes exist to impose the ZLB as a temporary alternative policy. Defaults to identity, and nothing will happen if this is the case.

  • tol::{<: Real}: Tolerance for the smallest permissible value for the nominal interest rate. Defaults to -1e-14.

Outputs

  • states::Matrix{S}: matrix of size nstates x horizon of forecasted states
  • obs::Matrix{S}: matrix of size nobs x horizon of forecasted observables
  • pseudo::Matrix{S}: matrix of size npseudo x horizon of forecasted pseudo-observables
  • shocks::Matrix{S}: matrix of size nshocks x horizon of shock innovations
source
DSGE.forecast_scenarioMethod
forecast_scenario(m, scen::Scenario; verbose = :low)

Simulate all draws of scen using the modal parameters of the model m. This function returns a Dict{Symbol, Array{Float64}.

source
DSGE.histforecastFunction
histforecast(var, hist, forecast;
    start_date = hist.means[1, :date], end_date = forecast.means[end, :date],
    names = Dict{Symbol, String}(), colors = Dict{Symbol, Any}(),
    alphas = Dict{Symbol, Float64}(), styles = Dict{Symbol, Symbol}(),
    bands_pcts = union(which_density_bands(hist, uniquify = true),
                       which_density_bands(forecast, uniquify = true)),
    bands_style = :fan, label_bands = false, transparent_bands = true,
    tick_size = 2)

User recipe called by plot_history_and_forecast.

Inputs

  • var::Symbol: e.g. obs_gdp
  • hist::MeansBands
  • forecast::MeansBands

Keyword Arguments

  • start_date::Date
  • end_date::Date
  • names::Dict{Symbol, String}: maps keys [:hist, :forecast, :bands] to labels. If a key is missing from names, a default value will be used
  • colors::Dict{Symbol, Any}: maps keys [:hist, :forecast, :bands] to colors
  • alphas::Dict{Symbol, Float64}: maps keys [:hist, :forecast, :bands] to transparency values (between 0.0 and 1.0)
  • styles::Dict{Symbol, Symbol}: maps keys [:hist, :forecast, :bands] to linestyles
  • bands_pcts::Vector{String}: which bands percentiles to plot
  • bands_style::Symbol: either :fan or :line
  • label_bands::Bool
  • transparent_bands::Bool
  • tick_size::Int: x-axis (time) tick size in units of years
  • add_new_model::Bool: Adding history from another model?
  • new_data::Array{Flaot64,1}: The new history to plot

Additionally, all Plots attributes (see docs.juliaplots.org/latest/attributes) are supported as keyword arguments.

source
DSGE.histforecast_vectorFunction
histforecast_vector(var, hist, forecast;
    start_date::Date = hf.args[2][1].means[1, :date],
    end_date::Date = hf.args[3][1].means[end, :date],
    names = Dict{Symbol, Vector{String}}(),
    colors = Dict{Symbol, Any}(),
    alphas = Dict{Symbol, Float64}(),
    styles = Dict{Symbol, Symbol}(),
    add_new_model::Bool = false, new_data::Array{Float64,1} = [],
    bands_pcts = union(which_density_bands(hf.args[2][1], uniquify = true),
        which_density_bands(hf.args[3][1], uniquify = true)),
    bands_style = :fan,
    label_bands = false,
    transparent_bands = true,
    add_trendline::Bool = false, trend_vals::Vector{Float64} = [1.0],
    tick_size = 2)

User recipe called by plot_history_and_forecast.

Inputs

  • var::Symbol: e.g. obs_gdp
  • hist::Vector{MeansBands}
  • forecast::Vector{MeansBands}

Keyword Arguments

  • start_date::Date
  • end_date::Date
  • names::Dict{Symbol, String}: maps keys [:hist, :forecast, :bands] to labels. If a key is missing from names, a default value will be used
  • colors::Dict{Symbol, Any}: maps keys [:hist, :forecast, :bands] to colors
  • alphas::Dict{Symbol, Float64}: maps keys [:hist, :forecast, :bands] to transparency values (between 0.0 and 1.0)
  • styles::Dict{Symbol, Symbol}: maps keys [:hist, :forecast, :bands] to linestyles
  • bands_pcts::Vector{String}: which bands percentiles to plot
  • bands_style::Symbol: either :fan or :line
  • label_bands::Bool
  • transparent_bands::Bool
  • tick_size::Int: x-axis (time) tick size in units of years
  • add_new_model::Bool: Adding history from another model?
  • new_data::Array{Float64,1}: The new history to plot
  • add_trendline::Bool: Whether add a trendline or not
  • trend_vals::Vector{Float64}: The values to multiply last historical value by for the trendline (should start with 0.0 to keep last historical value the same)
  • plot_all_histories::Bool: Whether to plot histories for each model (generally used for pseudo observables like Natural Rate.

Additionally, all Plots attributes (see docs.juliaplots.org/latest/attributes) are supported as keyword arguments.

source
DSGE.plot_history_and_forecastMethod
plot_history_and_forecast(m, var, class, input_type, cond_type;
    title = "", plot_handle = plot(), kwargs...)

plot_history_and_forecast(m, vars, class, input_type, cond_type;
    forecast_string = "", use_bdd = :unbdd,
    plotroot = figurespath(m, "forecast"), titles = [],
    plot_handles = fill(plot(), length(vars)), verbose = :low,
    add_trendline::Bool = false, trend_vals::Vector{Float64} = [1.0],
    kwargs...)

plot_history_and_forecast(ms::Vector, vars::Vector{Symbol}, class::Symbol,
    input_type::Symbol, cond_type::Symbol;
    forecast_string::String = "",
    use_bdd::Symbol = :bdd,
    modal_line::Bool = false, untrans::Bool = false,
    fourquarter::Bool = false,
    plotroot::String = figurespath(ms[1], "forecast"),
    titles::Vector{String} = String[],
    plot_handles::Vector{Plots.Plot} = Plots.Plot[plot() for i = 1:length(vars)],
    verbose::Symbol = :low, names = Dict{Symbol, Vector{String}}(),
    start_date::Date = Date(2019,3,31), end_date::Date = iterate_quarters(date_forecast_start(m), 12),
    add_new_model::Bool = false, new_model::AbstractDSGEModel = Model1002(),
    new_cond_type::Symbol = :full, outfile_end::String =  "",
    new_forecast_string::String = forecast_string,
    add_trendline::Bool = false, trend_vals::Vector{Float64} = [1.0],
    forecast_strings::Vector{String} = repeat([forecast_string], length(ms)),
    kwargs...)

Plot history and forecast for var or vars. If these correspond to a full-distribution forecast, you can specify the bands_style and bands_pcts.

Inputs

  • m::AbstractDSGEModel
  • var::Symbol or vars::Vector{Symbol}: variable(s) to be plotted, e.g. :obs_gdp or [:obs_gdp, :obs_nominalrate]
  • class::Symbol
  • input_type::Symbol
  • cond_type::Symbol

Keyword Arguments

  • forecast_string::String = ""
  • use_bdd::Symbol = :unbdd: specifies which combination of means and bands to use a. :bdd -> bounded bands and bounded means (from :bddforecastobs, etc.) b. :bdd_and_unbdd -> bounded bands (from :bddforecastobs, etc.) and unbounded means (from :forecastobs, etc.) c. :unbdd -> unbounded bands and unbounded means (from :forecastobs, etc.)
  • modal_line::Bool = false: if true, the modal line is plotted instead of the mean.
  • untrans::Bool = false: whether to plot untransformed (model units) history and forecast
  • fourquarter::Bool = false: whether to plot four-quarter history and forecast
  • plotroot::String = figurespath(m, "forecast"): if nonempty, plots will be saved in that directory
  • title::String = "" or titles::Vector{String} = []
  • plot_handle::Plot or plot_handles::Vector{Plot}: existing plot(s) on which to overlay new forecast plot(s)
  • verbose::Symbol = :low
  • names::Dict{Symbol, Vector{String}}: Legend names for lines
  • start_date::Date = Date(2019,3,31)`: Date plot should start at
  • end_date::Date = iteratequarters(dateforecast_start(m), 12)`: Date plot should end at
  • add_new_model::Bool = false: Whether history should be plotted using a new model
  • new_model::AbstractDSGEModel = Model1002(): The new model to use to plot history
  • new_cond_type::Symbol = :full: Cond_type of new model
  • outfile_end::String = "": String to append to the end of file name of plot
  • new_forecast_string::String = forecast_string: Forecast string for new model
  • add_trendline::Bool: Whether add a trendline or not
  • trend_vals::Vector{Float64}: The values to multiply last historical value by for the trendline (should start with 0.0 to keep last historical value the same)
  • forecast_strings::Vector{String}: Forecast string to use for each give model.

See ?histforecast or ?histforecast_vector' for additional keyword arguments, all of which can be passed intoplothistoryand_forecast`.

Output

  • p::Plot or plots::OrderedDict{Symbol, Plot}
source
DSGE.smoothMethod
smooth(m, df, system, s_0, P_0; cond_type = :none, draw_states = true,
    include_presample = false)

Computes and returns the smoothed values of states and shocks for the system system.

Inputs

  • m::AbstractDSGEModel: model object
  • df::DataFrame: data for observables. This should include the conditional period if cond_type in [:semi, :full]
  • system::System: System object representing the state-space system
  • s_0::Vector{S}: optional initial state vector
  • P_0::Matrix{S}: optional initial state covariance matrix

Keyword Arguments

  • cond_type: conditional case. See forecast_one for documentation of all cond_type options
  • draw_states: if using a simulation smoother (i.e. forecast_smoother(m) in [:carter_kohn, :durbin_koopman]), indicates whether to draw smoothed states from the distribution N(z_{t|T}, P_{t|T}) or to use the mean z_{t|T}. Defaults to false. If not using a simulation smoother, this flag has no effect (though the user will be warned if draw_states = true)
  • include_presample::Bool: indicates whether to include presample periods in the returned matrices. Defaults to false.
  • in_sample::Bool: indicates whether or not to discard out of sample rows in df_to_matrix call.

Outputs

  • states::Matrix{S}: array of size nstates x hist_periods of smoothed states (not including the presample)
  • shocks::Matrix{S}: array of size nshocks x hist_nperiods of smoothed shocks
  • pseudo::Matrix{S}: matrix of size npseudo x hist_periods of pseudo-observables computed from the smoothed states
  • initial_states::Vector{S}: vector of length nstates of the smoothed states in the last presample period. This is used as the initial state for computing the deterministic trend

Notes

states and shocks are returned from the smoother specified by forecast_smoother(m), which defaults to :durbin_koopman. This can be overridden by calling

m <= Setting(:forecast_smoother, :koopman_smoother))

before calling smooth.

source
DSGE.deterministic_trendsMethod
deterministic_trends(m, system, z0)

deterministic_trends(system, z0, nperiods, start_index, end_index)

deterministic_trends(m, system, z0, start_date, end_date)

deterministic_trends(m, system, z0, nperiods, start_index, end_index,
    regime_inds, regimes, cond_type)

deterministic_trends(m, system, old_system, z0, start_date, end_date)

deterministic_trends(m, system, old_system, z0, nperiods, start_index, end_index,
    regime_inds, regimes, cond_type)

Compute deterministic trend values of states, observables, and pseudo-observables, given a model object and system matrices. The deterministic trend for a single draw is simply the series that would be obtained by iterating the state-space system forward, beginning from a state vector z0 in the last presample period.

Inputs

  • m::AbstractDSGEModel: model object
  • system::System{S} or RegimeSwitchingSystem: state-space system matrices
  • z0::Vector{S}: initial state vector
  • start_date::Date: initial date for deterministic trends
  • end_date::Date: final date for deterministic trends
  • regime_inds::Vector{UnitRange}: indices of the data corresponding to each regime.
  • regimes::UnitRange: which regimes are involved in the date range for which we want to compute the deterministic trends
  • old_system::RegimeSwitchingSystem`: state-space system matrices for old model

where S<:AbstractFloat.

Outputs

  • states::Matrix{S}: matrix of size nstates x nperiods of state steady-state values
  • obs::Matrix{S}: matrix of size nobs x nperiods of observable steady-state values
  • pseudo::Matrix{S}: matrix of size npseudo x nperiods of pseudo-observable steady-state values

where nperiods is the number of quarters between date_shockdec_start(m) and date_shockdec_end(m), inclusive.

source
DSGE.shock_decompositionsMethod
shock_decompositions(m, system, histshocks)

shock_decompositions(system, forecast_horizons, histshocks, start_index,
    end_index)

shock_decompositions(m, system, histshocks, start_date, end_date)

shock_decompositions(m, system, forecast_horizons, histshocks, start_index,
    end_index, regime_inds, cond_type)

shock_decompositions(m, system, old_system,
                              histshocks, old_histshocks, start_date, end_date,
                              cond_type; full_shock_decomp)

shock_decompositions(m, system, old_system,
                              forecast_horizons, histshocks, old_histshocks,
                              start_index, end_index,
                              regime_inds, cond_type; full_shock_decomp)

Inputs

  • system::System{S} or RegimeSwitchingSystem{S}: state-space system matrices
  • forecast_horizons::Int: number of periods ahead to forecast
  • histshocks::Matrix{S}: matrix of size nshocks x hist_periods of historical smoothed shocks
  • start_index::Int: first index from which to return computed shock decompositions
  • end_index::Int: last index for which to return computed shock decompositions
  • start_date::Date: initial date for deterministic trends
  • end_date::Date: final date for deterministic trends
  • regime_inds::Vector{UnitRange}: indices of the data corresponding to each regime.
  • old_system::RegimeSwitchingSystem{S}: state-space system matrices for the old model
  • old_histshocks::Matrix{S}: matrix of size nshocks x hist_periods of historical smoothed shocks using the old system
  • full_shock_decomp::Bool: If true for old_system case, return difference in shockdecs between all shocks for the new and old systems. Else, return old system's shocks with 1 cumulative AIT shock added on.

where S<:AbstractFloat.

Outputs

  • states::Array{S, 3}: matrix of size nstates x nperiods x nshocks of state shock decompositions
  • obs::Array{S, 3}: matrix of size nobs x nperiods x nshocks of observable shock decompositions
  • pseudo::Array{S, 3}: matrix of size npseudo x nperiods x nshocks of pseudo-observable shock decompositions

where nperiods =endindex - startindex + 1`.

source
DSGE.trendsMethod
trends(system::System{S}) where {S<:AbstractFloat}
trends(system::RegimeSwitchingSystem{S}) where {S<:AbstractFloat}
trends(m::AbstractDSGEModel, system::RegimeSwitchingSystem{S},
                start_date::Dates.Date = date_presample_start(m),
                end_date::Dates.Date = prev_quarter(date_forecast_start(m)),
                cond_type::Symbol) where {S<:AbstractFloat}

Compute trend (steady-state) states, observables, and pseudo-observables. The trend is used for plotting shock decompositions. The first method applies to non-regime switching systems, the second to regime-switching systems that do not involve time variation in the CCC or DD vectors, and the third to regime-switching systems with time variation in CCC or DD.

source
DSGE.get_forecast_filenameMethod
get_forecast_filename(m, input_type, cond_type, output_var;
    pathfcn = rawpath, forecast_string = "", fileformat = :jld2)

get_forecast_filename(directory, filestring_base, input_type, cond_type,
    output_var; forecast_string = "", fileformat = :jld2)

Notes

  • If input_type == :subset, then the forecast_string is also appended to the

filenames. If in this case forecast_string is empty, get_forecast_filename throws an error.

  • In the second method, directory should be a string of the form "$saveroot/m990/ss2/forecast/raw/". (Note that a pathfcn is therefore not required.) filestring_base should be equivalent to the result of filestring_base(m).
source
DSGE.get_forecast_input_fileMethod
get_forecast_input_file(m, input_type)

Compute the appropriate forecast input filenames for model m and forecast input type input_type.

The default input files for each input_type can be overriden by adding entries to the Dict{Symbol, String} returned from forecast_input_file_overrides(m). For example:

overrides = forecast_input_file_overrides(m)
overrides[:mode] = "path/to/input/file.h5"
source
DSGE.get_forecast_output_filesMethod
get_forecast_output_files(m, input_type, cond_type, output_vars;
    forecast_string = "", fileformat = :jld2)

get_forecast_output_files(directory, filestring_base, input_type, cond_type,
    output_vars; forecast_string = "", fileformat = :jld2)

Compute the appropriate output filenames for model m, forecast input type input_type, and conditional type cond_type, for each output variable in output_vars. Returns a dictionary of file names with one entry for each output_var.

Arguments

  • See forecast_one for descriptions of other non-keyword arguments.

Keyword Arguments

  • forecast_string::String: subset identifier for when input_type = :subset
  • fileformat::Symbol: file extension, without a period. Defaults to :jld2, though :h5 is another common option.

Notes

See get_forecast_filename for more information.

source
DSGE.get_meansbands_input_fileMethod
get_meansbands_input_file(m, input_type, cond_type, output_var;
    forecast_string = "", fileformat = :jld2)
get_meansbands_input_file(directory, filestring_base, input_type, cond_type, output_var;
    forecast_string = "", fileformat = :jld2)

Returns a dictionary of raw forecast output files to read in to compute means and bands.

Inputs

Method 1:

  • m::AbstractDSGEModel

Method 2:

  • directory::String: directory location of input files to read
  • filestring_base::Vector{String}: a vector of strings to be added as a suffix. These usually come from model settings for which print = true. It should not include entries for cond_type and input_type (these will be added automatically).

Both methods:

  • input_type::Symbol: See ?forecast_one
  • cond_type::Symbol: See ?forecast_one
  • output_var::Symbol: See ?forecast_one
  • forecast_string::String: See ?forecast_one
  • fileformat: file extension of saved files
source
DSGE.get_meansbands_output_fileMethod
get_meansbands_output_file(m, input_type, cond_type, output_var;
    forecast_string = "", fileformat = :jld2)
get_meansbands_output_file(directory, filestring_base, input_type, cond_type, output_var;
    forecast_string = "", fileformat = :jld2)

Returns a dictionary of raw forecast output files in which to save computed means and bands.

Inputs

Method 1:

  • m::AbstractDSGEModel: Model object

Method 2:

  • directory::String: directory location of input files to read
  • filestring_base::Vector{String}: a vector of strings to be added as a suffix. These usually come from model settings for which print=true. It should not include entries for cond_type and input_type (these will be added automatically).

Both methods:

  • input_type::Symbol: See ?forecast_one
  • cond_type::Symbol: See ?forecast_one
  • output_vars::Symbol: See ?forecast_one
  • forecast_string::String: See ?forecast_one
  • fileformat: file extension of saved files
source
DSGE.get_scenario_filenameMethod
get_scenario_filename(m, scen::AbstractScenario, output_var;
    pathfcn = rawpath, fileformat = :jld2, directory = "")

Get scenario file name of the form pathfcn(m, "scenarios", output_var * filestring * string(fileformat)). If directory is provided (nonempty), then the same file name in that directory will be returned instead.

source
DSGE.get_scenario_mb_input_fileMethod
get_scenario_mb_input_file(m, scen::AbstractScenario, output_var)

Call get_scenario_filename while replacing forecastut and forecast4q in output_var with forecast.

source
DSGE.get_scenario_mb_output_fileMethod
get_scenario_mb_output_file(m, scen::AbstractScenario, output_var;
    directory = "")

Call get_scenario_filename while tacking on "mb" to the front of the base file name.

source
DSGE.get_scenario_output_filesMethod
get_scenario_output_files(m, scen::SingleScenario, output_vars)

Return a Dict{Symbol, String} mapping output_vars to the raw simulated scenario outputs for scen.

source
DSGE.plot_scenarioMethod
plot_scenario(m, var, class, scen; title = "", kwargs...)

plot_scenario(m, vars, class, scen; untrans = false, fourquarter = false,
    plotroot = figurespath(m, "scenarios"), titles = [], tick_size = 1,
    kwargs...)

Plot var or vars in deviations from baseline for the alternative scenario specified by key and vint.

Inputs

  • var::Symbol or vars::Vector{Symbol}: variable(s) to be plotted, e.g. :obs_gdp or [:obs_gdp, :obs_nominalrate]
  • class::Symbol
  • scen::AbstractScenario: scenario

Keyword Arguments

  • untrans::Bool: whether to plot untransformed (model units) forecast
  • fourquarter::Bool: whether to plot four-quarter forecast
  • plotroot::String: if nonempty, plots will be saved in that directory
  • title::String or titles::Vector{String}
  • tick_size::Int: x-axis (time) tick size in units of years
  • legend

See ?histforecast for additional keyword arguments, all of which can be passed into plot_scenario.

Output

  • p::Plot or plots::OrderedDict{Symbol, Plot}
source
DSGE.read_bdd_and_unbdd_mbMethod
read_bdd_and_unbdd_mb(bdd_fn::String, unbdd_fn::String; modal_line::Bool = false)

Read in the bounded and unbounded forecast MeansBands from bdd_fn and unbdd_fn. Create and return a MeansBands with the unbounded means and bounded bands. If modal_line is true, then the unbdd_fn is known to load in a modal forecast but should be treated as having the same input_type as the bounded forecast.

source
DSGE.read_mbMethod
read_mb(fn::String)

read_mb(fn1::String, fn2::String)

read_mb(m, input_type, cond_type, output_var; forecast_string = "",
    use_bdd = :unbdd, modal_line = false, directory = workpath(m, "forecast"))

Read in a MeansBands object saved in fn, or use the model object m to determine the file location.

The second method construct a MeansBands object with means from the modal object and bands, where fn1 is the file location of the bands and fn2 is the file location of the means.

If bdd_and_unbdd, then output_var must be either :forecast or :forecast4q. Then this function calls read_bdd_and_unbdd to return a MeansBands with unbounded means and bounded bands. If modal line is set to true, then the modal mean rather than the full-distribution mean is returned.

source
DSGE.read_scenario_mbMethod
read_scenario_mb(m, scen::AbstractScenario, output_var; directory = "")

Read in an alternative scenario MeansBands object.

source
DSGE.read_scenario_outputMethod
read_scenario_output(m, scen::SingleScenario, class, product, var_name)

read_scenario_output(m, agg::ScenarioAggregate, class, product, var_name)

Given either scen or agg, read in and return all draws of and the appropriate reverse transform for var_name.

The third function that takes in two models is used for when we have scenarios from two different models.

source
DSGE.write_means_tables_shockdecMethod
write_means_tables_shockdec(m, input_type, cond_type, class;
    forecast_string = "",
    read_dirname = workpath(m, "forecast"),
    write_dirname = tablespath(m, "forecast"),
    kwargs...)

write_means_tables_shockdec(write_dirname, filestring_base, mb_shockdec,
    mb_trend, mb_dettrend, mb_hist, mb_forecast; tablevars = get_variables(mb),
    columnvars = get_shocks(mb), groups = [])

Inputs

Method 1 only:

  • m::AbstractDSGEModel
  • input_type::Symbol
  • cond_type::Symbol
  • class::Symbol

Method 2 only:

  • write_dirname::String: directory to which tables are saved
  • filestring_base::Vector{String}: the result of filestring_base(m), typically ["vint=yymmdd"]`
  • mb_shockdec::MeansBands
  • mb_trend::MeansBands
  • mb_dettrend::MeansBands
  • mb_hist::MeansBands: optional
  • mb_forecast::MeansBands: optional

Keyword Arguments

  • tablevars::Vector{Symbol}: which series to write tables for
  • columnvars::Vector{Symbol}: which shocks to include as columns in the tables
  • groups::Vector{ShockGroup}: if provided, shocks will be grouped accordingly

Method 1 only:

  • forecast_string::String
  • use_bdd::Symbol: whether to use unbounded means and bounded bands. Applies only for class(output_var) in [:forecast, :forecast4q]
  • read_dirname::String: directory from which MeansBands are read in
  • write_dirname::String: directory to which tables are saved
source
DSGE.write_meansbands_tables_allMethod
write_meansbands_tables_all(m, input_type, cond_type, output_vars;
    forecast_string = "", dirname = tablespath(m, "forecast"),
    vars = [], shocks = [], shock_groups = [])

Write all output_vars corresponding to model m to tables in dirname.

Inputs

  • m::AbstractDSGEModel
  • input_type::Symbol: See ?forecast_one
  • cond_type::Symbol: See ?forecast_one
  • output_vars::Symbol: See ?forecast_one

Keyword Arguments

  • forecast_string::String: See ?forecast_one
  • vars::Vector{Symbol}: Vector of economic variables for which to print output_vars to tables. If omitted, all shocks will be printed.
  • shocks::Vector{Symbol}: Vector of shocks to print if output_vars contains a shock decomposition. If omitted, all shocks will be printed.
  • shock_groups::Vector{ShockGroup}: if provided, shocks will be grouped accordingly in shockdec tables
source
DSGE.write_meansbands_tables_timeseriesMethod
write_meansbands_tables_timeseries(m, input_type, cond_type, output_var;
    forecast_string = "", bdd_and_unbdd = false,
    read_dirname = workpath(m, "forecast"),
    write_dirname = tablespath(m, "forecast"), kwargs...)

write_meansbands_tables_timeseries(dirname, filestring_base, mb;
    tablevars = get_variables(mb))

Inputs

Method 1 only:

  • m::AbstractDSGEModel
  • input_type::Symbol
  • cond_type::Symbol
  • output_var::Symbol: class(output_var) must be one of [:hist, :histut, :hist4q, :forecast, :forecastut, :forecast4q, :bddforecast, :bddforecastut, :bddforecast4q, :trend, :dettrend, :histforecast, :histforecastut, :histforecast4q]
  • read_dirname::String: directory to which meansbands objects are read from

Method 2 only:

  • read_dirname::String: directory from which MeansBands are read in
  • write_dirname::String: directory to which tables are saved
  • filestring_base::Vector{String}: the result of filestring_base(m), typically ["vint=yymmdd"]`

Keyword Arguments

  • tablevars::Vector{Symbol}: which series to write tables for

Method 1 only:

  • forecast_string::String
  • bdd_and_unbdd::Bool: whether to use unbounded means and bounded bands. Applies only for class(output_var) in [:forecast, :forecast4q]
  • dirname::String: directory to which tables are saved
source
DSGE.write_ref_trialMethod
write_ref_trial(trial, trial_name)

Write a reference trial to a JLD file, to act as the standard that new trials are benchmarked against.

Arguments

  • trial::BenchmarkTools.Trial: The trial object that is being written.
  • trial_name::String: The name of the trial being written.
source
DSGE.decomposition_forecastMethod
decomposition_forecast(m, df, params, cond_type, keep_startdate, keep_enddate, shockdec_splitdate;
    outputs = [:forecast, :shockdec], check = false)

Equivalent of forecast_one_draw for forecast decomposition. keep_startdate = date_forecast_start(m_new) corresponds to time T+1, keep_enddate = date_forecast_end(m_old) to time T+H, and shockdec_splitdate = date_mainsample_end(m_old) to time T-k.

Returns out::Dict{Symbol, Array{Float64}}, which has keys determined as follows:

  • If :forecast in outputs or check = true:

    • :forecast<class>
  • If :shockdec in outputs:

    • :trend<class>
    • :dettrend<class>
    • :data<class>: like a shockdec, but only applying smoothed shocks up to shockdec_splitdate
    • :news<class>: like a shockdec, but only applying smoothed shocks after shockdec_splitdate
source
DSGE.decomposition_periodsMethod
decomposition_periods(m_new, m_old, df_new, df_old, cond_new, cond_old)

Returns T, k, and H, where:

  • New model has T periods of data
  • Old model has T-k periods of data
  • Old and new models both forecast up to T+H
source
DSGE.forecast_one_drawMethod
forecast_one_draw(m, input_type, cond_type, output_vars, params, df;
    verbose = :low, only_filter = false)

Compute output_vars for a single parameter draw, params. Called by forecast_one.

Inputs

  • m::AbstractDSGEModel{Float64}: model object
  • input_type::Symbol: See ?forecast_one.
  • cond_type::Symbol: See ?forecast_one.
  • output_vars::Vector{Symbol}: vector of desired output variables. See Outputs section
  • params::Vector{Float64}: parameter vector
  • df::DataFrame: historical data.
  • verbose::Symbol: desired frequency of function progress messages printed to standard out. One of :none, :low, or :high.

Output

  • forecast_outputs::Dict{Symbol, Array{Float64}}: dictionary of forecast outputs. Keys are output_vars, which is some subset of:
  - `:histstates`: `Matrix{Float64}` of smoothed historical states
  - `:histobs`: `Matrix{Float64}` of smoothed historical data
  - `:histpseudo`: `Matrix{Float64}` of smoothed historical
    pseudo-observables
  - `:histshocks`: `Matrix{Float64}` of smoothed historical shocks
  - `:forecaststates`: `Matrix{Float64}` of forecasted states
  - `:forecastobs`: `Matrix{Float64}` of forecasted observables
  - `:forecastpseudo`: `Matrix{Float64}` of forecasted pseudo-observables
  - `:forecastshocks`: `Matrix{Float64}` of forecasted shocks
  - `:bddforecaststates`, `:bddforecastobs`, `:bddforecastpseudo`, and
    `:bddforecastshocks`: `Matrix{Float64}`s of forecasts where we enforce
    the zero lower bound to be `forecast_zlb_value(m)`
  - `:shockdecstates`: `Array{Float64, 3}` of state shock decompositions
  - `:shockdecobs`: `Array{Float64, 3}` of observable shock decompositions
  - `:shockdecpseudo`: `Array{Float64, 3}` of pseudo-observable shock
    decompositions
  - `:dettrendstates`: `Matrix{Float64}` of state deterministic trends
  - `:dettrendobs`: `Matrix{Float64}` of observable deterministic trends
  - `:dettrendpseudo`: `Matrix{Float64}` of pseudo-observable deterministic
    trends
  - `:trendstates`: `Vector{Float64}` of state trends, i.e. the `CCC` vector
  - `:trendobs`: `Vector{Float64}` of observable trends, i.e. the `DD` vector
  - `:trendpseudo`: `Vector{Float64}` of pseudo-observable trends, i.e. the
    `DD_pseudo` vector
  - `:irfstates`: `Array{Float64, 3}` of state impulse responses
  - `:irfobs`: `Array{Float64, 3}` of observable impulse responses
  - `:irfpseudo`: `Array{Float64, 3}` of pseudo-observable impulse responses
source
DSGE.prepare_forecast_inputs!Method
prepare_forecast_inputs!(m, input_type, cond_type, output_vars;
    df = DataFrame(), verbose = :none)

Add required outputs using add_requisite_output_vars and load data if necessary.

Inputs

  • m::AbstractDSGEModel: model object
  • input_type::Symbol: See ?forecast_one.
  • cond_type::Symbol: See ?forecast_one.
  • output_vars::Vector{Symbol}: vector of desired output variables. See ?forecast_one_draw

Keyword Arguments

  • df::DataFrame: historical data. If cond_type in [:semi, :full], then the final row of df should be the period containing conditional data. If not provided, then df will be loaded using load_data with the appropriate cond_type
  • only_filter::Bool: do not run the smoother and only run the filter. This limits the number of output variables which can be calculated.
  • verbose::Symbol: desired frequency of function progress messages printed to standard out. One of :none, :low, or :high

Outputs

  • output_vars
  • df
source
DSGE.default_update_regime_eqcond_info!Method
default_update_regime_eqcond_info!(m::AbstractDSGEModel, eqcond_dict::AbstractDict{Int64, EqcondEntry},
                                   zlb_start_regime::Int64, liftoff_regime::Int64, liftoff_policy::AltPolicy)

is the default method for updating the setting regime_eqcond_info during the automatic endogenous ZLB enforcement as a temporary policy.

The input eqcond_dict should be a dictionary whose historical/conditional horizon regimes should already be set (and won't be affected by this function. The regimes which will be altered are those for which the temporary ZLB will apply (according to zlb_start_regime and liftoff_regime), and if there is time-varying credibility, forecast regimes past the ZLB.

In general, other implementations of the update_regime_eqcond_info! function should not change the EqcondEntry during historical/conditional horizon regimes, but any other regimes in the forecast horizon should be/can be set.

source
DSGE.forecast_scenario_drawMethod
forecast_scenario_draw(m, scen::Scenario, system, draw_index)

Filter shocks and use them to forecast the draw_indexth draw of scen.

source
DSGE.write_scenario_forecastsMethod
write_scenario_forecasts(m, scenario_output_files, forecast_output;
    verbose = :low)

Write scenario outputs in forecast_output to values(scenario_output_files).

source
DSGE.shock_decompositions_sequenceMethod
shock_decompositions_sequence(m, system, histshocks, n_back, back_shocks)

shock_decompositions_sequence(system, forecast_horizons, histshocks, start_index,
    end_index, n_back, back_shocks)

shock_decompositions_sequence(m, system, histshocks, start_date, end_date, n_back, back_shocks)

shock_decompositions_sequence(m, system, forecast_horizons, histshocks, start_index,
    end_index, regime_inds, cond_type, n_back, back_shocks)

shock_decompositions_sequence(m, system, old_system,
                              histshocks, old_histshocks, start_date, end_date,
                              cond_type; full_shock_decomp, n_back, back_shocks)

shock_decompositions_sequence(m, system, old_system,
                              forecast_horizons, histshocks, old_histshocks,
                              start_index, end_index,
                              regime_inds, cond_type; full_shock_decomp, n_back, back_shocks)

Inputs

  • system::System{S} or RegimeSwitchingSystem{S}: state-space system matrices
  • forecast_horizons::Int: number of periods ahead to forecast
  • histshocks::Matrix{S}: matrix of size nshocks x hist_periods of historical smoothed shocks
  • start_index::Int: first index from which to return computed shock decompositions
  • end_index::Int: last index for which to return computed shock decompositions
  • start_date::Date: initial date for deterministic trends
  • end_date::Date: final date for deterministic trends
  • regime_inds::Vector{UnitRange}: indices of the data corresponding to each regime.
  • old_system::RegimeSwitchingSystem{S}: state-space system matrices for the old model
  • old_histshocks::Matrix{S}: matrix of size nshocks x hist_periods of historical smoothed shocks using the old system
  • full_shock_decomp::Bool: If true for old_system case, return difference in shockdecs between all shocks for the new and old systems. Else, return old system's shocks with 1 cumulative AIT shock added on.

where S<:AbstractFloat.

Outputs

  • states::Array{S, 4}: matrix of size nhistperiods x nstates x nperiods x nshocks of state shock decompositions
  • obs::Array{S, 4}: matrix of size nhistperiods x nobs x nperiods x nshocks of observable shock decompositions
  • pseudo::Array{S, 4}: matrix of size nhistperiods x npseudo x nperiods x nshocks of pseudo-observable shock decompositions

where nperiods =endindex - startindex + 1andnhistperiods=size(histshocks,2)`

source
DSGE.combine_raw_forecast_output_and_metadataMethod
combine_raw_forecast_output_and_metadata(m, forecast_output_files; verbose = :low)

Writes the raw forecast output data (arr) saved in the temporary h5 file to the jld2 file containing the rest of the forecast metadata. The intermediary h5 step exists because jld2 does not support chunked memory assignment in the same way that jld and h5 permitted previously.

source
DSGE.count_scenario_draws!Method
count_scenario_draws!(m, scen::Scenario)

Return the number of draws for scen, determined using get_scenario_input_file(m, scen), and update the n_draws field of scen with this count.

source
DSGE.get_scenario_mb_metadataMethod
get_scenario_mb_metadata(m, scen::SingleScenario, output_var)

get_scenario_mb_metadata(m, agg::ScenarioAggregate, output_var)

Return the MeansBands metadata dictionary for scen.

source
DSGE.load_scenario_targets!Method
load_scenario_targets!(m, scen::Scenario, draw_index)

Add the targets from the draw_indexth draw of the raw scenario targets to scen.targets.

source
DSGE.read_forecast_metadataMethod
read_forecast_metadata(file::JLDFile)

Read metadata from forecast output files. This includes dictionaries mapping dates, as well as state, observable, pseudo-observable, and shock names, to their respective indices in the saved forecast output array. Depending on the output_var, the saved dictionaries might include:

  • date_indices::Dict{Date, Int}: not saved for IRFs
  • state_indices::Dict{Symbol, Int}
  • observable_indices::Dict{Symbol, Int}
  • pseudoobservable_indices::Dict{Symbol, Int}
  • shock_indices::Dict{Symbol, Int}
  • state_revtransforms::Dict{Symbol, Symbol}: states are not transformed, so all values are :identity
  • observable_revtransforms::Dict{Symbol, Symbol}
  • pseudoobservable_revtransforms::Dict{Symbol, Symbol}
  • shock_revtransforms::Dict{Symbol, Symbol}: shocks are not transformed, so all values are :identity
source
DSGE.read_forecast_seriesMethod
read_forecast_series(filepath, product, var_ind)
read_forecast_series(filepath, var_ind, shock_ind)

Read only the forecast output for a particular variable (e.g. for a particular observable) and possibly a particular shock. Result should be a matrix of size ndraws x nperiods.

source
DSGE.read_regime_switching_trendMethod
read_regime_switching_trend(filepath, var_ind)

Read only the trend output for a particular variable (e.g. for a particular observable). Result should be a matrix of size ndraws × n_regimes or ndraws × n_data_periods, depending on whether the state space system was time-varying or not.

source
DSGE.write_forecast_blockMethod
write_forecast_block(file::JLDFile, arr::Array, block_number::Int,
    block_inds::AbstractRange{Int64})

Writes arr to the subarray of file indicated by block_inds.

source
DSGE.write_forecast_metadataMethod
write_forecast_metadata(m::AbstractDSGEModel, file::JLDFile, var::Symbol)

Write metadata about the saved forecast output var to filepath.

Specifically, we save dictionaries mapping dates, as well as state, observable, pseudo-observable, and shock names, to their respective indices in the saved forecast output array. The saved dictionaries include:

  • date_indices::Dict{Date, Int}: saved for all forecast outputs except IRFs
  • state_names::Dict{Symbol, Int}: saved for var in [:histstates, :forecaststates, :shockdecstates]
  • observable_names::Dict{Symbol, Int}: saved for var in [:forecastobs, :shockdecobs]
  • observable_revtransforms::Dict{Symbol, Symbol}: saved identifiers for reverse transforms used for observables
  • pseudoobservable_names::Dict{Symbol, Int}: saved for var in [:histpseudo, :forecastpseudo, :shockdecpseudo]
  • pseudoobservable_revtransforms::Dict{Symbol, Symbol}: saved identifiers for reverse transforms used for pseudoobservables
  • shock_names::Dict{Symbol, Int}: saved for var in [:histshocks, :forecastshocks, :shockdecstates, :shockdecobs, :shockdecpseudo]

Note that we don't save dates or transformations for impulse response functions.

source
DSGE.write_forecast_outputsMethod
write_forecast_outputs(m, output_vars, forecast_output_files, forecast_output;
    df = DataFrame(), block_number = Nullable{Int64}(), block_inds = 1:0,
    verbose = :low)

Writes the elements of forecast_output indexed by output_vars to file, given forecast_output_files, which maps output_vars to file names.

If output_vars contains :histobs, data must be passed in as df.

source
DSGE.write_meansbands_tableMethod
write_meansbands_table(dirname, filestring_base, mb, df, tablevar)

Inputs

  • dirname::String: directory to which tables are saved. Defaults to tablespath(m, "forecast")
  • filestring_base::Vector{String}: the result of filestring_base(m), typically ["vint=yymmdd"]`
  • mb::MeansBands: used for computing the output file name
  • df::DataFrame: the result of calling one of prepare_meansbands_table_timeseries, prepare_means_table_shockdec, or prepare_means_table, irf
  • tablevar::Symbol: used for computing the base output file name
source