DynamicFactorModeling.jl

Documentation for DynamicFactorModeling.jl

DynamicFactorModeling.DFMMeansType
DFMMeans(F::Array{Float64}, B::Array{Float64}, S::Array{Float64}, P::Array{Float64}, P2::Array{Float64})

Description: HDFM Bayesian estimator-generated latent factor and hyperparameter sample means (expected values).

Inputs:

  • F = MCMC-generated latent factor sample mean.
  • B = MCMC-generated observation equation regression coefficient sample means.
  • S = MCMC-generated observable variable idiosyncratic error disturbance variance sample means.
  • P = MCMC-generated latent factor autoregressive coefficient sample means.
  • P2 = MCMC-generated idiosyncratic error autoregressive coefficient sample means.
source
DynamicFactorModeling.DFMResultsType
DFMResults(F::Array{Float64}, B::Array{Float64}, S::Array{Float64}, P::Array{Float64}, P2::Array{Float64}, means::DFMMeans)

Description: HDMF Bayesian estimator-generated MCMC posterior distribution samples and their means for latent factors and hyperparameters.

Inputs:

  • F = MCMC-generated latent factor sample.
  • B = MCMC-generated observation equation regression coefficient sample.
  • S = MCMC-generated observable variable idiosyncratic error disturbance variance sample.
  • P = MCMC-generated latent factor autoregressive coefficient sample.
  • P2 = MCMC-generated idiosyncratic error autoregressive coefficient sample.
  • means = HDFM Bayesian estimator-generated latent factor and hyperparameter sample means (expected values).
source
DynamicFactorModeling.DFMStructType
DFMStruct(factorlags::Int64, errorlags::Int64, ndraws::Int64, burnin::Int64)

Description: 1-level DFM lag structure specification and MCMC sample size for Bayesian estimation.

Inputs:

  • factorlags = Number of lags in the autoregressive specification of the latent factors.
  • errorlags = Number of lags in the autoregressive specification of the observable variable idiosyncratic errors.
  • ndraws = Number of MCMC draws used for posterior distributions.
  • burnin = Number of initial MCMC draws discarded.
source
DynamicFactorModeling.HDFMType
HDFM(
    nlevels::Int64                   
    nvar::Int64                     
    nfactors::Array{Int64,1}        
    fassign::Array{Int64,2}          
    flags::Array{Int64,1}         
    varlags::Array{Int64,1}        
    varcoefs::Array{Any,2}          
    varlagcoefs::Array{Any,2}    
    fcoefs::Array{Any,1}           
    fvars::Array{Any,1}             
    varvars::Array{Any,1}  
)

Creates an object of type HDFM that contains all parameters necessary to specify a multi-level linear dynamic factor data-generating process. This is a convenient alternative to specifying an HDFM directly in state-space form.

Inputs:

  • nlevels = number of levels in the multi-level model structure.
  • nvar = number of variables.
  • nfactors = number of factors for each level (vector of length nlevels).
  • fassign = determines which factor is assigned to which variable for each level (integer matrix of size nvar × nlevels).
  • flags = number of autoregressive lags for factors of each level (factors of the same level are restricted to having the same number of lags; vector of length nlevels).
  • varlags = number of observed variable error autoregressive lags (vector of length nvar).
  • varcoefs = vector of coefficients for each variable in the observation equation (length 1+nlevels, where first entry represents the intercept).
  • fcoefs = list of nlevels number of matrices, for which each row contains vectors of the autoregressive lag coefficients of the corresponding factor.
  • fvars = list of nlevels number of vectors, where each entry contains the disturbance variance of the corresponding factors.
  • varvars = vector of nvar number of entries, where each entry contains the innovation variance of the corresponding variable.
source
DynamicFactorModeling.HDFMStructType
HDFMStruct(nlevels::Int64, nfactors::Array{Int64,1}, factorassign::Array{Int64,2}, factorlags::Array{Int64,1}, errorlags::Array{Int64,1}, ndraws::Int64, burnin::Int64)

Description: Multi-level/hierarchical DFM (HDFM) level, factor assignment, and lag structure specification, and MCMC sample size for Bayesian estimation.

Inputs:

  • nlevels = Number of levels in the HDFM specification.
  • nvars = Number of observable variables in the HDFM specification.
  • nfactors = Number of factor per level in the HDFM specification.
  • factorassign = Factors assigned to each variable across all levels.
  • factorlags = Number of lags in the autoregressive specification of the latent factors.
  • errorlags = Number of lags in the autoregressive specification of the observable variable idiosyncratic errors.
  • ndraws = Number of MCMC draws used for posterior distributions.
  • burnin = Number of initial MCMC draws discarded.
