Computation of the steady state and dimensionality reduction

Note

Most of the code for this section lives in the submodule SteadyState (src/SubModules/SteadyState/). The function prepare_linearization is generated by the preprocessor using the template src/Preprocessor/template_fcns/prepare_linearization.jl. The function then lives in the submodule PerturbationSolution (src/SubModules/PerturbationSolution/).

The model features uninsured idiosyncratic income shocks h and, in the baseline configuration, two assets: liquid assets (bonds) b and illiquid assets (capital) k. Entrepreneurs (last income state) receive no labor income but firm profits, while workers additionally receive labor‑union profits.[BBL]

The steady‑state equilibrium consists of

  • the aggregate capital stock and other aggregate quantities,
  • marginal value functions V_b and V_k on the three‑dimensional grid ($b \times k \times h$), and
  • the ergodic joint distribution of households over these idiosyncratic states.

We reduce the dimensionality of the problem following [BL] by applying a Discrete Cosine Transform (DCT) to the marginal value functions and by approximating the joint distribution with a copula and state‑dependent marginals.

The models baseline household problem is described in detail in the Computational Notes and the description of the Household Problem. Yet, the model can also be configured as a one‑asset model or as a complete‑markets model. There are also other options for the household problem. Finally, different transition types for the distribution can be chosen.

At a high level, the workflow is:

  1. call_find_steadystate → solves the household steady state via BASEforHANK.SteadyState.find_steadystate.
  2. call_prepare_linearization → computes distributional statistics, performs the first‑stage DCT reduction, and prepares the objects needed for perturbation around the steady state (via BASEforHANK.PerturbationSolution.prepare_linearization).
  3. compute_steadystate is a convenience wrapper that executes both steps in sequence.

compute_steadystate

compute_steadystate(m_par) runs the full steady‑state pipeline in one call:

  1. Calls call_find_steadystate to solve the household steady state.
  2. Calls call_prepare_linearization to construct the state and control vectors, dimensionality‑reduction indices, and related bookkeeping structures.

It returns a SteadyResults struct that bundles all objects needed for subsequent linearization and estimation.

call_find_steadystate and find_steadystate

BASEforHANK.call_find_steadystateFunction
call_find_steadystate(m_par; n_par_kwargs::NamedTuple=NamedTuple())

Computes the steady state of the household problem and fills the SteadyStateStruct struct (without further steps of preparing the linearization).

Arguments

  • m_par::ModelParameters
  • n_par_kwargs::NamedTuple=NamedTuple(): Additional keyword arguments passed to NumericalParameters() when initializing the numerical parameters.

Returns

source
BASEforHANK.SteadyState.find_steadystateFunction
find_steadystate(m_par; n_par_kwargs::NamedTuple = NamedTuple())

Find the stationary equilibrium capital stock as well as the associated value functions and the stationary distribution of idiosyncratic states.

This function solves for the market clearing capital stock in the Aiyagari model with idiosyncratic income risk. It uses CustomBrent() to find the root of the excess capital demand function, which is defined in Kdiff(). The solution is obtained first on a coarse grid and then on the actual grid.

Arguments

  • m_par::ModelParameters
  • n_par_kwargs::NamedTuple=NamedTuple(): Additional keyword arguments passed to NumericalParameters() when initializing the numerical parameters.

Returns

  • KSS: Steady-state capital stock
  • vfSS: Value functions
  • distrSS::Array{Float64,3}: Steady-state distribution of idiosyncratic states
  • n_par::NumericalParameters, m_par::ModelParameters
source

call_find_steadystate(m_par; n_par_kwargs...) is the main entry point for solving the stationary household equilibrium. It returns a SteadyStateStruct(KSS, vfSS, distrSS, n_par) with

  • KSS: steady‑state aggregate capital (and, in two‑asset versions, the corresponding steady‑state bond position),
  • vfSS: steady‑state marginal value functions (ValueFunctions struct),
  • distrSS: stationary distribution over (b,k,h), and
  • n_par: the NumericalParameters used in the solution.

Internally, it calls BASEforHANK.SteadyState.find_steadystate, which implements the following steps:

Initialize numerical parameters

