Post-Estimation Analysis
The code for this page lives in src/SubModules/PostEstimation/ and is part of the PostEstimation submodule.
The PostEstimation tools help analyze model results after solving or estimating a model. This includes impulse responses (linear and distributional), variance decompositions, historical decompositions, and advanced diagnostics like GIRFs and unconditional moments.
1. Impulse Response Functions (IRFs)
Impulse responses trace the effect of a shock through time. In the linear case, these are computed from the state transition (hx) and state-to-control (gx) matrices derived from the linear solution.
compute_irfs: Given a set of shocks and selected variables, builds the impulse responses for a specified horizon using(State2Control, LOMstate). Accepts options for interval computation (via posterior draws) and returns both IRFs and optional bounds.plot_irfs: Plots a set of IRFs (possibly from multiple models) with legends, styles, and optional data overlays. Designed for quick comparison across shocks and variables.
BASEforHANK.PostEstimation.compute_irfs — Function
compute_irfs(
exovars,
gx,
hx,
XSS,
ids;
T = 1000,
init_val = fill(0.01, length(exovars)),
verbose = true,
irf_interval_options = nothing,
distribution = false,
transform_elements = nothing,
comp_ids = nothing,
n_par = nothing,
m_par = nothing,
)Computes impulse response functions (IRFs) for a given set of shocks to exogenous variables.
Arguments
exovars::Vector{Int64}: Vector of positional indices of exogenous variables to which shocks are applied.gx::Matrix{Float64}: Control matrix (states to controls equations)hx::Matrix{Float64}: Transition matrix for states (state transition equations)XSS::Vector{Float64}: Steady state values of the model variables.ids: Indexes of the model variables.
Keyword Arguments
T::Int64: Number of periods for which to compute IRFs.init_val::Vector{Float64}: Initial value of the shock to each exogenous variable, defaults to 0.01 for all shocks.verbose::Bool: Print progress to console.irf_interval_options: If provided, computes confidence intervals for IRFs based on parameter draws. Should be a dictionary with keys"draws"(matrix of parameter draws) and"e_set"(estimation settings struct).distribution::Bool: Compute distributional IRFs, defaults to false.comp_ids: Compression indices for the distribution, as created byprepare_linearization. Needed ifdistributionis true. Defaults to nothing.n_par: Model parameters, needed for distributional IRFs. Defaults to nothing.m_par: Model parameters, needed for distributional IRFs. Defaults to nothing.transform_elements: Transformation elements for the distributional IRFs. Defaults to nothing.
Returns
IRFs,IRFs_lvl: 3D array of (level) IRFs for each exogenous variable, with dimensions:- 1: States and controls
- 2: Time periods
- 3: Exogenous variables
IRFs_order: Names of the exogenous variables and their indices in the IRFs array.IRFs_dist: Dictionary containing the full distributional IRFs, if requested.
BASEforHANK.PostEstimation.plot_irfs — Function
plot_irfs(
shocks_to_plot,
vars_to_plot,
IRFs_to_plot,
IRFs_order,
ids;
intervals = false,
e_set = nothing,
plot_data_irfs::Bool = false,
horizon = 40,
factor = 100,
legend_on_all = false,
show_fig = true,
save_fig = false,
save_fig_indiv = false,
path = "",
suffix = "",
yscale = "standard",
style_options = (lw = 2, color = :auto, linestyle = :solid),
)Plots impulse response functions (IRFs) for specified shocks and variables, given a set of IRFs as produced by compute_irfs.
Arguments
shocks_to_plot::Vector{Tuple{Symbol,String}}: A vector of tuples, where each tuple contains a shock variable (as aSymbol) and its corresponding label (String).vars_to_plot::Vector{Tuple{Symbol,String}}: A vector of tuples, each containing a variable to plot (as aSymbol) and its corresponding label (String).IRFs_to_plot::Vector{Tuple{Array{Float64,3},String}}: A vector of tuples, where each tuple consists of a 3D array of IRFs (Array{Float64,3}) and a string representing the specification name.IRFs_order::Vector{Symbol}: A vector of symbols specifying the order of shocks in the IRF arrays.ids: A structure mapping variable symbols to their corresponding indices in the IRF arrays, must be identical for all IRFs versions.
Keyword Arguments
intervals::Vector{Tuple}: A vector of tuples, each tuple for each model specified. Each tuple contains the lower and upper bounds for the IRFs. Default is false (no intervals plotted).e_set::BASEforHANK.EconSet: An economic settings object containing relevant data and configurations. Needed for IRF matching plots.plot_data_irfs::Bool: Iftrue, overlays empirical IRFs from data on the plots. Default isfalse.horizon::Int64: The time horizon (number of periods) over which IRFs are plotted. Default is40.factor::Int64: Scaling factor for the IRFs (default:100).legend_on_all::Bool: Iftrue, includes the legend on all subplots; otherwise, only the first subplot has a legend. Default isfalse.show_fig::Bool: Iftrue, displays the plot. Default istrue.save_fig::Bool: Iftrue, saves the combined plot as a PDF. Default isfalse.save_fig_indiv::Bool: Iftrue, saves individual plots for each variable directly inside a shock-specific folder. Default isfalse.path::String: The directory path where the generated plots should be saved. Default is an empty string (no saving).suffix::String: A suffix to append to the saved plot filenames. Default is an empty string.yscale::Union{String, Tuple{Number,Number}, Dict{Symbol,Tuple{Number,Number}}}: Y-axis scaling specification. When set to"common", computes a common y-axis limit across all subplots from the data; when provided as a tuple, uses it as(ymin, ymax); when provided as a dictionary, applies specified y-axis limits for each variable. Default is"standard", which applies default scaling.style_options::NamedTuple: A named tuple specifying stylistic options for the plots, including line width (lw), color (default::auto), and linestyle (default::solid). Default is(lw = 2, color = :auto, linestyle = :solid).
Highlight: Distributional IRFs (heatmaps) for HANK
Distributional IRFs show how the cross-sectional distribution (e.g., of assets or income) shifts over time in response to a shock. This is crucial in HANK models where aggregate dynamics depend on the distribution.
BASEforHANK.PostEstimation.plot_distributional_irfs — Function
plot_distributional_irfs(
shocks_to_plot::Vector{Tuple{Symbol,String}},
vars_to_plot::Vector{Tuple{String,String}},
IRFs_to_plot::Dict{String,Array{Float64}},
IRFs_order::Vector{Symbol},
n_par;
bounds::Dict{String,Tuple{Float64,Float64}} = Dict{String,Tuple{Float64,Float64}}(),
horizon::Int64 = 40,
legend::Bool = false,
show_fig::Bool = true,
save_fig::Bool = false,
path::String = "",
fps::Int = 2,
suffix::String = ""
)Plots distributional impulse response functions (IRFs) for specified shocks and variables.
This function generates either static 3D surface plots for univariate distributions or animated 3D surface plots for bivariate distributions. It iterates through the specified shocks and variables, creating and optionally saving a plot for each combination.
Arguments
shocks_to_plot::Vector{Tuple{Symbol,String}}: A vector of tuples, where each tuple contains a shock variable (as aSymbol) and its corresponding label (String).vars_to_plot::Vector{Tuple{Symbol,String}}: A vector of tuples, each containing a variable to plot (as aSymbol) and its corresponding label (String).IRFs_to_plot::Dict{String,Array{Float64}}: A dictionary of distributional IRFs, where the keys are variable names and the values are arrays of IRFs.IRFs_order::Vector{Symbol}: A vector of symbols specifying the order of shocks in the IRF arrays.n_par: Model parameters, needed for distributional IRFs. Defaults to nothing.
Keyword Arguments
bounds::Dict{String,Tuple{Float64,Float64}}: A dictionary to set axis limits for specific variables. Keys can be "b", "k", "h". Example:Dict("k" => (0.0, 50.0)).horizon::Int64: The time horizon (number of periods) over which IRFs are plotted.legend::Bool: Toggles the display of the legend on plots. Default isfalse.show_fig::Bool: Iftrue, displays the plot (default:true).save_fig::Bool: Iftrue, saves the plot as a PDF. Default isfalse.path::String: The directory path where the generated plots should be saved. Default is an empty string (no saving).fps::Int: Frames per second for animated GIF outputs. Default is2.suffix::String: A suffix to append to the saved plot filenames. Default is an empty string.
Additional plotting helper:
plot_irfs_cat provides categorical plotting options to organize many IRFs by variables or shocks, useful for side-by-side comparisons and publication-ready figures.
BASEforHANK.PostEstimation.plot_irfs_cat — Function
plot_irfs_cat(
shock_categories,
vars_to_plot,
IRFs_to_plot,
IRFs_order,
ids;
horizon = 40,
factor = 100,
show_fig = true,
save_fig = false,
save_fig_indiv = false,
path = "",
suffix = "",
yscale = "standard",
style_options = (lw = 2, color = :auto, linestyle = :solid)
)Plots impulse response functions (IRFs) for specified shocks and variables, given IRFs as produced by compute_irfs, organized by shock categories.
Arguments
shock_categories::Dict{Tuple{String,String},Vector{Symbol}}: A dictionary where each key represents a category of shocks (with a label and a string for saving), and each value is a vector of symbols representing the shocks in that category.vars_to_plot::Vector{Tuple{Symbol,String}}: A vector of tuples, each containing a variable to plot (as aSymbol) and its corresponding label (String).IRFs_to_plot::Array{Float64,3}: A 3D array of IRFs.IRFs_order::Vector{Symbol}: A vector of symbols specifying the order of shocks in the IRF arrays.ids: A structure mapping variable symbols to their corresponding indices in the IRF arrays, which must be identical for all IRFs versions.
Keyword Arguments
horizon::Int64: The time horizon (number of periods) over which IRFs are plotted. Default is40.factor::Int64: Scaling factor for the IRFs (default:100).show_fig::Bool: Iftrue, displays the plot. Default istrue.save_fig::Bool: Iftrue, saves the combined plot as a PDF. Default isfalse.save_fig_indiv::Bool: Iftrue, saves individual plots for each variable or shock panel directly inside a category-specific folder. Default isfalse.path::String: The directory path where the generated plots should be saved. Default is an empty string (no saving).suffix::String: A suffix to append to the saved plot filenames. Default is an empty string.yscale::Union{String, Tuple{Number,Number}, Dict{Symbol,Tuple{Number,Number}}}: Y-axis scaling specification. When set to"common", computes a common y-axis limit across all subplots from the data; when provided as a tuple, uses it as(ymin, ymax); when provided as a dictionary, applies specified y-axis limits for each variable. Default is"standard", which applies default scaling.style_options::NamedTuple: A named tuple specifying stylistic options for the plots, including line width (lw), color (default::auto), and linestyle (default::solid). Default is(lw = 2, color = :auto, linestyle = :solid).
2. Variance Decomposition
Forecast Error Variance Decomposition (FEVD) attributes variance in forecast errors to different structural shocks over a given horizon. These tools compute and visualize the contribution of each shock.
compute_vardecomp: Uses the linear state-space representation to quantify how much each shock contributes to forecast errors of selected variables across horizons.plot_vardecomp: Visualizes FEVD results as stacked areas or bars by shock.
BASEforHANK.PostEstimation.compute_vardecomp — Function
compute_vardecomp(IRFs)Computes the variance decomposition of the impulse response functions (IRFs) of a model.
Arguments
IRFs::Array{Float64}: 3D array of IRFs for each exogenous variable, with dimensions, returned fromcompute_irfs:- 1: States and controls
- 2: Time periods
- 3: Exogenous variables
Returns
VDs: 3D array of variance decompositions for each exogenous variable, with dimensions:- 1: States and controls
- 2: Time periods
- 3: Exogenous variables
exovars_names: Names of the exogenous variables.
BASEforHANK.PostEstimation.plot_vardecomp — Function
plot_vardecomp(
vars_to_plot,
VDs_to_plot,
VDs_order,
ids;
shock_categories = Dict(),
VD_horizons = [4, 16, 100],
colorlist = [...],
show_fig = true,
save_fig = false,
path = "",
suffix = "",
factor = 100.0,
)Plots variance decompositions (VDs) for selected variables across multiple models and forecast horizons, given a set of VDs as produced by compute_vardecomp, potentially organized by shock categories.
Arguments
vars_to_plot::Vector{Tuple{Symbol,String}}: A vector of tuples, where each tuple consists of a variable to plot (as aSymbol) and its corresponding label (String).VDs_to_plot::Vector{Tuple{Array{Float64,3},String}}: A vector of tuples, where each tuple consists of a 3D array of VDs (Array{Float64,3}) and a string representing the specification name.VDs_order::Vector{Symbol}: A vector of symbols specifying the order of shocks in the VD arrays.ids: A structure mapping variable symbols to their corresponding indices in the VD arrays, must be identical for all VDs versions.
Keyword Arguments
shock_categories::Dict{Tuple{String,String},Vector{Symbol}}: A dictionary where each key represents a category of shocks (with a label and a string for saving), and each value is a vector of symbols representing the shocks in that category.VD_horizons::Vector{Int64}: Forecast horizons to consider. Default is[4, 16, 100].colorlist::Vector: A list of colors used to distinguish different shocks in the stacked bar plots. Default includes predefined hex colors and symbols.show_fig::Bool: Iftrue, displays the plot. Default istrue.save_fig::Bool: Iftrue, saves the plot as a PDF. Default isfalse.path::String: The directory path where the generated plots should be saved. Default is an empty string (no saving).suffix::String: A suffix to append to the saved plot filenames. Default is an empty string.factor::Float64 = 100.0: Scaling factor for variance decomposition values.
Business-cycle-frequency variant:
compute_vardecomp_bcfreq: Focuses FEVD on business-cycle frequency bands via spectral methods.plot_vardecomp_bcfreq: Plots the frequency-targeted variance contributions for clear comparison across shocks.
BASEforHANK.PostEstimation.compute_vardecomp_bcfreq — Function
compute_vardecomp_bcfreq(exovars_full, stds, gx, hx; passband = (6, 32), ngrid = 512)This function is designed to produce a variance decomposition at business cycle frequencies. It produces a variance decomposition of the linearized solution. It returns the variance decomposition at business cycle frequencies based on Uhlig (2001) and the unconditional variance.
Arguments
exovars_full::Vector{Int64}: Vector of positional indices of all exogenous variables to which shocks; note: here, no subset of exogenous variables is allowed.stds::Vector{Float64}: Vector of standard deviations of the shocks.gx::Matrix{Float64}: Control matrix (states to controls equations)hx::Matrix{Float64}: Transition matrix for states (state transition equations)
Keyword Arguments
passband::Tuple{Int64,Int64}: A tuple specifying the horizons associated with the business cycle. Default is (6, 32).ngrid::Int64: The number of grid points for the computation of the band pass filter. Default is 512.
Returns
VD: Variance decomposition at business cycle frequencies.UnconditionalVar: Unconditional variance.
BASEforHANK.PostEstimation.plot_vardecomp_bcfreq — Function
plot_vardecomp_bcfreq(
vars_to_plot,
VDbcs_to_plot,
VDbcs_order,
ids;
shock_categories = Dict(),
colorlist = [...],
show_fig = true,
save_fig = false,
path = "",
suffix = "",
factor = 100.0,
)Plots busincess cycle frequency variance decompositions (VDbcs) for selected variables across multiple models and forecast horizons, given a set of VDbcs as produced by compute_vardecomp_bcfreq, potentially organized by shock categories.
Arguments
vars_to_plot::Vector{Tuple{Symbol,String}}: A vector of tuples, where each tuple consists of a variable to plot (as aSymbol) and its corresponding label (String).VDbcs_to_plot::Vector{Tuple{Matrix{Float64},String}}: A vector of tuples, where each tuple consists of a matrix of VDbcs (Matrix{Float64}) and a string representing the specification name.VDbcs_order::Vector{Symbol}: A vector of symbols specifying the order of shocks in the VDbc arrays.ids: A structure mapping variable symbols to their corresponding indices in the VD arrays, must be identical for all VDbcs versions.
Keyword Arguments
shock_categories::Dict{Tuple{String,String},Vector{Symbol}}: A dictionary where each key represents a category of shocks (with a label and a string for saving), and each value is a vector of symbols representing the shocks in that category.colorlist::Vector: A list of colors used to distinguish different shocks in the stacked bar plots. Default includes predefined hex colors and symbols.show_fig::Bool: Iftrue, displays the plot. Default istrue.save_fig::Bool: Iftrue, saves the plot as a PDF. Default isfalse.path::String: The directory path where the generated plots should be saved. Default is an empty string (no saving).suffix::String: A suffix to append to the saved plot filenames. Default is an empty string.factor::Float64 = 100.0: Scaling factor for variance decomposition values.
3. Historical Decomposition
Historical decomposition recovers the contribution of each structural shock to the realized path of observables over a sample, consistent with the linear state-space representation.
compute_hist_decomp: Reconstructs each variable as a sum of shock-specific components plus initial conditions and measurement errors.plot_hist_decomp: Plots these components over time, highlighting the shocks driving historical movements.
BASEforHANK.PostEstimation.compute_hist_decomp — Function
compute_hist_decomp(exovars_full, gx, hx, smoother_output, ids)Computes the historical decomposition of the model's variables in response to a set of shocks, including the effect of initial conditions and each shock.
Arguments
exovars_full::Vector{Int64}: Vector of positional indices of all exogenous variables to which shocks; note: here, no subset of exogenous variables is allowed.gx::Matrix{Float64}: Control matrix (states to controls equations)hx::Matrix{Float64}: Transition matrix for states (state transition equations)smoother_output: The output from a smoother function.ids: Indexes of the model variables.
Returns
ShockContr::Array{Float64}: A 3D array representing the decomposition of the effects of each shock and the initial condition on states and controls over time. The dimensions are:- 1: States and controls
- 2: Time periods
- 3: Shocks (including an additional index for the initial condition)
ShockContr_order::Vector{Any}: A vector containing the names or identifiers of the shocks and the initial condition, in the order they appear in theShockContrarray.
BASEforHANK.PostEstimation.plot_hist_decomp — Function
plot_hist_decomp(
vars_to_plot,
HDs_to_plot,
HDs_order,
ids;
shock_categories = Dict(),
timeline = collect(1:size(HDs_to_plot, 2)),
colorlist = [...],
show_fig = true,
save_fig = false,
path = "",
suffix = ""
)Plots historical decompositions (HDs) for specified variables, showing the contribution of different shocks over time, given a set of HDs as produced by compute_hist_decomp, potentially organized by shock categories.
Arguments
vars_to_plot::Vector{Tuple{Symbol,String}}: A vector of tuples, where each tuple consists of a variable to plot (as aSymbol) and its corresponding label (String).HDs_to_plot::Array{Float64,3}: A 3D array containing historical decompositions, where dimensions correspond to variables, time periods, and shocks.HDs_order::Vector{Symbol}: A vector of symbols specifying the order of shocks in the HDs array.ids: A structure mapping variable symbols to their corresponding indices in the HDs array, ensuring consistency across all plotted variables.
Keyword Arguments
shock_categories::Dict{Tuple{String,String},Vector{Symbol}}: A dictionary where each key represents a category of shocks (with a label and a string for saving), and each value is a vector of symbols representing the shocks in that category.timeline::Vector: A vector specifying the time axis for the plots. Default is a sequence from1to the number of time periods inHDs_to_plot.colorlist::Vector: A list of colors used to distinguish different shocks in the stacked bar plots. Default includes predefined hex colors and symbols.show_fig::Bool: Iftrue, displays the plot. Default istrue.save_fig::Bool: Iftrue, saves the plot as a PDF. Default isfalse.path::String: The directory path where the generated plots should be saved. Default is an empty string (no saving).suffix::String: A suffix to append to the saved plot filenames. Default is an empty string.
4. Advanced Analysis (GIRF and Moments)
Generalized IRFs (GIRFs) capture potentially state-dependent or non-linear responses (e.g., sign-dependent shocks) when higher-order solutions are available. Unconditional moment tools summarize means and other moments implied by the solution.
GIRF_FO,GIRF_SO: Compute generalized IRFs using first- or second-order approximations, allowing for asymmetric or state-dependent responses.uncondFirstMoment_SO_analytical: Provides analytical unconditional first moments under the second-order solution, useful for checks and summaries.
BASEforHANK.PostEstimation.GIRF_FO — Function
GIRF_FO(l, ν, lr, sr; shock_names=shock_names, η=nothing)Compute the first-order Generalized Impulse Response Function (GIRF) up to horizon l for a given error vector ν.
BASEforHANK.PostEstimation.GIRF_SO — Function
GIRF_SO(l, ν, lr_full, lr_reduc, sr_reduc, sor_full, variances, xf_reduc; shock_names=shock_names, η_red=nothing, η_full=nothing)Compute the second-order Generalized Impulse Response Function (GIRF) at horizon l for a given error vector ν.
BASEforHANK.PostEstimation.uncondFirstMoment_SO_analytical — Function
uncondFirstMoment_SO_analytical(sr, lr, sor, covariance; shock_names=shock_names)Compute the unconditional first moments (mean) of states and controls using a second-order approximation.
Arguments
sr: steady-state resultslr: linear solution resultssor: second-order solution resultscovariance: variance-covariance matrix of shocksshock_names: names of the shocks (optional)
Returns
Ex: unconditional mean of statesEy: unconditional mean of controlsΣx: variance-covariance matrix of statesΣy: variance-covariance matrix of controls