SimSSMIVary {simStateSpace}R Documentation

Simulate Data from a State Space Model (Individual-Varying Parameters)

Description

This function simulates data using a state space model. It assumes that the parameters can vary across individuals.

Usage

SimSSMIVary(
  n,
  time,
  delta_t = 1,
  mu0,
  sigma0_l,
  alpha,
  beta,
  psi_l,
  nu,
  lambda,
  theta_l,
  type = 0,
  x = NULL,
  gamma = NULL,
  kappa = NULL
)

Arguments

n

Positive integer. Number of individuals.

time

Positive integer. Number of time points.

delta_t

Numeric. Time interval. The default value is 1.0 with an option to use a numeric value for the discretized state space model parameterization of the linear stochastic differential equation model.

mu0

List of numeric vectors. Each element of the list is the mean of initial latent variable values (\boldsymbol{\mu}_{\boldsymbol{\eta} \mid 0}).

sigma0_l

List of numeric matrices. Each element of the list is the Cholesky factorization (t(chol(sigma0))) of the covariance matrix of initial latent variable values (\boldsymbol{\Sigma}_{\boldsymbol{\eta} \mid 0}).

alpha

List of numeric vectors. Each element of the list is the vector of constant values for the dynamic model (\boldsymbol{\alpha}).

beta

List of numeric matrices. Each element of the list is the transition matrix relating the values of the latent variables at the previous to the current time point (\boldsymbol{\beta}).

psi_l

List of numeric matrices. Each element of the list is the Cholesky factorization (t(chol(psi))) of the covariance matrix of the process noise (\boldsymbol{\Psi}).

nu

List of numeric vectors. Each element of the list is the vector of intercept values for the measurement model (\boldsymbol{\nu}).

lambda

List of numeric matrices. Each element of the list is the factor loading matrix linking the latent variables to the observed variables (\boldsymbol{\Lambda}).

theta_l

List of numeric matrices. Each element of the list is the Cholesky factorization (t(chol(theta))) of the covariance matrix of the measurement error (\boldsymbol{\Theta}).

type

Integer. State space model type. See Details in SimSSMFixed() for more information.

x

List. Each element of the list is a matrix of covariates for each individual i in n. The number of columns in each matrix should be equal to time.

gamma

List of numeric matrices. Each element of the list is the matrix linking the covariates to the latent variables at current time point (\boldsymbol{\Gamma}).

kappa

List of numeric matrices. Each element of the list is the matrix linking the covariates to the observed variables at current time point (\boldsymbol{\kappa}).

Details

Parameters can vary across individuals by providing a list of parameter values. If the length of any of the parameters (mu0, sigma0_l, alpha, beta, psi_l, nu, lambda, theta_l, gamma, or kappa) is less the n, the function will cycle through the available values.

Value

Returns an object of class simstatespace which is a list with the following elements:

Author(s)

Ivan Jacob Agaloos Pesigan

References

Chow, S.-M., Ho, M. R., Hamaker, E. L., & Dolan, C. V. (2010). Equivalence and differences between structural equation modeling and state-space modeling techniques. Structural Equation Modeling: A Multidisciplinary Journal, 17(2), 303–332. doi:10.1080/10705511003661553

See Also

Other Simulation of State Space Models Data Functions: LinSDE2SSM(), SimSSMFixed(), SimSSMLinGrowth(), SimSSMLinGrowthIVary(), SimSSMLinSDEFixed(), SimSSMLinSDEIVary(), SimSSMOUFixed(), SimSSMOUIVary(), SimSSMVARFixed(), SimSSMVARIVary()

Examples

# prepare parameters
# In this example, beta varies across individuals.
set.seed(42)
## number of individuals
n <- 5
## time points
time <- 50
## dynamic structure
p <- 3
mu0 <- list(
  rep(x = 0, times = p)
)
sigma0 <- 0.001 * diag(p)
sigma0_l <- list(
  t(chol(sigma0))
)
alpha <- list(
  rep(x = 0, times = p)
)
beta <- list(
  0.1 * diag(p),
  0.2 * diag(p),
  0.3 * diag(p),
  0.4 * diag(p),
  0.5 * diag(p)
)
psi <- 0.001 * diag(p)
psi_l <- list(
  t(chol(psi))
)
## measurement model
k <- 3
nu <- list(
  rep(x = 0, times = k)
)
lambda <- list(
  diag(k)
)
theta <- 0.001 * diag(k)
theta_l <- list(
  t(chol(theta))
)
## covariates
j <- 2
x <- lapply(
  X = seq_len(n),
  FUN = function(i) {
    matrix(
      data = stats::rnorm(n = time * j),
      nrow = j,
      ncol = time
    )
  }
)
gamma <- list(
  diag(x = 0.10, nrow = p, ncol = j)
)
kappa <- list(
  diag(x = 0.10, nrow = k, ncol = j)
)

# Type 0
ssm <- SimSSMIVary(
  n = n,
  time = time,
  mu0 = mu0,
  sigma0_l = sigma0_l,
  alpha = alpha,
  beta = beta,
  psi_l = psi_l,
  nu = nu,
  lambda = lambda,
  theta_l = theta_l,
  type = 0
)

plot(ssm)

# Type 1
ssm <- SimSSMIVary(
  n = n,
  time = time,
  mu0 = mu0,
  sigma0_l = sigma0_l,
  alpha = alpha,
  beta = beta,
  psi_l = psi_l,
  nu = nu,
  lambda = lambda,
  theta_l = theta_l,
  type = 1,
  x = x,
  gamma = gamma
)

plot(ssm)

# Type 2
ssm <- SimSSMIVary(
  n = n,
  time = time,
  mu0 = mu0,
  sigma0_l = sigma0_l,
  alpha = alpha,
  beta = beta,
  psi_l = psi_l,
  nu = nu,
  lambda = lambda,
  theta_l = theta_l,
  type = 2,
  x = x,
  gamma = gamma,
  kappa = kappa
)

plot(ssm)


[Package simStateSpace version 1.2.1 Index]