Construct n_par = NumericalParameters(; m_par, n_par_kwargs...). This object contains, among others:

  • the model variant (CompleteMarkets, OneAsset, or TwoAsset),
  • the transition type for aggregate shocks and income (LinearTransition, NonLinearTransition, …),
  • grid sizes for assets and income (nb, nk, nh),
  • the copula discretization (nb_copula, nk_copula, nh_copula),
  • options for the steady‑state search (e.g. rKmin, rKmax, ϵ).

The income process grid and transition matrix are built using BASEforHANK.Tools.Tauchen, based on the persistence and variance implied by m_par. An entrepreneurial state is appended if n_par.entrepreneur == true.

Handle complete‑markets shortcut

If the model is of type CompleteMarkets, the code skips the full steady‑state search and uses an analytical complete‑markets interest rate and capital function (if provided by the user).

Search for the equilibrium capital stock

For incomplete‑markets models, find_steadystate uses a customized Brent root‑finding algorithm BASEforHANK.Tools.CustomBrent to find the steady‑state capital stock KSS. The equilibrium condition is expressed in terms of BASEforHANK.SteadyState.Kdiff, the excess‑capital function.

  • For each candidate capital stock KD, Kdiff computes the implied interest rate and wage.
  • Given these prices, it calls BASEforHANK.SteadyState.Ksupply to solve the household problem and obtain aggregate asset supplies.
  • It then returns the difference between supplied and demanded capital.

The search is done first on a coarse grid (controlled by coarse‑grid parameters in NumericalParameters) and then refined on the full grid.

Household problem and stationary distribution

Ksupply itself:

  • uses BASEforHANK.SteadyState.EGM_policyupdate to compute optimal policies given tomorrow’s marginal value functions, current prices and incomes, and
  • constructs the sparse transition matrix over the household state space using internal routines such as MakeWeights, MakeTransition, and MultipleDirectTransitions,
  • finds the stationary distribution as the (normalized) eigenvector associated with eigenvalue one, or by iterating on the transition.

The steady‑state marginal value functions and distribution are then stored as vfSS and distrSS.

Steady‑state helper functions

The steady‑state solver relies on a few core helper functions that compute policies, transition objects, and the stationary distribution. Key functions include:

  • Kdiff: excess-capital function used by the root finder to locate the equilibrium aggregate capital.
  • Ksupply: given prices and value functions, solves the household problem (via EGM) and returns aggregate savings and the stationary distribution.
  • find_ss_policies: solves for steady-state policy functions (EGM/value function iteration variants depending on model type).
  • find_transitions: builds sparse transition matrices (Γ, Γa, Γn) from policy functions.
  • find_ss_distribution: computes or iterates to the ergodic joint distribution (linear or nonlinear transition methods).
  • updateW (and EGM_policyupdate): low-level routines used inside policy iteration to update marginal value / policy objects.
BASEforHANK.SteadyState.KdiffFunction
Kdiff(
    KD,
    n_par,
    m_par,
    initialize = true,
    valueFunc_guess::AbstractArray,
    distr_guess = zeros(1, 1, 1)
)

This function is used to find the stationary equilibrium of the household block of the model, in particular, it is used in find_steadystate().

Calculate the difference between the capital stock that is demanded/assumed (KD) and the capital stock that prevails under that demanded capital stock's implied prices when households face idiosyncratic income risk (Aiyagari model).

Requires global functions from the IncomesETC module and Ksupply().

Arguments

  • KD::Float64: Assumed capital demand (guess)
  • n_par::NumericalParameters, m_par::ModelParameters
  • initialize::Bool = true: If true, initialize the marginal value functions and stationary distribution, otherwise use the provided guesses that follow. Providing the guesses can be used in CustomBrent() to speed up the solution since the results from the previous iteration can be used as starting values.
  • Wb_guess::AbstractArray = zeros(1, 1, 1): Guess for marginal value of liquid assets
  • Wk_guess::AbstractArray = zeros(1, 1, 1): Guess for marginal value of illiquid assets
  • distr_guess::AbstractArray = zeros(1, 1, 1): Guess for stationary distribution

Returns

  • diff::Float64: Difference between the demanded and supplied capital stock
  • Wb::AbstractArray: Marginal value of liquid assets, implied by capital demand
  • Wk::AbstractArray: Marginal value of illiquid assets, implied by capital demand
  • distr::AbstractArray: Stationary distribution of idiosyncratic states, implied by capital demand
