Standard Algorithms
Solving the Model
DSGE.gensys
— Functiongensys(Γ0, Γ1, c, Ψ, Π)
gensys(Γ0, Γ1, c, Ψ, Π, div)
gensys(F::LinearAlgebra.GeneralizedSchur, c, Ψ, Π)
gensys(F::LinearAlgebra.GeneralizedSchur, c, Ψ, Π, div)
Generate state-space solution to canonical-form DSGE model.
System given as
Γ0*y(t) = Γ1*y(t-1) + c + Ψ*z(t) + Π*η(t),
with z an exogenous variable process and η being endogenously determined one-step-ahead expectational errors.
Returned system is
y(t) = G1*y(t-1) + C + impact*z(t) + ywt*inv(I-fmat*inv(L))*fwt*z(t+1)
Returned values are
G1, C, impact, fmat, fwt, ywt, gev, eu, loose
If z(t)
is i.i.d., the last term drops out.
If div
is omitted from argument list, a div
>1 is calculated.
Return codes
eu[1] = 1
for existenceeu[2] = 1
for uniquenesseu[1] = -1
for existence only with not-s.c. zeu = [-2, -2]
for coincident zeroseu = [-3, -3]
if a LAPACKException is thrown while computing the Schur decomposition
Notes
We constrain Julia to use the complex version of the schurfact
routine regardless of the types of Γ0
and Γ1
, to match the behavior of Matlab. Matlab always uses the complex version of the Schur decomposition, even if the inputs are real numbers.
DSGE.gensys2
— Functiongensys2(m::AbstractDSGEModel, Γ0::Matrix{Float64}, Γ1::Matrix{Float64}, C::Vector{Float64},
Ψ::Matrix{Float64}, Π::Matrix{Float64}, TTT::Matrix{Float64}, RRR::Matrix{Float64},
CCC::Vector{Float64})
calculates the state space transition matrices when a temporary alternative policy applies. This function is not the same as the gensys2 written by Chris Sims, which computes a second-order perturbation.
Optimization
DSGE.csminwel
— Functioncsminwel(fcn::Function, grad::Function, x0::Vector, H0::Matrix=1e-5.*eye(length(x0)), args...;
xtol::Real=1e-32, ftol::Float64=1e-14, grtol::Real=1e-8, iterations::Int=1000,
store_trace::Bool = false, show_trace::Bool = false, extended_trace::Bool = false,
verbose::Symbol = :none, rng::AbstractRNG = MersenneTwister(0), kwargs...)
Minimizes fcn
using the csminwel algorithm.
Arguments
fcn::Function
: The objective functiongrad::Function
: The gradient of the objective function. This argument can be omitted if
an analytical gradient is not available, which will cause a numerical gradient to be calculated.
x0::Vector
: The starting guess for the optimizer
Optional Arguments
H0::Matrix
: An initial guess for the Hessian matrix – must be
positive definite. If none is given, then a scaled down identity matrix is used.
args...
: Other positional arguments to be passed tof
on each
function call
Keyword Arguments
ftol::{T<:Real}=1e-14
: Threshold for convergence in terms of change
in function value across iterations.
iterations::Int=100
: Maximum number of iterationskwargs...
: Other keyword arguments to be passed tof
on each
function call
DSGE.optimize!
— Functionoptimize!(m::Union{AbstractDSGEModel,AbstractVARModel}, data::Matrix;
method::Symbol = :csminwel,
xtol::Real = 1e-32, # default from Optim.jl
ftol::Float64 = 1e-14, # Default from csminwel
grtol::Real = 1e-8, # default from Optim.jl
iterations::Int = 1000,
store_trace::Bool = false,
show_trace::Bool = false,
extended_trace::Bool = false,
mle::Bool = false, # default from estimate.jl
step_size::Float64 = .01,
toggle::Bool = true, # default from estimate.jl
verbose::Symbol = :none)
Wrapper function to send a model to csminwel (or another optimization routine).
Hessian Approximation
DSGE.hessian!
— Methodhessian!(m::Union{AbstractDSGEModel,AbstractVARModel}, x::Vector{T}, data::AbstractArray;
check_neg_diag::Bool = true, toggle::Bool = false,
verbose::Symbol = :none) where {T<:AbstractFloat}
Compute Hessian of DSGE/VAR posterior function evaluated at x.
DSGE.hessizero
— Methodhessizero(fcn::Function, x::Vector{T};
check_neg_diag::Bool=false,
verbose::Symbol=:none,
distr::Bool=true) where T<:AbstractFloat
Compute Hessian of function fcn
evaluated at x
.
Arguments
check_neg_diag
: Throw an error if any negative diagonal elements are detected.verbose
: Print verbose outputdistr
: Use available parallel workers to increase performance.
Sampling
DSGE.metropolis_hastings
— Functionfunction metropolis_hastings(propdist::Distribution,
loglikelihood::Function,
parameters::ParameterVector{S},
data::Matrix{T},
cc0::T,
cc::T;
n_blocks::Int64 = 1,
n_param_blocks::Int64 = 1,
n_sim::Int64 = 100,
n_burn::Int64 = 0,
mhthin::Int64 = 1,
adaptive_accpt::Bool = false,
α::T = 1.0,
c::T = 0.5,
verbose::Symbol = :low,
savepath::String = "mhsave.h5",
rng::MersenneTwister = MersenneTwister(0),
regime_switching::Bool = false,
toggle::Bool = true,
testing::Bool = false) where {S<:Number, T<:AbstractFloat}
Implements the Metropolis-Hastings MCMC algorithm for sampling from the posterior distribution of the parameters.
Arguments
proposal_dist
: The proposal distribution that Metropolis-Hastings begins sampling from.m
: The model objectdata
: Data matrix for observablescc0
: Jump size for initializing Metropolis-Hastings.cc
: Jump size for the rest of Metropolis-Hastings.
Optional Arguments
n_blocks::Int = 1
: Number of blocks (for memory-management purposes)n_param_blocks::Int = 1
: Number of parameter blocksn_sim::Int = 100
: Number of simulations per save block. Note: # saved observations will ben_sim * n_param_blocks * (n_blocks - b_burn)
. The # of total simulations will ben_sim * n_param_blocks * n_blocks * mhthin
.n_burn::Int = 0
: Length of burn-in periodmhthin::Int = 1
: Thinning parameter (for mhthin = d, keep only every dth draw)adaptive_accpt::Bool = false
: Whether or not to adaptively adjust acceptance prob.α::T = 1.0
: Tuning parameter (step size) for proposal density computation in adaptive casec::T = 0.5
: Tuning parameter (mixture proportion) for proposal density computation in adaptive caseregime_switching::Bool = false
: do the parameters involve regime-switching?toggle
: if true, toggle the fields of any regime-switching parameters to regime 1.verbose::Bool
: The desired frequency of function progress messages printed to standard out. One of:
- `:none`: No status updates will be reported.
- `:low`: Status updates provided at each block.
- `:high`: Status updates provided at each draw.
savepath::String = "mhsave.h5"
: String specifying path to output filerng::MersenneTwister = MersenneTwister(0)
: Chosen seed (overridden if testing = true)testing::Bool = false
: Conditional for use when testing (determines fixed seeding)
metropolis_hastings(propdist::Distribution, m::Union{AbstractDSGEModel,AbstractVARModel},
data::Matrix{T}, cc0::T, cc::T; filestring_addl::Vector{String} = [],
regime_switching::Bool = false, toggle::Bool = false,
verbose::Symbol = :low) where {T<:AbstractFloat}
Wrapper function for DSGE models which calls Metropolis-Hastings MCMC algorithm for sampling from the posterior distribution of the parameters.
Arguments
propdist
: The proposal distribution that Metropolis-Hastings begins sampling from.m
: The model objectdata
: Data matrix for observablescc0
: Jump size for initializing Metropolis-Hastings.cc
: Jump size for the rest of Metropolis-Hastings.
Optional Arguments
verbose
: The desired frequency of function progress messages printed to standard out. One of:
- `:none`: No status updates will be reported.
- `:low`: Status updates provided at each block.
- `:high`: Status updates provided at each draw.
filestring_addl
: additional strings to add to the names of output filesregime_switching
: do the parameters involve regime-switching?toggle
: if true, toggle the fields of any regime-switching parameters to regime 1.
State Space Filters and Smoothers
Sequential Monte Carlo
See SMC.jl.