source
DynamicFactorModeling.SSModelType

SSModel( H::Array{Float64,2}, A::Array{Float64,2}, F::Array{Float64,2}, μ::Array{Float64,1}, R::Array{Float64,2}, Q::Array{Float64,2}, Z::Array{Float64,2} )

A type object containing all parameters necessary to specify a data-generating process in state-space form. Measurement Equation:

  • y{t} = H β{t} + A z{t} + e{t}

Transition Equation:

  • β{t} = μ + F β{t-1} + v_{t}
  • e_{t} ~ i.i.d.N(0,R)
  • v_{t} ~ i.i.d.N(0,Q)
  • z_{t} ~ i.i.d.N(0,Z)
  • E(e{t} v{s}') = 0

Inputs:

  • H = measurement equation state vector coefficient matrix.
  • A = measurement equation predetermined vector coefficient matrix.
  • F = state equation companion matrix.
  • μ = state equation intercept vector.
  • R = measurement equation error covariance matrix.
  • Q = state equation innovation covariance matrix.
  • Z = predetermined vector covariance matrix.
source
DynamicFactorModeling.KN1LevelEstimatorMethod
KNSingleFactorEstimator(data::Array{Float64,2}, dfm::DFMStruct)

Description: Estimate a single-factor DFM using the Kim-Nelson approach.

Inputs:

  • data = Matrix with each column being a data series.
  • dfm = Model structure specification.

Outputs:

  • results = HDMF Bayesian estimator-generated MCMC posterior distribution samples and their means for latent factors and hyperparameters.
source
DynamicFactorModeling.KN2LevelEstimatorMethod
KN2LevelEstimator(data::Array{Float64,2}, hdfm::HDFMStruct)

Description: Estimate a two-level HDFM using the Kim-Nelson approach. Both the latent factors and hyperparameters are estimated using the Bayesian approach outlined in Kim and Nelson (1999).

Inputs:

  • data = Matrix with each column being a data series.
  • hdfm = Model structure specification.

Outputs:

  • results = HDMF Bayesian estimator-generated MCMC posterior distribution samples and their means for latent factors and hyperparameters.
source
DynamicFactorModeling.KNFactorSamplerMethod
KNFactorSampler(data_y, ssmodel)

Description: Draw a sample series of dynamic factor from conditional distribution in Ch 8, Kim & Nelson (1999). Measurement Equation: y{t} = H{t} β{t} + A z{t} + e{t}. Transition Equation: β{t} = μ + F β{t-1} + v{t}; e{t} ~ i.i.d.N(0,R); v{t} ~ i.i.d.N(0,Q); z{t} ~ i.i.d.N(0,Z); E(et v_s') = 0.

Inputs:

  • data = observed data
  • H = measurement eq. state coef. matrix
  • A = measurement eq. exogenous coef. matrix
  • F = state eq. companion matrix
  • μ = state eq. intercept term
  • R = covariance matrix on measurement disturbance
  • Q = covariance matrix on state disturbance
  • Z = covariance matrix on predetermined var vector
source
DynamicFactorModeling.createSSforHDFMMethod
createSSforHDFM(hdfm::HDFM))

Description: Create state-space form coefficient and variance matrices for an HDFM object. Measurement Equation:

  • y{t} = H β{t} + A z{t} + e{t}

Transition Equation:

  • β{t} = μ + F β{t-1} + v_{t}
  • e_{t} ~ i.i.d.N(0,R)
  • v_{t} ~ i.i.d.N(0,Q)
  • z_{t} ~ i.i.d.N(0,Z)
  • E(e{t} v{s}') = 0

Inputs:

  • hdfm::HDFM

Output:

  • H = measurement equation state vector coefficient matrix.
  • A = measurement equation predetermined vector coefficient matrix.
  • F = state equation companion matrix.
  • μ = state equation intercept vector.
  • R = measurement equation error covariance matrix.
  • Q = state equation innovation covariance matrix.
  • Z = predetermined vector covariance matrix.
source
DynamicFactorModeling.kalmanFilterMethod
kalmanFilter(data, ssmodel)

Description: Apply Kalman filter to observed data. Measurement Equation: y{t} = H{t} β{t} + A z{t} + e{t} . Transition Equation: β{t} = μ + F β{t-1} + v{t}; e{t} ~ i.i.d.N(0,R); v{t} ~ i.i.d.N(0,Q); z{t} ~ i.i.d.N(0,Z); E(et v_s') = 0.

Inputs:

  • data = observed data
  • H = measurement eq. state coef. matrix
  • A = measurement eq. exogenous coef. matrix
  • F = state eq. companion matrix
  • μ = state eq. intercept term
  • R = covariance matrix on measurement disturbance
  • Q = covariance matrix on state disturbance
  • Z = covariance matrix on predetermined var vector
source
DynamicFactorModeling.kalmanSmootherMethod
kalmanSmoother(data, ssmodel)

Description: Apply Kalman smoother to observed data. Measurement Equation: y{t} = H{t} β{t} + A z{t} + e{t}. Transition Equation: β{t} = μ + F β{t-1} + v{t}; e{t} ~ i.i.d.N(0,R); v{t} ~ i.i.d.N(0,Q); z{t} ~ i.i.d.N(0,Z); E(et v_s') = 0.

Inputs:

  • data = observed data
  • H = measurement eq. state coef. matrix
  • A = measurement eq. exogenous coef. matrix
  • F = state eq. companion matrix
  • μ = state eq. intercept term
  • R = covariance matrix on measurement disturbance
  • Q = covariance matrix on state disturbance
  • Z = covariance matrix on predetermined var vector
source
DynamicFactorModeling.mvnMethod
mvn(μ, Σ, n::Int64)

Draw n number of observations from a multivariate normal distribution with mean vector μ and covariance matrix Σ. Use cholesky decomposition to generate n draws of X = Z Q + μ, where Z is (d × 1) N(0,1) vector, and Q is upper-triangular cholesky matrix. Cov. matrix Σ does not require non-degenerate random variables (nonzero diag.).

Inputs:

  • μ = mean vector
  • Σ = covariance matrix
  • n = number of draws

Output:

  • X::Array{Float64, 2} = simulated data matrix composed of n-number of draws of X ~ N(μ,Σ)
source
DynamicFactorModeling.mvnMethod
mvn(μ, Σ)

Draw from a multivariate normal distribution with mean vector μ and covariance matrix Σ. Use cholesky decomposition to generate X = Z Q + μ, where Z is (d × 1) N(0,1) vector, and Q is upper-triangular cholesky matrix. Cov. matrix Σ does not require non-degenerate random variables (nonzero diagonal).

Inputs:

  • μ = mean vector
  • Σ = covariance matrix

Outputs:

  • X::Array{Float64, 1} = observed draw of X ~ N(μ,Σ)
source
DynamicFactorModeling.simulateSSModelMethod
simulateSSModel(num_obs::Int64, ssmodel::SSModel)

Generate data from a DGP in state space form. Measurement Equation: y{t} = H β{t} + A z{t} + e{t} Transition Equation: β{t} = μ + F β{t-1} + v{t} e{t} ~ i.i.d.N(0,R) v{t} ~ i.i.d.N(0,Q) z{t} ~ i.i.d.N(0,Z) E(et vs') = 0

Inputs:

  • num_obs = number of observations
  • ssmodel::SSModel

Output:

  • data_y = simulated sample of observed vector
  • data_z = simulated sample of exogenous variables
  • data_β = simulated sample of state vector
source
DynamicFactorModeling.vardecomp2levelMethod
vardecomp2level(data::Array{Float64, 2}, factor::Array{Float64}, betas::Array{Float64}, factorassign::Array{Float64})

Description: Compute the portion of the variation of each observable series that may be attributed to their corresponding/assigned latent factors across all levels.

Inputs:

  • data = Matrix with each column representing a data series.
  • factor = Matrix containing latent factor estimates.
  • betas = Matrix containing observation equation coefficient parameter estimates.
  • factorassign = Matrix containing the indeces of factors across all levels (columns) assigned to each observable series (rows).

Output:

  • vardecomps = Matrix containing the variance contributions of factors across all levels (columns) corresponding to each observable series (rows).
source
DynamicFactorModeling.ΓinvMethod
Γinv(T::Int64, θ::Float64)

Gamma inverse distribution with T degrees of freedom and scale parameter θ.

Inputs:

  • T = degrees of freedom
  • θ = scale parameter

Output:

  • σ2 = draw of X ~ Γ_inverse(T,θ)
source