source
BASEforHANK.SteadyState.KsupplyFunction
Ksupply(
    args_hh_prob,
    n_par,
    m_par,
    vf,
    distr_guess,
    net_income,
    eff_int,
)

Calculate the aggregate savings when households face idiosyncratic income risk by first solving the household problem using the endogenous grid method (EGM) and then finding the stationary distribution of idiosyncratic states.

Idiosyncratic state is tuple (b,k,y), where:

  1. b: liquid assets,
  2. k: illiquid assets,
  3. y: labor income.

This function is used in find_steadystate() to find the stationary equilibrium, as an input to Kdiff(), and in BASEforHANK.PerturbationSolution.prepare_linearization() to prepare the linearization of the model.

Arguments

  • args_hh_prob: Vector of arguments to the household problem
  • n_par: Numerical parameters
  • m_par: Model parameters
  • vf: Initial guess for marginal value functions (ValueFunctions struct)
  • distr_guess: Initial guess for stationary distribution
  • net_income: Net incomes, output of functions from the IncomesETC module
  • eff_int: Effective interest rates, output of functions from the IncomesETC module

Returns

  • K: Aggregate saving in illiquid assets
  • B: Aggregate saving in liquid assets
  • transition_matrices: Transition matrices struct containing sparse transition matrices
  • pf: Optimal policy functions (PolicyFunctions struct)
  • vf: Updated marginal value functions (ValueFunctions struct)
  • distr: Ergodic stationary distribution
source
BASEforHANK.SteadyState.find_ss_policiesFunction
find_ss_policies(args_hh_prob, n_par, m_par, vf::ValueFunctionsCompleteMarkets, net_income, eff_int)

Find steady-state policy functions for the complete markets model.

For complete markets, no policy iteration is needed.

Arguments

  • args_hh_prob: Vector of arguments to the household problem (unused for complete markets)
  • n_par: Numerical parameters
  • m_par: Model parameters
  • vf: Initial value functions (unused for complete markets)
  • net_income: Net incomes (unused for complete markets)
  • eff_int: Effective interest rates

Returns

  • pf: Empty policy functions struct for complete markets
  • vf: Value functions struct with effective interest rates
source
find_ss_policies(args_hh_prob, n_par, m_par, vf::ValueFunctionsOneAsset, net_income, eff_int)

Find steady-state policy functions for the one-asset model using value function iteration.

Solves the household problem for agents with only liquid assets using the endogenous grid method (EGM) until the marginal value functions converge.

Arguments

  • args_hh_prob: Vector of arguments to the household problem
  • n_par: Numerical parameters
  • m_par: Model parameters
  • vf: Initial guess for marginal value functions
  • net_income: Net incomes from the IncomesETC module
  • eff_int: Effective interest rates

Returns

  • pf: Optimal policy functions for the one-asset model
  • vf: Converged marginal value functions
source
find_ss_policies(args_hh_prob, n_par, m_par, vf::ValueFunctionsTwoAssets, net_income, eff_int)

Find steady-state policy functions for the two-asset model using value function iteration.

Solves the household problem for agents with both liquid and illiquid assets using the endogenous grid method (EGM) until the marginal value functions converge. Handles both adjustment and non-adjustment cases for illiquid assets.

Arguments

  • args_hh_prob: Vector of arguments to the household problem
  • n_par: Numerical parameters
  • m_par: Model parameters
  • vf: Initial guess for marginal value functions (liquid and illiquid assets)
  • net_income: Net incomes from the IncomesETC module
  • eff_int: Effective interest rates

Returns

  • pf: Optimal policy functions for the two-asset model
  • vf: Converged marginal value functions for both assets
source
BASEforHANK.SteadyState.find_transitionsFunction
find_transitions(pf::PolicyFunctionsTwoAssets, n_par, m_par)

Construct transition matrices for the two-asset model based on optimal policy functions.

Creates sparse transition matrices for both adjustment and non-adjustment cases of illiquid assets, then combines them using the adjustment probability λ.

Arguments

  • pf: Policy functions for the two-asset model
  • n_par: Numerical parameters
  • m_par: Model parameters (including adjustment probability λ)

Returns

  • TransitionMatricesTwoAssets: Struct containing joint transition matrix Γ, adjustment transition matrix Γa, and non-adjustment transition matrix Γn
