Impulse Responses
We provide many different types of impulse responses for DSGEs, VARs, VECMs, DSGE-VARs, and DSGE-VECMs. The forecast step allows the user to automatically compute "structural" impulse responses specifically for DSGEs, but for some purposes, a user may just want impulse responses without having to compute any other quantities. We provide this functionality with the impulse_responses
function. See the end of this page for the docstrings of all available impulse response functions.
We overload impulse_responses to cover specific use cases. For any AbstractDSGEModel
, we can compute the impulse responses over a specified horizon for all endogenous state variables, observables, and pseudo-observables by running
m = AnSchorfheide()
system = compute_system(m)
horizon = 40
states_irf, obs_irf, pseudo_irf = impulse_response(system, horizon)
For an AbstractRepModel
(a subtype of AbstractDSGEModel
for representative agent models), we can also grab the impulse responses by running
states_irf, obs_irf, pseudo_irf = impulse_response(m, system)
This use case requires the user to add a setting under the key :impulse_response_horizons
, which is set by default to 40.
If a user wants to specify a subset of the exogenous shocks and the size of those shocks, the user can run
shock_names = [:g_sh, :b_sh]
shock_values = [1.0, 1.0]
impulse_responses(m, system, horizon, shock_names, shock_values)
For the response of an endogenous state or observable to a specific shock,
shock_name = :g_sh
var_name = :obs_gdp
var_value = 0.
impulse_responses(m, system, horizon, shock_name , var_name, var_value)
DSGE Impulse Responses
There are two categories of impulse responses for DSGEs provided by DSGE.jl. It is easy to distinguish them by examining the state space form of a DSGE model (see Solving):
Impulse responses in the first category are "structural" impulse responses, which are the response of states and observables to the exogenous structural shocks $\epsilon_t$.
Impulse responses in the second category are "observables-identified" impulse responses. First, we may suppose that the measurement equation generically follows
where $F(\cdot)$ is some function of the unknown states $s_t$, and $\eta_t$ are random innovation to the observables $y_t$. By innovations, we mean that these random variables are potentially endogenous shocks, i.e. shocks which do not have a causal interpretation. An "observables-identified" impulse response specifies a certain response of $y_t$ to the innovations $\eta_t$, and uses this response to identify the underlying structural shocks which are consistent with these innovations. A DSGE identifies these innovations using the state space form of a DSGE.
We provide three types of "observables-identified" impulse responses for DSGEs.
- Short-Run Cholesky
- Long-Run Cholesky
- Maximizing Explained Cyclical Variance
We document the details of the identification in the docstrings of these impulse response functions. For the first two types of impulse responses, search the docstrings at the end of the page for
function impulse_responses(system::System{S}, horizon::Int, permute_mat::AbstractMatrix{T},
shocks::AbstractVector{S} = Vector{S}(undef, 0);
restriction::Symbol = :short_run, flip_shocks::Bool = false,
get_shocks::Bool = false) where {S <: Real, T <: Number}
For the third type of impulse response, search for
function impulse_responses(system::System{S}, horizon::Int, frequency_band::Tuple{S,S},
n_obs_shock::Int; flip_shocks::Bool = false,
get_shocks::Bool = false) where {S <: Real}
VAR Impulse Responses
While we have not yet implemented a VAR model, we do have impulse responses often used on VARs because of DSGE-VARs. Consider the VAR
where $X_t$ is a matrix of the lags of $y_t$, $\beta$ are the VAR coefficients, and $\epsilon_t \sim N(0, \Omega)$ are the innovations to observables.
We provide three types of impulse responses, each of which provide a different way of identifying orthogonalized shocks from the innovations.
- Short-Run Cholesky
- Long-Run Cholesky
- Maximizing Explained Cyclical Variance
These impulse responses are named the same as the observables-identified impulse responses for DSGEs because they are considering the same response of observables to the innovations. However, the treatment of identification is different when using a VAR because the mathematical structure of a VAR is not the same as a DSGE's. As a result, there are slight differences between these impulse responses and the observables-identified DSGE impulse responses.
impulse_responses(β, Σ, n_obs_shock, horizon, shock_size = 1;
method = :cholesky, flip_shocks = false, use_intercept = true,
frequency_band = (2π/32, 2π/6)) where {S<:Real}
DSGE-VAR Impulse Responses
There are two types of impulse responses we can compute for a DSGE-VAR. For both types, we first draw from the posterior distributions of the $VAR(p)$ coefficients and innovations variance-covariance matrix, where $p$ is the number of lags. With these draws, we can do one of two things:
- Compute the VAR impulse response implied by the draws.
- Use the DSGE's structural impact response (i.e. the first period of an impulse response) to identify a mapping from the (endogenous) innovations in the VAR to the structural shocks of the DSGE.
The first type of impulse response uses the same code as the VAR impulse responses once we compute the coefficients and innovations variance-covariance matrix. We call the second type of impulse responses "DSGE-VAR rotation impulse responses" because we effectively use the DSGE to identify a rotation matrix.
For the first type of impulse response, see
function impulse_responses(m::AbstractDSGEVARModel{S}, data::AbstractArray{S}, method::Symbol,
n_obs_shock::Int; horizon::Int = 0 ,use_intercept::Bool = false,
flip_shocks::Bool = false, verbose::Symbol = :none) where {S <: Real}
function impulse_responses(m::AbstractDSGEVARModel{S}, method::Symbol,
n_obs_shock::Int; horizon::Int = 0 ,use_intercept::Bool = false,
flip_shocks::Bool = false, verbose::Symbol = :none) where {S <: Real}
The second function is for the specific case when $\lambda = \infty$, where the data does not matter for the impulse response.
For the second type of impulse responses, see
function impulse_responses(m::AbstractDSGEVARModel{S}, data::AbstractArray{S},
X̂::Matrix{S} = Matrix{S}(undef, 0, 0);
horizon::Int = 0, MM::Matrix{S} = Matrix{S}(undef, 0, 0),
flip_shocks::Bool = false, draw_shocks::Bool = false,
deviations::Bool = false,
verbose::Symbol = :none) where {S <: Real}
VECM Impulse Responses
While we have not yet implemented a VECM model, we do have impulse responses often used on VECMs because of DSGE-VECMs. Consider the VECM
where $e_{t - 1}$ are cointegrating relationships (the error correction terms); $X_t$ is a matrix of the lags of $\Delta y_t$; $\beta_e$ and $\beta_v$ are the VECM coefficients; and $\epsilon_t \sim N(0, \Omega)$ are the innovations to observables. We identify orthogonalized shocks for VECMs from the innovations using the short-run Cholesky method. Other methods have yet to be implemented, hence passing keywords specific to these methods will not do anything (namely frequency_band
). Find the docstring of the following function for details.
impulse_responses(β, Σ, coint_mat, n_obs_shock, horizon, shock_size = 1;
method = :cholesky, flip_shocks = false, use_intercept = true,
frequency_band = (2π/32, 2π/6)) where {S<:Real}
DSGE-VECM Impulse Responses
Most of the impulse responses for DSGE-VARs have been implemented for DSGE-VECMs. The two impulse responses that have not been implemented, due to the differences in VECMs and VARs, are the maxBC
and cholesky_long_run
impulse responses. For the first type of impulse responses, which use the VECM impulse response code, see the functions
function impulse_responses(m::AbstractDSGEVECMModel{S}, data::AbstractArray{S},
coint_mat::AbstractMatrix{S}, method::Symbol,
n_obs_shock::Int; horizon::Int = 0 ,use_intercept::Bool = false,
flip_shocks::Bool = false, verbose::Symbol = :none) where {S <: Real}
function impulse_responses(m::AbstractDSGEVECMModel{S}, coint_mat::AbstractMatrix{S}, method::Symbol,
n_obs_shock::Int; horizon::Int = 0 ,use_intercept::Bool = false,
flip_shocks::Bool = false, verbose::Symbol = :none) where {S <: Real}
The second function is for the specific case when $\lambda = \infty$, where the data does not matter for the impulse response.
For the second type of impulse responses, which we call rotation impulse responses, see
function impulse_responses(m::AbstractDSGEVECMModel{S}, data::AbstractArray{S},
coint_mat::AbstractMatrix{S},
X̂::Matrix{S} = Matrix{S}(undef, 0, 0);
horizon::Int = 0, MM::Matrix{S} = Matrix{S}(undef, 0, 0),
flip_shocks::Bool = false, draw_shocks::Bool = false,
deviations::Bool = false,
verbose::Symbol = :none) where {S <: Real}
Wrappers for Impulse Response Functions
The forecast_one
function provides a wrapper for computing structural DSGE impulse responses when drawing from a distribution of parameters and for saving these impulse responses as MeansBands
objects (see Computing Means and Bands).
However, forecast_one
currently cannot compute observables-identified DSGE impulse responses, VAR impulse responses, or DSGE-VAR impulse responses, and we do not plan on modifying forecast_one
to make it possible to do so. Instead, we provide three wrapper functions specifically for computing means and bands for these impulse responses. Please see their docstrings for details. The first wrapper is for observables-identified DSGE impulse responses, and the second two are for DSGE-VAR impulse responses. The first one applies generically to a DSGE-VAR with any prior weight $\lambda$, but the second one is a convenience wrapper for the case of $\lambda = \infty$, which is equivalent to computing the impulse responses of the VAR approximation to a DSGE.
No wrappers for DSGE-VECM impulse responses have been implemented because we have not constructed a DSGE model that can be interfaced with the DSGEVECM
type yet. As a result, there are no explicit test cases for these wrappers. We have decided against implementing wrappers for DSGE-VECM impulse responses until we have test cases to guarantee the wrappers do not have bugs.
For observables-identified DSGE impulse responses, find
function impulse_responses(m::AbstractDSGEModel, paras::Matrix{S},
input_type::Symbol, method::Symbol, n_obs_shock::Int,
output_vars::Vector{Symbol} =
[:irfstates, :irfobs, :irfpseudo]; parallel::Bool = false,
permute_mat::Matrix{S} = Matrix{Float64}(undef,0,0),
frequency_band::Tuple{S,S} = (2*π/32, 2*π/6),
flip_shocks::Bool = false,
density_bands::Vector{Float64} = [.5, .6, .7, .8, .9],
create_meansbands::Bool = false, test_meansbands::Bool = false,
minimize::Bool = true,
forecast_string::String = "",
do_rev_transform::Bool = false,
verbose::Symbol = :high) where {S<:Real}
For DSGE-VAR impulse responses, find
function impulse_responses(m::AbstractDSGEVARModel{S}, paras::Matrix{S},
data::Matrix{S}, input_type::Symbol, method::Symbol;
parallel::Bool = false,
frequency_band::Tuple{S,S} = (2*π/32, 2*π/6),
n_obs_shock::Int = 1, draw_shocks::Bool = false,
flip_shocks::Bool = false,
X̂::AbstractMatrix{S} = Matrix{S}(undef, 0, 0),
deviations::Bool = false
density_bands::Vector{Float64} = [.5, .6, .7, .8, .9],
create_meansbands::Bool = false, test_meansbands::Bool = false,
minimize::Bool = true,
forecast_string::String = "",
verbose::Symbol = :high) where {S<:Real}
For the impulse response of a VAR approximation to a DSGE, find
function impulse_responses(m::AbstractDSGEModel, paras::Union{Vector{S}, Matrix{S}},
input_type::Symbol, method::Symbol,
lags::Int, observables::Vector{Symbol},
shocks::Vector{Symbol},
n_obs_shock::Int; parallel::Bool = false,
frequency_band::Tuple{S,S} = (2*π/32, 2*π/6),
flip_shocks::Bool = false,
use_intercept::Bool = false,
density_bands::Vector{Float64} = [.5, .6, .7, .8, .9],
create_meansbands::Bool = false,
minimize::Bool = true,
forecast_string::String = "",
verbose::Symbol = :high) where {S<:Real}
Docstrings
DSGE.impulse_responses
— Functionimpulse_responses(m, system)
impulse_responses(system, horizon)
Compute impulse responses for a single draw.
Inputs
m::AbstractDSGEModel
: model objectsystem::System{S}
: state-space system matriceshorizon::Int
: number of periods ahead to forecastflip_shocks::Bool
: Whether to compute IRFs in response to a positive shock (by default the shock magnitude is a negative 1 std. shock)
where S<:AbstractFloat
Outputs
states::Array{S, 3}
: matrix of sizenstates
xhorizon
xnshocks
of state impulse response functionsobs::Array{S, 3}
: matrix of sizenobs
xhorizon
xnshocks
of observable impulse response functionspseudo::Array{S, 3}
: matrix of sizenpseudo
xhorizon
xnshocks
of pseudo-observable impulse response functions
function impulse_responses(system::System{S}, horizon::Int, permute_mat::AbstractMatrix{T},
shocks::AbstractVector{S} = Vector{S}(undef, 0);
restriction::Symbol = :short_run, flip_shocks::Bool = false,
get_shocks::Bool = false) where {S <: Real, T <: Number}
computes impulse responses using a Cholesky-identified shock to observables, with the ability to apply a permutation matrix.
Inputs
system::System{S}
: state-space system matriceshorizon::Int
: number of periods ahead to forecastpermute_mat::AbstractMatrix{T}
: permutation matrixshocks::AbstractVector{S}
: vector of orthogonal shocks to observables. See the section Restriction Types below for a mathematical explanation. If no vector is provided, then we assume the user wants a shock vector whose first entry is one and all other entries are zero (a 1 standard deviation shock to the innovation affecting the first variable).restriction::Symbol
: type of restriction for the Cholesky identification. Can be either:short_run
or:long_run
. We also let:cholesky
refer specifically to the:short_run
restriction and:choleskyLR
or:cholesky_long_run
refer to the:long_run
restriction.flip_shocks::Bool
: Whether to compute IRFs in response to a positive shock (by default the shock magnitude is a negative 1 std. shock)get_shocks::Bool
: Whether to return the Cholesky-identified shocks and the underlying structural shocks
Outputs
states::Matrix{S}
: matrix of sizenstates
xhorizon
xnshocks
of state impulse response functionsobs::Matrix{S}
: matrix of sizenobs
xhorizon
xnshocks
of observable impulse response functionspseudo::Matrix{S}
: matrix of sizenpseudo
xhorizon
xnshocks
of pseudo-observable impulse response functionscholesky_shock::Vector{S}
: Cholesky-identified shock on impact (only ifget_shocks = true
)structural_shock::Vector{S}
: structural shocks causing the Cholesky-identified shock (only ifget_shocks = true
)
Restriction Types
Consider a state space system
sₜ = TTT * sₜ₋₁ + RRR * sqrt(QQ) * ϵₜ
yₜ = ZZ * sₜ
where ϵₜ ∼ 𝒩 (0, I)
and QQ
is a diagonal matrix specifying the variances of the structural shocks ϵₜ
. Under this specification, the units of ϵₜ
are standard deviations.
A Cholesky identification searches for a lower triangular matrix M
such that M * M'
equals the covariance of observables. The restriction is short-run if we take the short-run covariance and is long run if we take the long-run covariance (i.e. the covariance of observables when shocks occur to the stationary level of sₜ). Explicitly, the short- and long-run restrictions are, respectively,
Covₛᵣ = (ZZ * RRR) * QQ * (ZZ * RRR)',
Covₗᵣ = (ZZ * (I - TTT)⁻¹ * RRR) * QQ * (ZZ * (I - TTT)⁻¹ * RRR)'.
Taking the Cholesky decomposition of the appropriate covariance matrix yields the lower triangular matrix M
.
Using M
, we may re-write the measurement equation of the state space as
yₜ = ZZ * TTT * sₜ₋₁ + M * uₜ,
where uₜ ∼ 𝒩 (0, I)
are orthogonal shocks identified from the covariance of yₜ
using M
. The orthogonal shocks uₜ
are precisely the input argument shocks
. We can map uₜ
to ϵₜ
, i.e. identify structural shocks, by solving the linear system
(ZZ * RRR * sqrt(QQ)) * ϵₜ = M * uₜ
when using the short-run restriction and
(ZZ * (I - TTT)⁻¹ * RRR * sqrt(QQ)) * ϵₜ = M * uₜ
when using the long-run restriction.
function impulse_responses(system::System{S}, horizon::Int, frequency_band::Tuple{S,S},
n_obs_shock::Int; flip_shocks::Bool = false,
get_shocks::Bool = false) where {S<:Real}
computes impulse responses by identifying the shock which maximizes the variance explained of a chosen observable within a certain frequency band.
For typical business-cycle frequences, we recommend setting frequency_band = (2 * π / 32, 2 * π / 6)
.
Inputs
system::System{S}
: state-space system matriceshorizon::Int
: number of periods ahead to forecastfrequency_band::Tuple{S, S}
: the frequencies at which we want to maximize the variance explained. The first frequency must be less than the second one.n_obs_shock::Int
: the index of the observable whose variance we want to maximally explain in the state space model implied bysystem
.flip_shocks::Bool
: Whether to compute IRFs in response to a positive shock (by default the shock magnitude is a negative 1 std. shock)get_shocks::Bool
: Whether to return the structural shocks.
where S <: Real
Outputs
states::Matrix{S}
: matrix of sizenstates
xhorizon
xnshocks
of state impulse response functionsobs::Matrix{S}
: matrix of sizenobs
xhorizon
xnshocks
of observable impulse response functionspseudo::Matrix{S}
: matrix of sizenpseudo
xhorizon
xnshocks
of pseudo-observable impulse response functionsstructural_shock::Vector{S}
: structural shocks causing the Cholesky-identified shock (only ifget_shocks = true
)
function impulse_responses(m::AbstractDSGEModel, paras::Union{Vector{S}, Matrix{S}},
input_type::Symbol, method::Symbol,
lags::Int, observables::Vector{Symbol},
shocks::Vector{Symbol},
n_obs_shock::Int; parallel::Bool = false,
frequency_band::Tuple{S,S} = (2*π/32, 2*π/6),
flip_shocks::Bool = false,
use_intercept::Bool = false,
density_bands::Vector{Float64} = [.5, .6, .7, .8, .9],
create_meansbands::Bool = false,
minimize::Bool = true,
forecast_string::String = "",
verbose::Symbol = :high) where {S<:Real}
computes the impulse responses of a VAR(p) approximation to a DSGE.
Inputs
m::Union{AbstractDSGEModel,AbstractDSGEVARModel}
: DSGE/DSGEVAR model objectparas::Matrix{S}
orparas::Vector{S}
: parameters to calibrate the modelinput_type::Symbol
::mode
specifies a modal impulse response, and:full
specifies a full-distribution forecast ifparas
is not given. This argument is also used to construct the file names of computedMeansBands
.method::Symbol
: type of impulse response to compute. The options are:cholesky
,:maximum_business_cycle_variance
or:maxBC
, and:cholesky_long_run
or:choleskyLR
. See?cholesky_shock
,?maxBC_shock
, and?choleskyLR_shock
.lags::Int
: number of lags in the VAR(p) approximation, i.e. p = lagsobservables::Vector{Symbol}
: observables to be used in the VAR. These can be any of the observables or pseudo-observables inm
.shocks::Vector{Symbol}
: (structural) exogenous shocks to be used in the DSGE-VAR. These shocks must be inm
.n_obs_shock::Int
: the index of the observable corresponding to the orthogonalized shock causing the impulse response.
Keywords
parallel::Bool
: use parallel workers or notfrequency_band::Tuple{S,S}
: See?maxBC_shock
.flip_shocks::Bool
: impulse response shocks are negative by default. Set totrue
for a positive signed shock.density_bands::Vector{Float64}
: bands for full-distribution IRF computationscreate_meansbands::Bool
: set totrue
to save output as aMeansBands
object.minimize::Bool
: choose shortest interval if true, otherwise just chop off lowest and highest (percent/2)forecast_string::String
: string tag for identifying this impulse responseverbose::Symbol
: quantity of output desired
function impulse_responses(m::AbstractDSGEVARModel{S}, paras::Matrix{S},
data::Matrix{S}, input_type::Symbol, method::Symbol;
parallel::Bool = false,
frequency_band::Tuple{S,S} = (2*π/32, 2*π/6),
n_obs_shock::Int = 1, draw_shocks::Bool = false,
flip_shocks::Bool = false,
X̂::AbstractMatrix{S} = Matrix{S}(undef, 0, 0),
deviations::Bool = false, normalize_rotation::Bool = false,
density_bands::Vector{Float64} = [.5, .6, .7, .8, .9],
create_meansbands::Bool = false, test_meansbands::Bool = false,
minimize::Bool = true,
forecast_string::String = "",
verbose::Symbol = :high) where {S<:Real}
computes the VAR impulse responses identified by the state space system
sₜ = TTT × sₜ₋₁ + RRR × impact[:, i],
yₜ = ZZ × sₜ + DD + MM × impact[:, i],
where impact[:, i]
is a linear combination of (orthogonal) structural shocks ϵₜ ∼ 𝒩 (0, I)
, and MM × impact[:, i]
are the correlated measurement errors.
The VAR impulse responses are computed according to
ŷₜ₊₁ = X̂ₜ₊₁β + uₜ₊₁,
where X̂ₜ₊₁
are the lags of observables in period t + 1
, i.e. yₜ, yₜ₋₁, ..., yₜ₋ₚ
.
If the method is :rotation
, the shock uₜ₊₁
is identified via
Σᵤ = 𝔼[u × u'] = chol(Σᵤ) × Ω × ϵₜ,
where the rotation matrix Ω
is the Q
matrix from a QR decomposition of the impact response matrix corresponding to the state space system, i.e.
Ω, _ = qr(∂yₜ / ∂ϵₜ').
Otherwise, we draw a β and Σᵤ from the posterior implied by the DSGE and data, and we then compute normal VAR impulse responses given those coefficients and innovations variance-covariance matrix.
NOTE: this function generally involves taking random draws from probability distributions, so seeds need to be set to achieve reproducibility. ****
Inputs
m::Union{AbstractDSGEModel,AbstractDSGEVARModel}
: DSGE/DSGEVAR model objectparas::Matrix{S}
orparas::Vector{S}
: parameters to calibrate the modelinput_type::Symbol
::mode
specifies a modal impulse response, and:full
specifies a full-distribution forecast ifparas
is not given. This argument is also used to construct the file names of computedMeansBands
.method::Symbol
: type of impulse response to compute. The options are:cholesky
and:rotation
. For the first, see?cholesky_shock
, and for the latter, we use the DSGE model to identify the rotation matrix which maps the DSGE's structural shocks to the innovations in the VAR's observables.lags::Int
: number of lags in the VAR(p) approximation, i.e. p = lagsobservables::Vector{Symbol}
: observables to be used in the VAR. These can be any of the observables or pseudo-observables inm
.shocks::Vector{Symbol}
: (structural) exogenous shocks to be used in the DSGE-VAR. These shocks must be inm
.n_obs_shock::Int
: the index of the observable corresponding to the orthogonalized shock causing the impulse response.
Keywords
parallel::Bool
: use parallel workers or notfrequency_band::Tuple{S,S}
: See?maxBC_shock
.n_obs_shock::Int
: Index of observable to be shocked when using a Cholesky-based impulse responsedraw_shocks::Bool
: true if you want to draw shocks along the entire horizonflip_shocks::Bool
: impulse response shocks are negative by default. Set totrue
for a positive signed shock.X̂::AbstractMatrix{S}
: matrix stacking the intercept and lags of the data for rotation IRFs. Set to a vector of zeros with length1 + n_observables * p
to compute the rotation IRFs in deviations from the baseline forecast.deviations::Bool
: set true to compute the impulse response in deviations rather than as a forecast. Mechnically, we ignoreX̂
(treated as zeros) and the intercept term.normalize_rotation::Bool
: set to true to normalize the rotation so that rows have the correct sign. This requires as many structural shocks as there are observables in the DSGE-VAR.density_bands::Vector{Float64}
: bands for full-distribution IRF computationscreate_meansbands::Bool
: set totrue
to save output as aMeansBands
object.minimize::Bool
: choose shortest interval if true, otherwise just chop off lowest and highest (percent/2)forecast_string::String
: string tag for identifying this impulse responseverbose::Symbol
: quantity of output desired
impulse_responses(β, Σ, n_obs_shock, horizon, shock_size = 1;
method = :cholesky, flip_shocks = false, use_intercept = true,
frequency_band = (2π/32, 2π/6)) where {S<:Real}
computes the impulse responses of a VAR system represented in the form
yₜ = Xₜβ + ϵₜ,
where Xₜ
stacks the lags of yₜ (with dimensions nobservables x nregressors), and
ϵₜ ∼ 𝒩 (0, Σ).
Inputs
β::AbstractMatrix{S}
: coefficient matrixΣ::AbstractMatrix{S}
: innovations variance-covariance matrixn_obs_shock::Int
: index of the observable corresponding to the orthogonalized shock causing the impulse response.shock_size::S
: number of standard deviations of the shock
Keywords
method::Symbol
: type of impulse response to compute. The available options are:cholesky
(default),:maximum_business_cycle_variance
or:maxBC
, and:cholesky_long_run
or:choleskyLR
. See?cholesky_shock
,?maxBC_shock
, and?cholesky_long_run_shock
.flip_shocks::Bool
: by default, we compute the impulse responses to a negative shock. Setflip_shocks = true
to obtain a positive shock.use_intercept::Bool
:impulse_responses
assumesβ
has constant term(s). If there are no such terms, thenuse_intercept
must be set tofalse
.frequency_band::Tuple{S,S}
: See?maxBC_shock
.
Outputs
Y::AbstractMatrix
: Impulse response matrix with dimensions horizons x n_observables
impulse_responses(β, Σ, coint_mat, n_obs_shock, horizon, shock_size = 1;
method = :cholesky, flip_shocks = false, use_intercept = true,
frequency_band = (2π/32, 2π/6)) where {S<:Real}
computes the impulse responses of a VECM system represented in the form
Δyₜ = eₜ₋₁ βₑ + Xₜ βᵥ + ϵₜ,
where βₑ
are the coefficients for the error correction terms; eₜ
are the error correction terms specifying the cointegrating relationships; βᵥ
are the coefficients for the VAR terms (including the intercept); Xₜ
are the lags of observables in period t
, i.e. yₜ₋₁, yₜ₋2, ..., yₜ₋ₚ
; and ϵₜ ∼ 𝒩 (0, Σ)
. We assume that β = [βₑ; βᵥ]
.
Inputs
β::AbstractMatrix{S}
: coefficient matrixΣ::AbstractMatrix{S}
: innovations variance-covariance matrixcoint_mat::AbstractMatrix{S}
: matrix specifying the cointegrating relationships. Multiplyingcoint_mat * data
, wheredata
is ann_observables × T
matrix, should yield ann_coint × T
matrix, wheren_coint
are the number of cointegrating relationships andT
are the number of periods of data.n_obs_shock::Int
: index of the observable corresponding to the orthogonalized shock causing the impulse response.shock_size::S
: number of standard deviations of the shock
Keywords
method::Symbol
: type of impulse response to compute. The available option is:cholesky
(default).flip_shocks::Bool
: by default, we compute the impulse responses to a negative shock. Setflip_shocks = true
to obtain a positive shock.use_intercept::Bool
:impulse_responses
assumesβ
has constant term(s). If there are no such terms, thenuse_intercept
must be set tofalse
.frequency_band::Tuple{S,S}
: See?maxBC_shock
.
Outputs
Y::AbstractMatrix
: Impulse response matrix with dimensions horizons x n_observables
function impulse_responses(m::AbstractDSGEVARModel{S}, data::AbstractArray{S}, method::Symbol,
n_obs_shock::Int; horizon::Int = 0,
flip_shocks::Bool = false, verbose::Symbol = :none) where {S <: Real}
function impulse_responses(m::AbstractDSGEVARModel{S}, method::Symbol,
n_obs_shock::Int; horizon::Int = 0, use_intercept::Bool = false,
flip_shocks::Bool = false, verbose::Symbol = :none) where {S <: Real}
computes the VAR impulse responses identified by the DSGE
sₜ = TTT × sₜ₋₁ + RRR × impact[:, i],
yₜ = ZZ × sₜ + DD + MM × impact[:, i],
where impact[:, i]
is a linear combination of (orthogonal) structural shocks ϵₜ ∼ 𝒩 (0, I)
, and MM × impact[:, i]
are the correlated measurement errors.
We draw a β and Σᵤ from the posterior implied by the DSGE and data, and we then compute normal VAR impulse responses given those coefficients and innovations variance-covariance matrix. The weight placed on the DSGE is encoded by the field λ
of the DSGEVAR object m
.
Given β, Σᵤ, we compute impulse responses to the VAR system
ŷₜ₊₁ = X̂ₜ₊₁β + uₜ₊₁,
where X̂ₜ₊₁
are the lags of observables in period t + 1
, i.e. yₜ, yₜ₋₁, ..., yₜ₋ₚ₊₁
using one of the available identification methods for VARs
If the second function is used (where data
is not an input), then we assume the user wants to compute the VAR approximation of the DSGE, regardless of the λ
value in m
. Note that this function will not update the value of λ
in m
(even though we are computing the DSGE-VAR(∞) approximation).
Inputs
method::Symbol
: The available methods are:cholesky
,:maxBC
, and:choleskyLR
. See the docstringsimpulse_responses
for VARs specifically.n_obs_shock::Int
: The index of the observable corresponding to the orthogonalized shock causing the impulse response.
Keyword Arguments
horizon::Int
: the desired horizon of the impulse responses.use_intercept::Bool
: use an intercept term for the VAR approximationflip_shocks::Bool
: default is a "negative" impulse response on impact. Set totrue
for the positive impulse response.
function impulse_responses(m::AbstractDSGEVARModel{S}, data::AbstractArray{S},
X̂::Matrix{S} = Matrix{S}(undef, 0, 0);
horizon::Int = 0, MM::Matrix{S} = Matrix{S}(undef, 0, 0),
flip_shocks::Bool = false, draw_shocks::Bool = false,
deviations::Bool = false, normalize_rotation::Bool = false,
verbose::Symbol = :none) where {S <: Real}
computes the VAR impulse responses identified by the DSGE
sₜ = TTT × sₜ₋₁ + RRR × impact[:, i],
yₜ = ZZ × sₜ + DD + MM × impact[:, i],
where impact[:, i]
is a linear combination of (orthogonal) structural shocks ϵₜ ∼ 𝒩 (0, I)
, and MM × impact[:, i]
are the correlated measurement errors.
The VAR impulse responses are computed according to
ŷₜ₊₁ = X̂ₜ₊₁β + uₜ₊₁,
where X̂ₜ₊₁
are the lags of observables in period t + 1
, i.e. yₜ, yₜ₋₁, ..., yₜ₋ₚ
. Note these impulse responses are not computed in deviations from the baseline forecast ŷₜ₊₁ = X̂ₜ₊₁β
. To compute these impulse responses, use the keyword deviations
.
The shock uₜ₊₁
is identified by assuming
Σᵤ = 𝔼[u × u'] = chol(Σᵤ) × Ω × ϵₜ,
where the rotation matrix Ω
is the Q
matrix from a QR decomposition of the impact response matrix corresponding to the state space system, i.e.
Ω, _ = qr(∂yₜ / ∂ϵₜ').
For reference, see Del Negro and Schorfheide (2004), Del Negro and Schorfheide (2006), and Del Negro and Schorfheide (2009).
Inputs
X̂::Matrix{S}
: covariates for the first "forecast" period of the impulse response, i.e. if we have a VAR withp
lags, then
X̂ = [1, ŷₜ, ŷₜ₋₁, ..., ŷₜ₋ₚ₊₁]
so that, when β is the vector of VAR coefficients, then
𝔼[ŷₜ₊₁] = kron(I, X̂') * β.
Internally, we do equivalent matrix operations to avoid allocating the Kronecker product.
NOTE: this function generally involves taking random draws from probability distributions, so seeds need to be set to achieve reproducibility. ****
Keywords
horizon::Int
: horizon of impulse responsesflip_shocks::Bool
: impulse response shocks are negative by default. Set totrue
for a positive signed shock.draw_shocks::Bool
: true if you want to draw shocks along the entire horizondeviations::Bool
: set true to compute the impulse response in deviations rather than as a forecast. Mechnically, we ignoreX̂
(treated as zeros) and the intercept term.normalize_rotation::Bool
: set to true to normalize the rotation so that rows have the correct sign. Requires same number of structural shocks as observables.verbose::Symbol
: quantity of output desired
function impulse_responses(TTT::Matrix{S}, RRR::Matrix{S}, ZZ::Matrix{S},
DD::Vector{S}, MM::Matrix{S}, impact::Matrix{S}, horizon::Int;
accumulate::Bool = false,
cum_inds::Union{Int,UnitRange{Int},Vector{Int}} = 0) where {S<:Real}
computes the impulse responses of the linear state space system to linear combinations of (orthogonal) structural shocks specified by impact
. Measurement error that is correlated with the impact matrix is allowed. We also include the option to accumulate certain observables.
This state space model takes the form
sₜ = TTT × sₜ₋₁ + RRR × impact[:, i],
yₜ = ZZ × sₜ + DD + MM × impact[:, i],
where impact[:, i]
is a linear combination of orthogonal structural shocks with mean zero and identity covariance, and MM × impact[:, i]
are the correlated measurement errors.
The impact
matrix must be nshock × nirf
, where nshock
is the number of structural shocks and nirf
is the number of desired impulse responses. For each row of impact
, we compute the corresponding impulse responses.
A standard choice for impact
is a square diagonal matrix. In this case, we compute the impulse response of observables to each structural shock, scaled by the desired size of the shock.
Keywords
accumulate
: set to true if an observable should be accumulated.cum_inds
: specifies the indices of variables to accumulated.
Outputs
irf_results::Matrix{S}
: anobs x horizon × nirf
matrix, wherenobs
is the number of observables.
function impulse_responses(TTT::Matrix{S}, RRR::Matrix{S}, ZZ::Matrix{S},
DD::Vector{S}, MM::Matrix{S}, QQ::Matrix{S},
k::Int, β::Matrix{S}, Σ::Matrix{S},
horizon::Int, X̂::Matrix{S} = zeros(S, k);
flip_shocks::Bool = false, draw_shocks::Bool = false,
deviations::Bool = false, normalize_rotation::Bool = false,
test_shocks::Matrix{S} =
Matrix{S}(undef, 0, 0)) where {S<:Real}
computes the VAR impulse responses identified by the state space system
sₜ = TTT × sₜ₋₁ + RRR × ϵₜ
yₜ = ZZ × sₜ + DD + MM × ϵₜ
where ϵₜ ∼ 𝒩 (0, QQ)
and MM × ϵₜ
are the correlated measurement errors.
The VAR impulse responses are computed according to
ŷₜ₊₁ = X̂ₜ₊₁β + uₜ₊₁,
where X̂ₜ₊₁
are the lags of observables in period t + 1
, i.e. yₜ, yₜ₋₁, ..., yₜ₋ₚ
. Note these impulse responses are not computed in deviations from the baseline forecast ŷₜ₊₁ = X̂ₜ₊₁β
. To compute these impulse responses, set the keyword deviations = true
.
The shock uₜ₊₁
is identified via
Σᵤ = 𝔼[u × u'] = chol(Σᵤ) × Ω × ϵₜ,
where the rotation matrix Ω
is the Q
matrix from a QR decomposition of the impact response matrix corresponding to the state space system, i.e.
Ω, _ = qr(∂yₜ / ∂ϵₜ').
To normalize the rotation so that rows have the correct sign, set normalize_rotation = true
. This requires that there are as many structural shocks as observables.
For reference, see Del Negro and Schorfheide (2004), Del Negro and Schorfheide (2006), and Del Negro and Schorfheide (2009).
function impulse_responses(m::AbstractDSGEVECMModel{S}, data::AbstractArray{S},
coint_mat::AbstractMatrix{S}, method::Symbol,
n_obs_shock::Int; horizon::Int = 0,
flip_shocks::Bool = false, verbose::Symbol = :none) where {S <: Real}
function impulse_responses(m::AbstractDSGEVECMModel{S}, coint_mat::AbstractMatrix{S}, method::Symbol,
n_obs_shock::Int; horizon::Int = 0, use_intercept::Bool = false,
flip_shocks::Bool = false, verbose::Symbol = :none) where {S <: Real}
computes the VECM impulse responses identified by the DSGE
sₜ = TTT × sₜ₋₁ + RRR × impact[:, i],
yₜ = ZZ × sₜ + DD + MM × impact[:, i],
where impact[:, i]
is a linear combination of (orthogonal) structural shocks ϵₜ ∼ 𝒩 (0, I)
, and MM × impact[:, i]
are the correlated measurement errors.
We draw a β and Σᵤ from the posterior implied by the DSGE and data, and we then compute normal VECM impulse responses given those coefficients, innovations variance-covariance matrix, and the matrix specifying cointegrating relationships in observables. The weight placed on the DSGE is encoded by the field λ
of the DSGEVECM object m
.
Given β, Σᵤ, we compute impulse responses using one of the available identifiction strategies to the VECM system
Δŷₜ₊₁ = eₜβₑ + X̂ₜ₊₁βᵥ + uₜ₊₁,
where βₑ
are the coefficients for the error correction terms; eₜ₊₁
are the error correction terms specifying the cointegrating relationships; βᵥ
are the coefficients for the VAR terms; X̂ₜ₊₁
are the lags of observables in period t + 1
, i.e. yₜ, yₜ₋₁, ..., yₜ₋ₚ
, and uₜ₊₁ ∼ 𝒩 (0, Σ)
.
If the second function is used (where data
is not an input), then we assume the user wants to compute the VECM approximation of the DSGE, regardless of the λ
value in m
. Note that this function will not update the value of λ
in m
(even though we are computing the DSGE-VECM(∞) approximation).
Inputs
coint_mat::AbstractMatrix{S}
: matrix specifying the cointegrating relationships in observables. Given a matrixdata
with dimensionsn_observables × T
, multiplyingcoint_mat * data
should yield an_coint × T
matrix, wheren_coint
is the number of cointegrating relationships andT
is the number of periods of data.method::Symbol
: The available methods are:cholesky
,:maxBC
, and:choleskyLR
. See the docstringsimpulse_responses
for VECMs specifically.n_obs_shock::Int
: The index of the observable corresponding to the orthogonalized shock causing the impulse response.
Keyword Arguments
horizon::Int
: the desired horizon of the impulse responses.use_intercept::Bool
: use an intercept term for the VECM approximationflip_shocks::Bool
: default is a "negative" impulse response on impact. Set totrue
for the positive impulse response.
function impulse_responses(m::AbstractDSGEVECMModel{S}, data::AbstractArray{S},
X̂::Matrix{S} = Matrix{S}(undef, 0, 0);
horizon::Int = 0, MM::Matrix{S} = Matrix{S}(undef, 0, 0),
flip_shocks::Bool = false, draw_shocks::Bool = false,
verbose::Symbol = :none) where {S <: Real}
computes the VECM impulse responses identified by the DSGE
sₜ = TTT × sₜ₋₁ + RRR × impact[:, i],
yₜ = ZZ × sₜ + DD + MM × impact[:, i],
where impact[:, i]
is a linear combination of (orthogonal) structural shocks ϵₜ ∼ 𝒩 (0, I)
, and MM × impact[:, i]
are the correlated measurement errors.
The VECM impulse responses are computed according to
Δŷₜ₊₁ = eₜβₑ + X̂ₜ₊₁βᵥ + uₜ₊₁,
where βₑ
are the coefficients for the error correction terms; eₜ
are the error correction terms specifying the cointegrating relationships; βᵥ
are the coefficients for the VAR terms; X̂ₜ₊₁
are the lags of observables in period t + 1
, i.e. yₜ, yₜ₋₁, ..., yₜ₋ₚ
, and uₜ₊₁ ∼ 𝒩 (0, Σ)
. Note these impulses responses are not computed in deviations from the baseline forecast Δŷₜ₊₁ = eₜ₊₁βₑ + X̂ₜ₊₁βᵥ
. To compute these impulse responses, use the keyword deviations
.
The shock uₜ₊₁
is identified by assuming
Σᵤ = 𝔼[u × u'] = chol(Σᵤ) × Ω × ϵₜ,
where the rotation matrix Ω
is the Q
matrix from a QR decomposition of the impact response matrix corresponding to the state space system, i.e.
Ω, _ = qr(∂yₜ / ∂ϵₜ').
The impact response matrix is constructed using only the stationary component of the state space system and ignores the cointegration components of ZZ
and DD
.
For reference, see Del Negro and Schorfheide (2004), Del Negro and Schorfheide (2006), and Del Negro and Schorfheide (2009).
Inputs
coint_mat::Matrix{S}
: matrix specifying the cointegrating relationships in the actualdata
matrix. Evaluatingcoint_mat * data
should yield a time series of the cointegrating relationships.X̂::Matrix{S}
: covariates for the first "forecast" period of the impulse response, i.e. if we have a VECM withp
lags, then
X̂ = [eₜ, 1, ŷₜ, ŷₜ₋₁, ..., ŷₜ₋ₚ₊₁]
so that, when β is the vector of VECM coefficients, then
𝔼[ŷₜ₊₁] = kron(I, X̂') * β.
Internally, we do equivalent matrix operations to avoid allocating the Kronecker product.
NOTE: this function generally involves taking random draws from probability distributions, so seeds need to be set to achieve reproducibility. ****
Keywords
horizon::Int
: horizon of impulse responsesflip_shocks::Bool
: impulse response shocks are negative by default. Set totrue
for a positive signed shock.draw_shocks::Bool
: true if you want to draw shocks along the entire horizondeviations::Bool
: set true to compute the impulse response in deviations rather than as a forecast. Mechnically, we ignoreX̂
(treated as zeros) and the intercept term.verbose::Symbol
: quantity of output desired
function impulse_responses(TTT::Matrix{S}, RRR::Matrix{S}, ZZ::Matrix{S},
DD::Vector{S}, MM::Matrix{S}, QQ::Matrix{S},
k::Int, n_obs::Int, n_coint::Int, β::Matrix{S}, Σ::Matrix{S},
coint_mat::Matrix{S}, horizon::Int, X̂::Matrix{S} = zeros(S, k);
flip_shocks::Bool = false, draw_shocks::Bool = false,
deviations::Bool = false,
test_shocks::Matrix{S} =
Matrix{S}(undef, 0, 0)) where {S<:Real}
computes the VECM impulse responses identified by the state space system
sₜ = TTT × sₜ₋₁ + RRR × ϵₜ
yₜ = ZZ × sₜ + DD + MM × ϵₜ
where ϵₜ ∼ 𝒩 (0, QQ)
and MM × ϵₜ
are the correlated measurement errors.
Consider the VECM
Δŷₜ₊₁ = eₜβₑ + X̂ₜ₊₁βᵥ + uₜ₊₁,
where βₑ
are the coefficients for the error correction terms; eₜ
are the error correction terms specifying the cointegrating relationships; βᵥ
are the coefficients for the VAR terms (including the intecept)o; X̂ₜ₊₁
are the lags of observables in period t + 1
, i.e. yₜ, yₜ₋₁, ..., yₜ₋ₚ₊₁
; and uₜ₊₁ ∼ 𝒩 (0, Σ)
. Note these impulses responses are not computed in deviations from the baseline forecast Δŷₜ₊₁ = eₜ₊₁βₑ + X̂ₜ₊₁βᵥ
. To compute these impulse responses, set the keyword deviations = true
.
The shock uₜ₊₁
is identified via
Σᵤ = 𝔼[u × u'] = chol(Σᵤ) × Ω × ϵₜ,
where the rotation matrix Ω
is the Q
matrix from a QR decomposition of the impact response matrix corresponding to the state space system, i.e.
Ω, _ = qr(∂yₜ / ∂ϵₜ').
The impact response matrix is constructed using only the stationary component of the state space system and ignores the cointegration components of ZZ
and DD
.
The data are assumed to have dimensions n_obs × T
, and the cointegration relationships in the data are given by coint_mat * data
, where coint_mat
has dimensions n_coint × n_obs
. The variable k
is the number of total regressors in the VECM, including cointegration terms.
For reference, see Del Negro and Schorfheide (2004), Del Negro and Schorfheide (2006), and Del Negro, Schorfheide, Smets, and Wouters (2007).