source
find_transitions(pf::PolicyFunctionsOneAsset, n_par, m_par)

Construct transition matrix for the one-asset model based on optimal policy functions.

Creates a sparse transition matrix based on the liquid asset policy functions.

Arguments

  • pf: Policy functions for the one-asset model
  • n_par: Numerical parameters
  • m_par: Model parameters

Returns

  • TransitionMatricesOneAsset: Struct containing the transition matrix Γ
source
BASEforHANK.SteadyState.find_ss_distributionFunction
find_ss_distribution(pf, distr_guess, m_par, n_par)

Find the stationary distribution for one-asset or two-asset models.

Computes the steady-state distribution by finding the left-hand unit eigenvector of the transition matrix.

Arguments

  • pf: Policy functions (OneAsset or TwoAssets)
  • distr_guess: Initial guess for the distribution
  • m_par: Model parameters
  • n_par: Numerical parameters

Returns

  • distr: Stationary distribution normalized to sum to 1
  • transition_matrices: Transition matrices struct
source
find_ss_distribution_splines( m_n_star::AbstractArray, cdf_guess_intial::AbstractArray,
    n_par::NumericalParameters, mbar::AbstractArray )

Find steady state disribution given steady state policies and guess for the interest
    rate through iteration. Iteration of distribution is done using monotonic spline
    interpolation to bring the next periods cdf's back to the reference grid.

# Arguments
- `m_n_star::AbstractArray`: optimal savings choice given fixed grid of assets and
  income shocks
- `cdf_guess_intial::AbstractArray`: initial guess for stationary cumulative joint
  distribution (in m) `cdf_guess[m,y]
- `n_par::NumericalParameters`


# Returns
- `pdf_ss`: stationary distribution over m and y
source
find_ss_distribution(pf::PolicyFunctionsCompleteMarkets, distr_guess, m_par, n_par)

Find the stationary distribution for the complete markets model. Set to income state PDF.

Arguments

  • pf: Policy functions for complete markets (unused)
  • distr_guess: Initial guess for distribution (unused)
  • m_par: Model parameters (unused)
  • n_par: Numerical parameters (contains income transition matrix Π)

Returns

  • distr: Stationary distribution of the income process
  • transition_matrices: Empty transition matrices struct for complete markets
source
BASEforHANK.SteadyState.updateWFunction
updateW(
    Evf::Array,
    pf,
    args_hh_prob,
    m_par,
    n_par,
    beta_factor::Number
)

Update the marginal values for the liquid and illiquid assets, given the current continuation value for the illiquid asset EWkPrime, the policy functions [x_a_star, x_n_star, b_n_star], and today's prices [args_hh_prob].

Consult the document 'ComputationalNotes.md', specifically section 3 with the title 'Update the continuation values (CV1) and (CV2) using the Envelope conditions' for an explanation of the function's code.

Arguments

  • Evf: Expected continuation value functions
  • pf: Policy functions
  • args_hh_prob: Vector of household problem arguments
  • m_par: Model parameters
  • n_par: Numerical parameters
  • beta_factor: Discount factor shock multiplier

Returns

  • vf: Updated value functions
source
BASEforHANK.SteadyState.EGM_policyupdateFunction
EGM_policyupdate(
    vfPrime::ValueFunctions,
    args_hh_prob::Vector,
    net_income::Array,
    n_par,
    m_par,
    warnme::Bool,
    model::Union{CompleteMarkets,OneAsset,TwoAsset},
)

Find optimal policies, given tomorrow's marginal continuation values in vfPrime, today's prices [args_hh_prob], and income [net_income], using the Endogenous Grid Method.

Optimal policies are defined over the exogenously fixed grid, while values of optimal policies (b and k) can have off-grid values. Please refer to the subsection with the title 'Update the optimal policy functions' of the document 'Computational Notes.md', for a detailed explanation of the function's code. The FOC's mentioned in the code are the Euler Equations as in the 'documentation of the household problem'.

Arguments

  • vfPrime: Tomorrow's marginal value functions (ValueFunctions struct)
  • args_hh_prob: Vector of arguments to the household problem
  • net_income: Incomes, output of functions from the IncomesETC module
  • n_par: Numerical parameters
  • m_par: Model parameters
  • warnme: If true, warns about non-monotonicity in resources or liquid asset choices
  • model: Model type, either CompleteMarkets, OneAsset, or TwoAsset

Returns

  • pf: Optimal policy functions struct containing all policy arrays
source

These functions live in src/SubModules/SteadyState/IM_fcns/ and src/SubModules/SteadyState/EGM/ and together implement the policy solution, transition construction, and stationary-distribution routines used by find_steadystate.

call_prepare_linearization and dimensionality reduction

BASEforHANK.PerturbationSolution.prepare_linearizationFunction
prepare_linearization(K, Wb, Wk, distr, n_par, m_par)

Given the stationary equilibrium of the household side, computed in find_steadystate(), this function performs several steps:

  • Step 1: compute the stationary equilibrium.
  • Step 2: perform the dimensionality reduction of the marginal value functions as well as the distribution.
  • Step 3: compute the aggregate steady state from input_aggregate_steady_state.mod.
  • Step 4: produce indexes to access the variables in the linearized model.
  • Step 5: return the results.

Arguments

  • K::Float64: steady-state capital stock
  • Wb::Array{Float64,3}, Wk::Array{Float64,3}: steady-state marginal value functions
  • distr::Array{Float64,3}: steady-state distribution of idiosyncratic states
  • n_par::NumericalParameters,m_par::ModelParameters

Returns

  • XSS::Array{Float64,1}, XSSaggr::Array{Float64,1}: steady state vectors produced by @writeXSS()
  • indexes, indexes_aggr: structs for accessing XSS,XSSaggr by variable names, produced by @make_fn(), @make_fnaggr()
  • compressionIndexes::Array{Array{Int,1},1}: indexes for compressed marginal value functions (V_m and V_k)
  • n_par::NumericalParameters, m_par::ModelParameters: updated parameters
  • CDFSS, CDF_bSS, CDF_kSS, CDF_hSS: cumulative distribution functions (joint and marginals)
  • distrSS::Array{Float64,3}: steady state distribution of idiosyncratic states, computed by Ksupply()
source

call_prepare_linearization(ss, m_par) takes the SteadyStateStruct coming from call_find_steadystate and prepares all objects needed for the perturbation solution. Internally it calls BASEforHANK.PerturbationSolution.prepare_linearization (from the PerturbationSolution submodule) with

  • the steady‑state capital stock and prices,
  • the steady‑state marginal value functions vfSS,
  • the stationary distribution distrSS,
  • the numerical parameters n_par, and
  • the model parameters m_par.

The main steps are:

Compute additional steady‑state objects and distributional statistics

Using BASEforHANK.SteadyState.distrSummaries, the code computes a set of distributional statistics such as top‑10% wealth and income shares, Gini coefficients, and the standard deviation of log labor earnings.

First‑stage DCT reduction of marginal value functions

The function BASEforHANK.SteadyState.first_stage_reduction selects DCT coefficients that summarize fluctuations of V_b and V_k around the steady state:

  • It computes derivatives of the marginal value functions with respect to all prices that enter the household problem.
  • It transforms these derivatives using a DCT (FFTW.jl) into polynomial‑like coefficients.
  • Based on n_par.reduc_marginal_value and n_par.reduc_value, it keeps the coefficients that explain most of the “energy” (variance) in those coefficients and in the levels of V_b and V_k.

The selected coefficient indices are stored in compressionIndexes[1].

Copula mesh and node selection

For the copula representation, the preprocessor constructs a grid over (b,k,h) such that each bin contains (approximately) the same mass of the respective variable. The corresponding indices in state space are stored in compressionIndexes[2].

Construction of state and control vectors

Finally, prepare_linearization collects all steady‑state values in the vector XSS and its aggregate counterpart XSSaggr. It also builds the index structs indexes and indexes_aggr (using the parsing utilities in BASEforHANK.Parsing) that map named model variables to positions in the state and control vectors and allow keeping track of steady‑state levels and deviations separately.

  • State variables consist of marginal distributions over b, k and h and aggregate state variables (listed in state_names).
  • Control variables consist of the steady‑state marginal value functions (stored on the full grid) as well as aggregate controls (control_names).

Deviations from the steady state are stored in reduced form:

  • only the selected DCT coefficients are kept for marginal value functions,
  • deviations of marginal distributions have one element fewer than the grid, because they must sum to one.

The result is returned as a SteadyResults struct:

BASEforHANK.Parsing.SteadyStateStructType
SteadyStateStruct

Container for steady-state solution components.

Fields

  • KSS: steady-state aggregate capital (or placeholder type).
  • vfSS: value functions at steady state (ValueFunctions).
  • distrSS: steady-state distribution object (type varies by model).
  • n_par: numeric parameters related to solution dimensions.
source
BASEforHANK.Parsing.SteadyResultsType
SteadyResults

Holds outputs produced when computing steady-state results.

Fields

  • XSS: full steady-state variable collection.
  • XSSaggr: aggregate steady-state variables.
  • indexes, indexes_r, indexes_aggr: index maps for states/controls.
  • compressionIndexes: indexes used when compressing state space.
  • n_par, m_par: numeric parameters and model dimensions.
  • distrSS: steady-state distribution values (DistributionValues).
  • state_names, control_names: lists of state/control variable names.
source

Distributional summaries and first‑stage reduction

BASEforHANK.SteadyState.distrSummariesFunction
distrSummaries(distr, q, pf, n_par, net_income, gross_income, m_par)

Compute distributional summary statistics for a household distribution.

This file provides two dispatched implementations:

  • distrSummaries(distr::RepAgent, ...) — placeholder implementation that currently returns small numeric placeholders (eps()).
  • distrSummaries(distr::Union{CDF,CopulaOneAsset,CopulaTwoAssets}, ...) — concrete implementation that returns the tuple (TOP10Wshare, TOP10Ishare, TOP10Inetshare, giniwealth, giniconsumption, sdlogy).

Arguments

  • distr: distribution object (RepAgent, CDF, or Copula*).
  • q::Real: price of illiquid asset (used to compute total wealth when applicable).
  • pf::PolicyFunctions: policy functions (used to compute consumption summaries).
  • n_par::NumericalParameters: numerical/grid parameters and metadata.
  • net_income, gross_income: arrays following the codebase's internal indexing convention (see distr_summaries_incomes docstring for index usage).
  • m_par::ModelParameters: model parameters container.

Returns

  • For the CDF/Copula dispatch, a 6-tuple:

    • TOP10Wshare::Float64 — top 10% wealth share.
    • TOP10Ishare::Float64 — top 10% gross income share.
    • TOP10Inetshare::Float64 — top 10% net income share.
    • giniwealth::Float64 — Gini coefficient of the wealth distribution.
    • giniconsumption::Float64 — Gini coefficient of consumption.
    • sdlogy::Float64 — standard deviation of log labor earnings.

Notes

  • The RepAgent dispatch is currently a stub and should be implemented if used.
source
BASEforHANK.SteadyState.first_stage_reductionFunction
first_stage_reduction(
    vfSS,
    transition_matricesSS,
    pfSS,
    args_hh_prob,
    include_list_idx,
    n_par,
    m_par,
)

First-stage reduction for selecting DCT (Discrete Cosine Transform) coefficients to perturb, following Appendix C of BBL. Uses the steady-state value functions vfSS, policy functions pfSS, and transition matrices to compute the Jacobian of marginal values with respect to selected aggregate inputs, transform these via the inverse marginal utility, and select a small set of DCT coefficients that best represent the response shapes.

  • include_list_idx selects which elements of args_hh_prob are perturbed; others remain at steady state.
  • Works for one-asset and two-asset models via multiple dispatch on vfSS and transition_matricesSS.
  • Complete-markets variant is a no-op selection (representative-agent shortcut) and returns an empty selection.

Returns a pair (ind, J) where ind contains selected DCT indices by asset (liquid/illiquid when applicable) and J is the Jacobian of marginal values with respect to the selected inputs.

Arguments

  • vfSS::ValueFunctions: steady-state value functions.
  • transition_matricesSS::TransitionMatrices: steady-state exogenous transition matrices.
  • pfSS::PolicyFunctions: steady-state policy functions.
  • args_hh_prob::AbstractVector{Float64}: vector of household problem arguments at steady state (in original scale).
  • include_list_idx::AbstractVector{Int}: indexes of args_hh_prob to include in the reduction step.
  • n_par::NumericalParameters: numerical parameters (provides DCT sizes, etc.).
  • m_par::ModelParameters: model parameters (provides utility parameters, etc.).

Returns

- `ind::Union{Array{Array{Int}}, Array{Int}}`: selected DCT coefficient indices by asset.
- `J::AbstractMatrix{Float64}`: Jacobian of marginal values with respect to selected inputs.
source

distrSummaries computes distributional statistics from the stationary distribution and associated policy functions and incomes. first_stage_reduction implements the DCT‑based dimensionality reduction used in the perturbation solution.

Calibration of the steady state

The SteadyState module also provides a helper for calibrating the model to match empirical moments before running the steady‑state and linearization pipeline. The calibration routine is in detail described in the section Calibration.

Parameters

The economic parameters of the model are collected in the struct ModelParameters, which also stores priors for estimation and the stochastic process parameters for aggregate shocks.

BASEforHANK.Parsing.ModelParametersType
ModelParameters()

A structure to collect all model parameters, including calibrated values and prior distributions for estimation.

Overview

  • This struct is designed for macroeconomic models and includes parameters related to household preferences, income processes, technological factors, monetary policy, fiscal policy, and exogenous shocks.
  • The parameters are annotated with metadata such as names (both ASCII and LaTeX), prior distributions, and a boolean flag indicating whether they are estimated.
  • Uses the Parameters, FieldMetadata, and Flatten packages to facilitate parameter management.

Fields

Each field follows the structure:

parameter::T = value | "ascii_name" | L"latex_name" | prior_distribution | estimated
  • parameter: Internal model parameter name.
  • value: Default numerical value.
  • ascii_name: Human-readable name used in output or logging.
  • latex_name: Corresponding LaTeX notation for use in reports and documentation.
  • prior_distribution: Prior distribution for Bayesian estimation (if applicable).
  • estimated: Boolean indicating whether the parameter is estimated.
source

The numerical parameters are collected in the struct NumericalParameters. They contain, among others:

  • the grids on which the stationary equilibrium is solved,
  • discretization objects (e.g. the income transition matrix and joint distribution),
  • parameters governing numerical solution accuracy (e.g. ϵ, reduc_value, reduc_marginal_value, further_compress),
  • options for the steady-state search (rKmin, rKmax, coarse-grid settings), and
  • algorithmic choices such as sol_algo.

Some of the most important switches are:

  • model: chooses CompleteMarkets, OneAsset, or TwoAsset, which determines how many assets the household problem features and which policy/value functions are allocated inside n_par.
  • transition_type: selects LinearTransition versus NonLinearTransition, i.e. whether distributions evolve linearly around steady state, using the Young-method, or through the nonlinear law of motion using the DEGM-method.[BLWW]
  • distribution_states: toggles between CopulaStates and CDFStates (or RepAgentIndexes in complete-markets mode), thereby deciding whether linearization is performed on copula coefficients or cumulative distributions.

These knobs live in n_par because they control the shape of the solution space. See examples/baseline/Model/input_parameters.jl for default choices and NumericalParameters for the full keyword list.

In particular, nh, nk, and nb control the resolution for the income, illiquid‑asset, and liquid‑asset grids. The resolution of the copula used in the linearization does not need to coincide with these grids and is controlled by nh_copula, nk_copula, and nb_copula. The copula resolution should not exceed the underlying grid size.

Aggregation functions

Quick reference for computing aggregates and population shares from a distribution object (implementations in src/SubModules/SteadyState/IM_fcns/fcn_aggregation.jl).

Main helpers

  • aggregate_asset(...) — compute the aggregate of an asset. Accepts either a marginal vector (PDF or CDF) plus the asset grid and transition_type, or a distribution container (CDF, Copula*, RepAgent) together with n_par (which provides grid_<asset> and transition_type).
  • eval_cdf(distr, asset, n_par; bound) — returns the cumulative mass up to bound for asset (useful for population shares).
BASEforHANK.SteadyState.integrate_assetFunction
integrate_asset(cdf::AbstractVector, grid::AbstractVector, bound)

Compute aggregate asset holdings E[K] = ∫ k f(k) dk from a savings CDF using integration by parts.

Arguments:

  • cdf::AbstractVector: cumulative distribution F(k) defined on grid
  • grid::AbstractVector: increasing asset grid
  • bound: optional upper integration bound; if nothing integrates over the full grid. Methods exist for bound::Nothing and bound::Real (or Float64).

Details:

  • Uses an interpolator for the CDF and computes E[K] = k*F(k) |{kmin}^{bound} - ∫{kmin}^{bound} F(k) dk
source
BASEforHANK.SteadyState.aggregateFunction
aggregate(pdf::AbstractVector, grid::AbstractVector, bound)

Compute aggregate asset holdings by integrating over the PDF.

Arguments:

  • pdf::AbstractVector: Probability density function over asset grid
  • grid::AbstractVector: Asset grid
  • bound: Optional upper bound for integration; if nothing integrates over full grid.

Returns:

  • Aggregate asset holdings (dot product over grid or truncated grid)
source
BASEforHANK.SteadyState.aggregate_asset_helperFunction
aggregate_asset_helper(distr, grid, transition_type; bound=nothing, pdf_input=true)

Compute aggregate asset holdings from a marginal distribution for either linear or nonlinear transition schemes.

Arguments:

  • distr::AbstractVector: marginal distribution values (PDF if pdf_input=true, CDF if pdf_input=false)

  • grid::AbstractVector: asset grid

  • transition_type::LinearTransition | NonLinearTransition: determines aggregation method

    • LinearTransition: uses pointwise aggregation over the grid.

      • If distr is a PDF, calls aggregate(distr, grid, bound).
      • If distr is a CDF, converts to a PDF using cdf_to_pdf(distr) and calls aggregate.
    • NonLinearTransition: uses integration-by-parts on the CDF.

      • If distr is a PDF, first converts to a CDF via cumsum(distr) and calls integrate_asset.
      • If distr is a CDF, calls integrate_asset directly.
  • bound: optional upper integration bound (default nothing)

  • pdf_input::Bool: true if distr is a PDF, false if distr is a CDF

Returns:

  • Aggregate asset holdings (scalar)
source
BASEforHANK.SteadyState.aggregate_assetFunction
aggregate_asset(distr, asset, n_par, bound=nothing)

Aggregate asset holdings across the distribution of agents.

Arguments

  • distr::Union{CDF, CopulaOneAsset, CopulaTwoAssets, RepAgent}: The distribution of agents
  • asset::Symbol: The asset to aggregate (:b for bonds, :k for capital, :h for productivity)
  • n_par::NumericalParameters: Numerical parameters containing asset grids (uses grid_<asset> and transition_type)
  • bound: optional upper integration bound (default nothing = full grid)

Returns

  • The aggregate level of the specified asset across all agents.
source
BASEforHANK.SteadyState.aggregate_B_KFunction
aggregate_B_K(distr, m_par, n_par, model)

Calculate aggregate asset holdings depending on the model type.

Arguments

  • distr: Distribution (unused for complete markets)
  • m_par::ModelParameters: Model parameters (contains β for interest rate calculation, ψ for portfolio share)
  • n_par::NumericalParameters: Numerical parameters (grids etc.)
  • model::AbstractMacroModel: Model type (CompleteMarkets, OneAsset, TwoAsset)

Returns

  • (B, K): Aggregate bonds and aggregate capital (ordering as returned by each method)
source
BASEforHANK.SteadyState.get_asset_distrFunction
get_asset_distr(distr, asset)

Extract the marginal distribution for a specific asset from various distribution types.

Arguments:

  • distr::Union{CDF, CopulaOneAsset, CopulaTwoAssets, RepAgent}: The distribution of agents
  • asset::Symbol: The asset to extract (:b, :k, or :h)

Returns:

  • Vector representing the marginal distribution for the specified asset.
source
BASEforHANK.SteadyState.eval_cdfFunction
eval_cdf(distr, asset, n_par, bound=nothing)

Evaluate the cumulative distribution function given either a PDF or CDF distribution object for a specific asset.

Arguments:

  • distr::DistributionValues (or related Copula/CDF types): distribution of agents
  • asset::Symbol: asset to evaluate (:b, :k, or :h)
  • n_par::NumericalParameters: numerical parameters with grids
  • bound: optional upper bound; if nothing evaluates over full grid

Returns:

  • Value of the CDF up to bound (uses PDF->CDF conversion if input is a PDF)
source

References