Overview

The Stochastic Process Model (SPM) was developed several decades ago Yashin, Arbeev, Akushevich, et al. (2007), and applied for analyses of clinical, demographic, epidemiologic longitudinal data as well as in many other studies that relate stochastic dynamics of repeated measures to the probability of end-points (outcomes). SPM links the dynamic of stochastical variables with a hazard rate as a quadratic function of the state variables (Yashin, Arbeev, Akushevich, et al. 2007). The R-package, “stpm,” is a set of utilities to estimate parameters of stochastic process and modeling survival trajectories and time-to-event outcomes observed from longitudinal studies. It is a general framework for studying and modeling survival (censored) traits depending on random trajectories (stochastic paths) of variables.

Installation

Stable version from CRAN

install.packages("stpm")

Most-recent version from GitHub

require(devtools)
devtools::install_github("izhbannikov/stpm")

Data description

Data represents a typical longitudinal data in form of two datasets: longitudinal dataset (follow-up studies), in which one record represents a single observation, and vital (survival) statistics, where one record represents all information about the subject. Longitudinal dataset cat contain a subject ID (identification number), status (event(1)/censored(0)), time and measurements across the variables.

Below there is an example of clinical data that can be used in stpm and we will discuss the fields later.

Longitudinal table:

##   ID IndicatorDeath Age      DBP      BMI
## 1  1              0  30 80.00000 25.00000
## 2  1              0  32 80.51659 26.61245
## 3  1              0  34 77.78412 29.16790
## 4  1              0  36 77.86665 32.40359
## 5  1              0  38 96.55673 31.92014
## 6  1              0  40 94.48616 32.89139

Description of data fields

  • ID - subject unique identificatin number.
  • IndicatorDeath - 0/1, indicates death of a subject.
  • Age - current age of subject at observation.
  • DBP, BMI - covariates, here “DBP” represents a diastolic blood pressure, “BMI” a body-mass index.

“Short” and “Long” longitudinal data formats

The packate accepts longitudinal data in two formats: “short” and “long.”

“Short” format
##   id xi  t          y1
## 1  1  0 30   36.133394
## 2  1  0 31   -5.539881
## 3  1  0 32  -45.129493
## 4  1  0 33  -82.739624
## 5  1  0 34 -118.469249
## 6  1  0 35 -152.412392
“Long” format
##   id xi t1 t2          y1     y1.next
## 1  1  0 30 31   36.133394   -5.539881
## 2  1  0 31 32   -5.539881  -45.129493
## 3  1  0 32 33  -45.129493  -82.739624
## 4  1  0 33 34  -82.739624 -118.469249
## 5  1  0 34 35 -118.469249 -152.412392
## 6  1  0 35 36 -152.412392 -184.658379

Discrete- and continuous-time models

There are two main SPM types in the package: discrete-time model (Akushevich, Kulminski, and Manton 2005) and continuous-time model (Yashin, Arbeev, Akushevich, et al. 2007). Discrete model assumes equal intervals between follow-up observations. The example of discrete dataset is given below.

library(stpm)
data <- simdata_discr(N=10) # simulate data for 10 individuals, "long" format (default)
head(data)
##   id xi t1 t2          y1     y1.next
## 1  1  0 30 31   36.133394   -5.539881
## 2  1  0 31 32   -5.539881  -45.129493
## 3  1  0 32 33  -45.129493  -82.739624
## 4  1  0 33 34  -82.739624 -118.469249
## 5  1  0 34 35 -118.469249 -152.412392
## 6  1  0 35 36 -152.412392 -184.658379

In this case there are equal intervals between \(t_1\) and \(t_2\).

In the continuous-time SPM, in which intervals between observations are not equal (arbitrary or random). The example of such dataset is shown below:

library(stpm)
data <- simdata_cont(N=5, format="short") # simulate data for 5 individuals, "short" format
head(data)
##   id xi  t         y1
## 1  0  0 30   71.22668
## 2  0  0 31   28.86212
## 3  0  0 32  -11.43627
## 4  0  0 33  -49.76927
## 5  0  0 34  -86.23273
## 6  0  0 35 -120.91783

Discrete-time model

The discrete model assumes fixed time intervals between consecutive observations. In this model, \(\mathbf{Y}(t)\) (a \(k \times 1\) matrix of the values of covariates, where \(k\) is the number of considered covariates) and \(\mu(t, \mathbf{Y}(t))\) (the hazard rate) have the following form:

\(\mathbf{Y}(t+1) = \mathbf{u} + \mathbf{R} \mathbf{Y}(t) + \mathbf{\epsilon}\)

\(\mu (t, \mathbf{Y}(t)) = [\mu_0 + \mathbf{b} \mathbf{Y}(t) + \mathbf{Y}(t)^* \mathbf{Q} \mathbf{Y}(t)] e^{\theta t}\)

Coefficients \(\mathbf{u}\) (a \(k \times 1\) matrix, where \(k\) is a number of covariates), \(\mathbf{R}\) (a \(k \times k\) matrix), \(\mu_0\), \(\mathbf{b}\) (a \(1 \times k\) matrix), \(\mathbf{Q}\) (a \(k \times k\) matrix) are assumed to be constant in the particular implementation of this model in the R-package stpm. \(\mathbf{\epsilon}\) are normally-distributed random residuals, \(k \times 1\) matrix. A symbol ’*’ denotes transpose operation. \(\theta\) is a parameter to be estimated along with other parameters (\(\mathbf{u}\), \(\mathbf{R}\), \(\mathbf{\mu_0}\), \(\mathbf{b}\), \(\mathbf{Q}\)).

Example

library(stpm)
#Data simulation (200 individuals)
data <- simdata_discr(N=100)
#Estimation of parameters
pars <- spm_discrete(data)
pars
## $dmodel
## $dmodel$theta
## [1] 0.2
## 
## $dmodel$mu0
## [1] 9.136222991e-06
## 
## $dmodel$b
## [1] 2.477520974e-08
## 
## $dmodel$Q
##                 [,1]
## [1,] 1.678638383e-11
## 
## $dmodel$u
## [1] -39.86660584
## 
## $dmodel$u.std.err
##     (Intercept) 
## 1.072800413e-13 
## 
## $dmodel$R
##      [,1]
## [1,] 0.95
## 
## $dmodel$R.std.err
##                 y1_1
## [1,] 1.830749633e-16
## 
## $dmodel$Sigma
## [1] 3.194195589e-12
## 
## 
## $cmodel
## $cmodel$a
##       [,1]
## [1,] -0.05
## 
## $cmodel$f1
##              [,1]
## [1,] -797.3321169
## 
## $cmodel$Q
##                 [,1]
## [1,] 1.678638383e-11
## 
## $cmodel$f
##              [,1]
## [1,] -737.9555354
## 
## $cmodel$b
##                 [,1]
## [1,] 3.194195589e-12
## 
## $cmodel$mu0
##                  [,1]
## [1,] -5.278592871e-09
## 
## $cmodel$theta
## [1] 0.2
## 
## 
## attr(,"class")
## [1] "spm.discrete"

Continuous-time model

In the specification of the SPM described in 2007 paper by Yashin and collegaues (Yashin, Arbeev, Akushevich, et al. 2007) the stochastic differential equation describing the age dynamics of a covariate is:

\(d\mathbf{Y}(t)= \mathbf{a}(t)(\mathbf{Y}(t) -\mathbf{f}_1(t))dt + \mathbf{b}(t)d\mathbf{W}(t), \mathbf{Y}(t=t_0)\)

In this equation, \(\mathbf{Y}(t)\) (a \(k \times 1\) matrix) is the value of a particular covariate at a time (age) \(t\). \(\mathbf{f}_1(t)\) (a \(k \times 1\) matrix) corresponds to the long-term mean value of the stochastic process \(\mathbf{Y}(t)\), which describes a trajectory of individual covariate influenced by different factors represented by a random Wiener process \(\mathbf{W}(t)\). Coefficient \(\mathbf{a}(t)\) (a \(k \times k\) matrix) is a negative feedback coefficient, which characterizes the rate at which the process reverts to its mean. In the area of research on aging, \(\mathbf{f}_1(t)\) represents the mean allostatic trajectory and \(\mathbf{a}(t)\) represents the adaptive capacity of the organism. Coefficient \(\mathbf{b}(t)\) (a \(k \times 1\) matrix) characterizes a strength of the random disturbances from Wiener process \(\mathbf{W}(t)\).

The following function \(\mu(t, \mathbf{Y}(t))\) represents a hazard rate:

\(\mu(t, \mathbf{Y}(t)) = \mu_0(t) + (\mathbf{Y}(t) - \mathbf{f}(t))^* \mathbf{Q}(t) (\mathbf{Y}(t) - \mathbf{f}(t))\)

here \(\mu_0(t)\) is the baseline hazard, which represents a risk when \(\mathbf{Y}(t)\) follows its optimal trajectory; \(\mathbf{f}(t)\) (a \(k \times 1\) matrix) represents the optimal trajectory that minimizes the risk and \(\mathbf{Q}(t)\) (\(k \times k\) matrix) represents a sensitivity of risk function to deviation from the norm.

Example

library(stpm)
#Simulate some data for 50 individuals
data <- simdata_cont(N=50)
head(data)
##   id xi t1 t2            y1       y1.next
## 1  0  0 30 31   71.22667883   28.86212066
## 2  0  0 31 32   28.86212066  -11.43627413
## 3  0  0 32 33  -11.43627413  -49.76927447
## 4  0  0 33 34  -49.76927447  -86.23273469
## 5  0  0 34 35  -86.23273469 -120.91783418
## 6  0  0 35 36 -120.91783418 -153.91130546
#Estimate parameters
# a=-0.05, f1=80, Q=2e-8, f=80, b=5, mu0=2e-5, theta=0.08 are starting values for estimation procedure
pars <- spm_continuous(dat=data,a=-0.05, f1=80, Q=2e-8, f=80, b=5, mu0=2e-5, theta=0.08)
## Parameter a achieved lower/upper bound.
## 0 
## Parameter Q achieved lower/upper bound.
## 0 
## Parameter b achieved lower/upper bound.
## 5.5 
## Parameter mu0 achieved lower/upper bound.
## 1.8e-05
pars
## $a
##      [,1]
## [1,]    0
## 
## $f1
##             [,1]
## [1,] 87.20242295
## 
## $Q
##      [,1]
## [1,]    0
## 
## $f
##             [,1]
## [1,] 74.69160644
## 
## $b
##      [,1]
## [1,]  5.5
## 
## $mu0
## [1] 1.8e-05
## 
## $theta
## [1] 0.08763157956
## 
## $status
## [1] 5
## 
## $LogLik
## [1] -25284.7227
## 
## $objective
## [1] 25284.7227
## 
## $message
## [1] "NLOPT_MAXEVAL_REACHED: Optimization stopped because maxeval (above) was reached."
## 
## $limit
## [1] TRUE
## 
## attr(,"class")
## [1] "spm.continuous"

Coefficient conversion between continuous- and discrete-time models

The coefficient conversion between continuous- and discrete-time models is as follows (‘c’ and ‘d’ denote continuous- and discrete-time models respectively; note: these equations can be used if intervals between consecutive observations of discrete- and continuous-time models are equal; it also required that matrices \(\mathbf{a}_c\) and \(\mathbf{Q}_{c,d}\) must be full-rank matrices):

\(\mathbf{Q}_c = \mathbf{Q}_d\)

\(\mathbf{a}_c = \mathbf{R}_d - I(k)\)

\(\mathbf{b}_c = \mathbf{\Sigma}\)

\({\mathbf{f}_1}_c = -\mathbf{a}_c^{-1} \times \mathbf{u}_d\)

\(\mathbf{f}_c = -0.5 \mathbf{b}_d \times \mathbf{Q}^{-1}_d\)

\({\mu_0}_c = {\mu _0}_d - \mathbf{f}_c \times \mathbf{Q_c} \times \mathbf{f}_c^*\)

\(\theta_c = \theta_d\)

where \(k\) is a number of covariates, which is equal to model’s dimension and ’*’ denotes transpose operation; \(\mathbf{\Sigma}\) is a \(k \times 1\) matrix which contains s.d.s of corresponding residuals (residuals of a linear regression \(\mathbf{Y}(t+1) = \mathbf{u} + \mathbf{R}\mathbf{Y}(t) + \mathbf{\epsilon}\); s.d. is a standard deviation), \(I(k)\) is an identity \(k \times k\) matrix.

Model with time-dependent coefficients

In previous models, we assumed that coefficients is sort of time-dependant: we multiplied them on to \(e^{\theta t}\). In general, this may not be the case (Yashin, Arbeev, Kulminski, et al. 2007). We extend this to a general case, i.e. (we consider one-dimensional case):

\(\mathbf{a(t)} = \mathbf{par}_1 t + \mathbf{par}_2\) - linear function.

The corresponding equations will be equivalent to one-dimensional continuous case described above.

Example

library(stpm)
#Data preparation:
n <- 10
data <- simdata_time_dep(N=n)
# Estimation:
opt.par <- spm_time_dep(data, 
                        start = list(a = -0.05, f1 = 80, Q = 2e-08, f = 80, b = 5, mu0 = 0.001), 
                        frm = list(at = "a", f1t = "f1", Qt = "Q", ft = "f", bt = "b", mu0t= "mu0"))
opt.par
## $a
## [1] -0.05
## 
## $f1
## [1] 80
## 
## $Q
## [1] 2e-08
## 
## $f
## [1] 80
## 
## $b
## [1] 5
## 
## $mu0
## [1] 0.001

Setting lower and upper boundaries of the model parameters

Lower and upper boundaries can be set up with parameters \(lb\) and \(ub\), which represents simple numeric vectors. Note: lengths of \(lb\) and \(ub\) must be the same as the total length of the parameters. Lower and upper boundaries can be set for continuous-time and time-dependent models only.

Setting lb and ub for continuous-time model

One covariate

Below we show the example of setting up \(lb\) and \(ub\) when we have a single covariate:

library(stpm)
data <- simdata_cont(N=10, ystart = 80, a = -0.1, Q = 1e-06, mu0 = 1e-5, theta = 0.08, f1 = 80, f=80, b=1, dt=1, sd0=5)
ans <- spm_continuous(dat=data,
                      a = -0.1,
                      f1 = 82, 
                      Q = 1.4e-6,
                      f = 77,
                      b = 1,
                      mu0 = 1.6e-5,
                      theta = 0.1,
                      stopifbound = FALSE,
                      lb=c(-0.2, 60, 0.1e-6, 60, 0.1, 0.1e-5, 0.01), 
                      ub=c(0, 140, 5e-06, 140, 3, 5e-5, 0.20))
ans
## $a
##               [,1]
## [1,] -0.1021716235
## 
## $f1
##             [,1]
## [1,] 80.44006779
## 
## $Q
##                 [,1]
## [1,] 3.858234307e-07
## 
## $f
##             [,1]
## [1,] 72.10982571
## 
## $b
##              [,1]
## [1,] 0.9480323775
## 
## $mu0
## [1] 2.149143207e-05
## 
## $theta
## [1] 0.146372252
## 
## $status
## [1] 5
## 
## $LogLik
## [1] -660.5948954
## 
## $objective
## [1] 660.5604444
## 
## $message
## [1] "NLOPT_MAXEVAL_REACHED: Optimization stopped because maxeval (above) was reached."
## 
## $limit
## [1] FALSE
## 
## attr(,"class")
## [1] "spm.continuous"

Two covariates

This is an example for two physiological variables (covariates).

library(stpm)

data <- simdata_cont(N=10, 
                     a=matrix(c(-0.1,  0.001, 0.001, -0.1), nrow = 2, ncol = 2, byrow = T),
                     f1=t(matrix(c(100, 200), nrow = 2, ncol = 1, byrow = F)),
                     Q=matrix(c(1e-06, 1e-7, 1e-7,  1e-06), nrow = 2, ncol = 2, byrow = T),
                     f=t(matrix(c(100, 200), nrow = 2, ncol = 1, byrow = F)),
                     b=matrix(c(1, 2), nrow = 2, ncol = 1, byrow = F),
                     mu0=1e-4,
                     theta=0.08,
                     ystart = c(100,200), sd0=c(5, 10), dt=1)

a.d <- matrix(c(-0.15,  0.002, 0.002, -0.15), nrow = 2, ncol = 2, byrow = T)
f1.d <- t(matrix(c(95, 195), nrow = 2, ncol = 1, byrow = F))
Q.d <- matrix(c(1.2e-06, 1.2e-7, 1.2e-7,  1.2e-06), nrow = 2, ncol = 2, byrow = T)
f.d <- t(matrix(c(105, 205), nrow = 2, ncol = 1, byrow = F))
b.d <- matrix(c(1, 2), nrow = 2, ncol = 1, byrow = F)
mu0.d <- 1.1e-4
theta.d <- 0.07

ans <- spm_continuous(dat=data,
                      a = a.d, 
                      f1 = f1.d,
                      Q = Q.d,
                      f = f.d,
                      b = b.d,
                      mu0 = mu0.d,
                      theta = theta.d,
                      lb=c(-0.5, ifelse(a.d[2,1] > 0, a.d[2,1]-0.5*a.d[2,1], a.d[2,1]+0.5*a.d[2,1]), ifelse(a.d[1,2] > 0, a.d[1,2]-0.5*a.d[1,2], a.d[1,2]+0.5*a.d[1,2]), -0.5,  
                           80, 100, 
                           Q.d[1,1]-0.5*Q.d[1,1], ifelse(Q.d[2,1] > 0, Q.d[2,1]-0.5*Q.d[2,1], Q.d[2,1]+0.5*Q.d[2,1]), ifelse(Q.d[1,2] > 0, Q.d[1,2]-0.5*Q.d[1,2], Q.d[1,2]+0.5*Q.d[1,2]), Q.d[2,2]-0.5*Q.d[2,2],
                           80, 100,
                           0.1, 0.5,
                           0.1e-4,
                           0.01),
                      ub=c(-0.08,  0.002,  0.002, -0.08,  
                           110, 220, 
                           Q.d[1,1]+0.1*Q.d[1,1], ifelse(Q.d[2,1] > 0, Q.d[2,1]+0.1*Q.d[2,1], Q.d[2,1]-0.1*Q.d[2,1]), ifelse(Q.d[1,2] > 0, Q.d[1,2]+0.1*Q.d[1,2], Q.d[1,2]-0.1*Q.d[1,2]), Q.d[2,2]+0.1*Q.d[2,2],
                           110, 220,
                           1.5, 2.5,
                           1.2e-4,
                           0.10))
ans
## $a
##                [,1]            [,2]
## [1,] -0.15093725674  0.001664644568
## [2,]  0.00197532678 -0.148792208813
## 
## $f1
##             [,1]
## [1,] 106.8774537
## [2,] 194.8962017
## 
## $Q
##                 [,1]            [,2]
## [1,] 1.226075476e-06 1.226281156e-07
## [2,] 1.226793880e-07 1.227885451e-06
## 
## $f
##             [,1]
## [1,] 105.0787389
## [2,] 204.6303638
## 
## $b
##             [,1]
## [1,] 1.140562931
## [2,] 1.912680217
## 
## $mu0
## [1] 0.0001111659552
## 
## $theta
## [1] 0.07457707686
## 
## $status
## [1] 5
## 
## $LogLik
## [1] 1348.805053
## 
## $objective
## [1] -1951.919043
## 
## $message
## [1] "NLOPT_MAXEVAL_REACHED: Optimization stopped because maxeval (above) was reached."
## 
## $limit
## [1] FALSE
## 
## attr(,"class")
## [1] "spm.continuous"

Setting lb and ub for model with time-dependent coefficients

This model uses only one covariate, therefore setting-up model parameters is easy:

n <- 10
data <- simdata_time_dep(N=n)
# Estimation:
opt.par <- spm_time_dep(data, start=list(a=-0.05, f1=80, Q=2e-08, f=80, b=5, mu0=0.001), 
                        lb=c(-1, 30, 1e-8, 30, 1, 1e-6), ub=c(0, 120, 5e-8, 130, 10, 1e-2))
opt.par
## $a
## [1] -0.05
## 
## $f1
## [1] 80
## 
## $Q
## [1] 2e-08
## 
## $f
## [1] 80
## 
## $b
## [1] 5
## 
## $mu0
## [1] 0.001

Special case when some model parameter functions are equal to zero

Imagine a situation when one parameter function you want to be equal to zero: \(f=0\). Let’s emulate this case:

library(stpm)
n <- 10
data <- simdata_time_dep(N=n)
# Estimation:
opt.par <- spm_time_dep(data, frm = list(at="a", f1t="f1", Qt="Q", ft="0", bt="b", mu0t="mu0"))
opt.par
## $a
## [1] -0.05
## 
## $f1
## [1] 80
## 
## $Q
## [1] 2e-08
## 
## $b
## [1] 80
## 
## $mu0
## [1] 5
## 
## $<NA>
## <NA> 
##   NA

As you can see, there is no parameter \(f\) in \(opt.par\). This is because we set \(f=0\) in \(frm\)!

Then, is you want to set the constraints, you must not specify the starting value (parameter \(start\)) and \(lb\)/\(ub\) for the parameter \(f\) (otherwise, the function raises an error):

n <- 10
data <- simdata_time_dep(N=n)
# Temporarily commented below
# Estimation:
opt.par <- spm_time_dep(data, frm = list(at="a", f1t="f1", Qt="Q", ft="0", bt="b", mu0t="mu0"), 
                        start=list(a=-0.05, f1=80, Q=2e-08, b=5, mu0=0.001), 
                        lb=c(-1, 30, 1e-8, 1, 1e-6), ub=c(0, 120, 5e-8, 10, 1e-2))
opt.par
## $a
## [1] -0.05
## 
## $f1
## [1] 80
## 
## $Q
## [1] 2e-08
## 
## $b
## [1] 5
## 
## $mu0
## [1] 0.001

You can do the same manner if you want two or more parameters to be equal to zero.

Fast parameter estimating for one-dimentional model using spm_con_1d

Function spm_con_1d(...) allows for very fast parameter estimating for one-dimensional model. This function implements a analytical solution to estimate the parameters in the continuous SPM model by assuming all the parameters are constants. Below there is an example.

library(stpm) 
dat <- simdata_cont(N=500)
colnames(dat) <- c("id", "xi", "t1", "t2", "y", "y.next")
res <- spm_con_1d(as.data.frame(dat), a=-0.05, b=2, q=1e-8, f=80, f1=90, mu0=1e-3, theta=0.08)
## [1] "Initial values:"
## [1] -5e-02  2e+00  1e-08  8e+01  9e+01  1e-03  8e-02
## [1] "Lower bounds:"
## [1] -2.5e-01  2.0e-01  1.0e-10  4.0e+01  4.5e+01  1.0e-04  8.0e-02
## [1] "Upper bounds:"
## [1] -2.5e-03  1.0e+01  1.0e-07  1.6e+02  1.8e+02  1.0e-02  8.0e-02
res
## $est
##               Coeff. Std. Err.   z p value
## a     -0.05042051891       NaN NaN     NaN
## b      4.95368016985       NaN NaN     NaN
## q      0.00000000010       NaN NaN     NaN
## f     79.99999993052       NaN NaN     NaN
## f1    79.67871456091       NaN NaN     NaN
## mu0    0.00010000000       NaN NaN     NaN
## theta  0.08000000000       NaN NaN     NaN
## 
## $hessian
##                      a                b   q                f               f1
## a      2.675297952e+05  6.283305907e+03 NaN -1.562887391e-06  3.194937715e+00
## b      6.283305907e+03  1.734549016e+03 NaN -2.128573103e-11  5.229535555e-09
## q                  NaN              NaN NaN              NaN              NaN
## f     -1.562887391e-06 -2.128573103e-11 NaN  5.783794625e-06  4.346227451e-08
## f1     3.194937715e+00  5.229535555e-09 NaN  4.346227451e-08  3.309451871e+00
## mu0    1.123993491e+00  2.189275792e-02 NaN  5.435604246e-03 -4.434279628e-04
## theta  1.585146614e-02  3.100999887e-04 NaN  4.904923279e-05 -4.183547439e-06
##                    mu0            theta
## a      1.123993491e+00  1.585146614e-02
## b      2.189275792e-02  3.100999887e-04
## q                  NaN              NaN
## f      5.435604246e-03  4.904923279e-05
## f1    -4.434279628e-04 -4.183547439e-06
## mu0   -2.131313146e+08  1.848322797e+09
## theta  1.848322797e+09  1.704572385e+07
## 
## $lik
## [1] 50471.63876
## 
## $con
## [1] 4
## 
## $message
## [1] "NLOPT_XTOL_REACHED: Optimization stopped because xtol_rel or xtol_abs (above) was reached."

Simulation (individual trajectory projection, also known as microsimulations)

We added one- and multi- dimensional simulation to be able to generate test data for hyphotesis testing. Data, which can be simulated can be discrete (equal intervals between observations) and continuous (with arbitrary intervals).

Discrete-time simulation

The corresponding function is (k - a number of variables(covariates), equal to model’s dimension):

simdata_discr(N=100, a=-0.05, f1=80, Q=2e-8, f=80, b=5, mu0=1e-5, theta=0.08, ystart=80, tstart=30, tend=105, dt=1)

Here:

N - Number of individuals

a - A matrix of kxk, which characterize the rate of the adaptive response

f1 - A particular state, which if a deviation from the normal (or optimal). This is a vector with length of k

Q - A matrix of k by k, which is a non-negative-definite symmetric matrix

f - A vector-function (with length k) of the normal (or optimal) state

b - A diffusion coefficient, k by k matrix

mu0 - mortality at start period of time (baseline hazard)

theta - A displacement coefficient of the Gompertz function

ystart - A vector with length equal to number of dimensions used, defines starting values of covariates

tstart - A number that defines a start time (30 by default). Can be a number (30 by default) or a vector of two numbers: c(a, b) - in this case, starting value of time is simulated via uniform(a,b) distribution.

tend - A number, defines a final time (105 by default)

dt - A time interval between observations.

This function returns a table with simulated data, as shown in example below:

library(stpm)
data <- simdata_discr(N=10)
head(data)
##   id xi t1 t2          y1      y1.next
## 1  1  0 30 31 78.37568207  80.14572528
## 2  1  0 31 32 80.14572528  90.89021177
## 3  1  0 32 33 90.89021177  85.35114798
## 4  1  0 33 34 85.35114798  89.35097104
## 5  1  0 34 35 89.35097104  98.87706896
## 6  1  0 35 36 98.87706896 100.91052534

Continuous-time simulation

The corresponding function is (k - a number of variables(covariates), equal to model’s dimension):

simdata_cont(N=100, a=-0.05, f1=80, Q=2e-07, f=80, b=5, mu0=2e-05, theta=0.08, ystart=80, tstart=c(30,50), tend=105)

Here:

N - Number of individuals

a - A matrix of kxk, which characterize the rate of the adaptive response

f1 - A particular state, which if a deviation from the normal (or optimal). This is a vector with length of k

Q - A matrix of k by k, which is a non-negative-definite symmetric matrix

f - A vector-function (with length k) of the normal (or optimal) state

b - A diffusion coefficient, k by k matrix

mu0 - mortality at start period of time (baseline hazard)

theta - A displacement coefficient of the Gompertz function

ystart - A vector with length equal to number of dimensions used, defines starting values of covariates

tstart - A number that defines a start time (30 by default). Can be a number (30 by default) or a vector of two numbers: c(a, b) - in this case, starting value of time is simulated via uniform(a,b) distribution.

tend - A number, defines a final time (105 by default)

This function returns a table with simulated data, as shown in example below:

library(stpm)
data <- simdata_cont(N=10)
head(data)
##   id xi          t1          t2          y1     y1.next
## 1  0  0 37.32465961 39.21089196 78.83923290 78.60497011
## 2  0  0 39.21089196 40.73386605 78.60497011 78.82752695
## 3  0  0 40.73386605 41.78966292 78.82752695 78.45734388
## 4  0  0 41.78966292 43.49786515 78.45734388 76.33148060
## 5  0  0 43.49786515 45.01146966 76.33148060 72.53230450
## 6  0  0 45.01146966 46.32764222 72.53230450 63.88171410

SPM with partially observed covariates

Stochastic Process Model has many applications in analysis of longitudinal biodemographic data. Such data contain various physiological variables (known as covariates). Data can also potentially contain genetic information available for all or a part of participants. Taking advantage from both genetic and non-genetic information can provide future insights into a broad range of processes describing aging-related changes in the organism.

Method

In this package, SPM with partially observed covariates is implemented in form of GenSPM (Genetic SPM), presented in (Arbeev et al. 2009) and further advanced in (Arbeev et al. 2014), further elaborates the basic stochastic process model conception by introducing a categorical variable, \(Z\), which may be a specific value of a genetic marker or, in general, any categorical variable. Currently, \(Z\) has two gradations: 0 or 1 in a genetic group of interest, assuming that \(P(Z=1) = p\), \(p \in [0, 1]\), were \(p\) is the proportion of carriers and non-carriers of an allele in a population. Example of longitudinal data with genetic component \(Z\) is provided below.

library(stpm)
data <- sim_pobs(N=10)
head(data)
##   id xi          t1          t2 Z          y1     y1.next
## 1  0  0 59.31488727 60.34278590 0 78.48691974 85.37612372
## 2  0  0 60.34278590 61.29360276 0 85.37612372 74.07415518
## 3  0  0 61.29360276 62.31176143 0 74.07415518 82.04855154
## 4  0  0 62.31176143 63.30138311 0 82.04855154 71.59308147
## 5  0  0 63.30138311 64.38932967 0 71.59308147 73.10168234
## 6  0  0 64.38932967 65.39889616 0 73.10168234 66.49071358

In the specification of the SPM described in 2007 paper by Yashin and colleagues (Yashin, Arbeev, Akushevich, et al. 2007) the stochastic differential equation describing the age dynamics of a physiological variable (a dynamic component of the model) is:

\(dY(t) = a(Z, t)(Y(t) - f1(Z, t))dt + b(Z, t)dW(t), Y(t = t_0)\)

Here in this equation, \(Y(t)\) is a \(k \times 1\) matrix, where \(k\) is a number of covariates, which is a model dimension) describing the value of a physiological variable at a time (e.g. age) t. \(f_1(Z,t)\) is a \(k \times 1\) matrix that corresponds to the long-term average value of the stochastic process \(Y(t)\), which describes a trajectory of individual variable influenced by different factors represented by a random Wiener process \(W(t)\). The negative feedback coefficient \(a(Z,t)\) (\(k \times k\) matrix) characterizes the rate at which the stochastic process goes to its mean. In research on aging and well-being, \(f_1(Z,t)\) represents the average allostatic trajectory and \(a(t)\) in this case represents the adaptive capacity of the organism. Coefficient \(b(Z,t)\) (\(k \times 1\) matrix) characterizes a strength of the random disturbances from Wiener process \(W(t)\). All of these parameters depend on \(Z\) (a genetic marker having values 1 or 0). The following function \(\mu(t,Y(t))\) represents a hazard rate:

\(\mu(t,Y(t)) = \mu_0(t) + (Y(t) - f(Z, t))^*Q(Z, t)(Y(t) - f(Z, t))\)

In this equation: \(\mu_0(t)\) is the baseline hazard, which represents a risk when \(Y(t)\) follows its optimal trajectory; f(t) (\(k \times 1\) matrix) represents the optimal trajectory that minimizes the risk and \(Q(Z, t)\) (\(k \times k\) matrix) represents a sensitivity of risk function to deviation from the norm. In general, model coefficients \(a(Z, t)\), \(f1(Z, t)\), \(Q(Z, t)\), \(f(Z, t)\), \(b(Z, t)\) and \(\mu_0(t)\) are time(age)-dependent. Once we have data, we then can run analysis, i.e. estimate coefficients (they are assumed to be time-independent and data here is simulated):

library(stpm)
#Generating data:
data <- sim_pobs(N=10)
head(data)
##   id xi          t1          t2 Z          y1     y1.next
## 1  0  0 91.15401210 92.21965967 0 81.36255177 79.33769707
## 2  0  0 92.21965967 93.16039286 0 79.33769707 84.88844550
## 3  0  0 93.16039286 94.11000198 0 84.88844550 84.66125367
## 4  0  0 94.11000198 95.17582943 0 84.66125367 88.37073761
## 5  0  0 95.17582943 96.19965652 0 88.37073761 92.67507043
## 6  0  0 96.19965652 97.25163375 0 92.67507043 97.88009840
#Parameters estimation:
pars <- spm_pobs(x=data)
## Parameter mu0L achieved lower/upper bound.
## 9e-06
pars
## $aH
##                [,1]
## [1,] -0.04795259246
## 
## $aL
##                [,1]
## [1,] -0.01099943278
## 
## $f1H
##             [,1]
## [1,] 61.23900359
## 
## $f1L
##             [,1]
## [1,] 72.00524129
## 
## $QH
##                 [,1]
## [1,] 1.680364984e-08
## 
## $QL
##                 [,1]
## [1,] 2.240941098e-08
## 
## $fH
##             [,1]
## [1,] 64.05396857
## 
## $fL
##             [,1]
## [1,] 76.44837581
## 
## $bH
##             [,1]
## [1,] 4.024384937
## 
## $bL
##             [,1]
## [1,] 4.892744211
## 
## $mu0H
## [1] 7.282797918e-06
## 
## $mu0L
## [1] 9e-06
## 
## $thetaH
## [1] 0.07235009977
## 
## $thetaL
## [1] 0.09
## 
## $p
## [1] 0.2295535296
## 
## $limit
## [1] TRUE
## 
## attr(,"class")
## [1] "pobs.spm"

Here and represents parameters when \(Z\) = 1 (H) and 0 (L).

###Joint analysis of two datasets: first dataset with genetic and second dataset with non-genetic component

library(stpm)
data.genetic <- sim_pobs(N=5, mode='observed')
head(data.genetic)
##   id xi          t1          t2 Z          y1     y1.next
## 1  0  0 52.84443870 53.85097898 0 80.89977188 83.19630112
## 2  0  0 53.85097898 54.78242384 0 83.19630112 80.25782131
## 3  0  0 54.78242384 55.72679338 0 80.25782131 86.44334846
## 4  0  0 55.72679338 56.79405094 0 86.44334846 97.88669379
## 5  0  0 56.79405094 57.69972449 0 97.88669379 94.65498786
## 6  0  0 57.69972449 58.72380119 0 94.65498786 99.17243895
data.nongenetic <- sim_pobs(N=10, mode='unobserved')
head(data.nongenetic)
##   id xi          t1          t2          y1     y1.next
## 1  0  0 58.97322021 59.92391381 81.25664113 79.44120006
## 2  0  0 59.92391381 60.86036779 79.44120006 69.59782970
## 3  0  0 60.86036779 61.94182249 69.59782970 61.99190512
## 4  0  0 61.94182249 62.97351575 61.99190512 60.87646177
## 5  0  0 62.97351575 63.92779874 60.87646177 64.14649756
## 6  0  0 63.92779874 64.95188836 64.14649756 53.97920321
#Parameters estimation:
pars <- spm_pobs(x=data.genetic, y = data.nongenetic, mode='combined')
## Parameter fH achieved lower/upper bound.
## 66
pars
## $aH
##                [,1]
## [1,] -0.05038474709
## 
## $aL
##                [,1]
## [1,] -0.01023813814
## 
## $f1H
##             [,1]
## [1,] 56.48646989
## 
## $f1L
##             [,1]
## [1,] 87.91714119
## 
## $QH
##                 [,1]
## [1,] 1.651970204e-08
## 
## $QL
##                 [,1]
## [1,] 2.606859538e-08
## 
## $fH
##      [,1]
## [1,]   66
## 
## $fL
##             [,1]
## [1,] 83.61105956
## 
## $bH
##            [,1]
## [1,] 4.38852204
## 
## $bL
##             [,1]
## [1,] 5.133226046
## 
## $mu0H
## [1] 7.611760398e-06
## 
## $mu0L
## [1] 9.008979879e-06
## 
## $thetaH
## [1] 0.07288271682
## 
## $thetaL
## [1] 0.09006074495
## 
## $p
## [1] 0.2455866161
## 
## $limit
## [1] TRUE
## 
## attr(,"class")
## [1] "pobs.spm"

Here mode ‘observed’ is used for simlation of data with genetic component \(Z\) and ‘unobserved’ - without genetic component.

Genetic SPM ‘GSPM’

This type of SPM also uses genetic component by analogy from the previous chapters but uses explicit gradient function which speeds up computations significantly. See (He et al. 2017) for details. Below we provide examples of usage:

library(stpm) 
data(ex_spmcon1dg)
head(ex_data$spm_data)
##   id xi t1 t2           y      y.next
## 1  1  0 30 31 2.000000000 2.024328135
## 2  1  0 31 32 2.024328135 1.927486318
## 3  1  0 32 33 1.927486318 1.899083801
## 4  1  0 33 34 1.899083801 2.061574385
## 5  1  0 34 35 2.061574385 2.034558435
## 6  1  0 35 36 2.034558435 2.114382051
head(ex_data$gene_data)
##   id geno
## 1  1    1
## 2  2    1
## 3  3    0
## 4  4    0
## 5  5    1
## 6  6    0
res <- spm_con_1d_g(spm_data=ex_data$spm_data, 
                    gene_data=ex_data$gene_data, 
                    a = -0.02, b=0.2, q=0.01, f=3, f1=3, mu0=0.01, theta=1e-05, 
                    upper=c(-0.01,3,0.1,10,10,0.1,1e-05), lower=c(-1,0.01,0.00001,1,1,0.001,1e-07), 
                    effect=c('q'), method = "tnewton")
## [1] "Initial values:"
##  [1] -2e-02 -2e-02  2e-01  2e-01  1e-02  1e-02  3e+00  3e+00  3e+00  3e+00
## [11]  1e-02  1e-02  1e-05
## [1] "Lower bounds:"
##  [1] -1e+00 -1e+00  1e-02  1e-02  1e-05  1e-05  1e+00  1e+00  1e+00  1e+00
## [11]  1e-03  1e-03  1e-07
## [1] "Upper bounds:"
##  [1] -1e-02 -1e-02  3e+00  3e+00  1e-01  1e-01  1e+01  1e+01  1e+01  1e+01
## [11]  1e-01  1e-01  1e-05
res
## $est
##                Coeff.       Std. Err.                z      p value
## a     -0.031238526396 9.063948538e-04 -3.446458932e+01 0.0000000000
## b      0.101329776535 2.785757812e-04  3.637422323e+02 0.0000000000
## q_0    0.000010000000 3.124724202e-03  3.200282442e-03 0.9974465484
## q_2    0.003608789876 5.278940497e-03  6.836201087e-01 0.4942150834
## f      2.999731527110 3.052250225e-02  9.827934494e+01 0.0000000000
## f1     3.000056947322 1.716252597e-02  1.748027622e+02 0.0000000000
## mu0    0.001000000000 1.552500478e+04  6.441221851e-08 0.9999999486
## theta  0.000000100000 8.987112435e+01  1.112704450e-09 0.9999999991
## 
## $lik
## [1] -121546.402
## 
## $con
## [1] -1
## 
## $message
## [1] "NLOPT_FAILURE: Generic failure code."
## 
## $hessian
##                      a                b             q_0              q_2
## a      2.354652039e+06  6.786921795e+05   158.058062204    82.3853495314
## b      6.786921795e+05  1.325446527e+07  1517.123161715   775.1749357615
## q_0    1.580580622e+02  1.517123162e+03   -39.193641214   -13.8346645758
## q_2    8.238534953e+01  7.751749358e+02   -13.834664576   -13.8346698668
## f      3.489487410e-01 -4.196232148e-03 36528.531980456 19000.8175924782
## f1    -8.446522178e+04  5.065699701e+01    -2.840122726    -1.4581382745
## mu0    2.228043741e-03  8.307251846e-04     1.053761684     0.6444883984
## theta  4.230830818e-01  1.744572073e-01   194.433867000   117.7827687934
##                      f               f1              mu0            theta
## a      3.489487410e-01 -8.446522178e+04  2.228043741e-03   0.423083081841
## b     -4.196232148e-03  5.065699701e+01  8.307251846e-04   0.174457207322
## q_0    3.652853198e+04 -2.840122726e+00  1.053761684e+00 194.433867000043
## q_2    1.900081759e+04 -1.458138275e+00  6.444883984e-01 117.782768793404
## f      1.670661873e+02 -1.202720644e-02  3.587745596e-03   0.465787947178
## f1    -1.202720644e-02  6.470562259e+03 -5.237961886e-05  -0.006925314665
## mu0    3.587745596e-03 -5.237961886e-05  0.000000000e+00   0.000000000000
## theta  4.657879472e-01 -6.925314665e-03  0.000000000e+00   0.000000000000
## 
## $beta
##                  Coeff.      Std. Err.         Chi. Sq      p value
## beta_a               NA             NA              NA           NA
## beta_b               NA             NA              NA           NA
## beta_q   0.001799394938 0.004201832138 0.0007705738917 0.9778541969
## beta_f               NA             NA              NA           NA
## beta_f1              NA             NA              NA           NA
## beta_mu0             NA             NA              NA           NA

Here: spm_data - A dataset for the SPM model. See the STPM package for more details about the format.

gene_data - A two column dataset containing the genotypes for the individuals in spm_data. The first column id is the ID of the individuals in dataset spm_data, and the second column geno is the genotype.

a - The initial value for the paramter . The initial value will be predicted if not specified.

b - The initial value for the paramter . The initial value will be predicted if not specified.

q - The initial value for the paramter . The initial value will be predicted if not specified.

f - The initial value for the paramter . The initial value will be predicted if not specified.

f1 - The initial value for the paramter . The initial value will be predicted if not specified.

mu0 - The initial value for the paramter in the baseline hazard. The initial value will be predicted if not specified.

theta - The initial value for the paramter in the baseline hazard. The initial value will be predicted if not specified.

lower - A vector of the lower bound of the parameters.

upper - A vector of the upper bound of the parameters.

effect - A character vector of the parameters that are linked to genotypes. The vector can contain any combination of , , , , .

control - A list of the control parameters for the optimization paramters.

global - A logical variable indicating whether the MLSL (TRUE) or the L-BFGS (FALSE) algorithm is used for the optimization.

verbose - A logical variable indicating whether initial information is printed.

ahessian - A logical variable indicating whether the approximate (FALSE) or analytical (TRUE) Hessian is returned.

est - The estimates of the parameters.

hessian - The Hessian matrix of the estimates.

lik - The minus log-likelihood.

con - A number indicating the convergence. See the ‘nloptr’ package for more details.

message - Extra message about the convergence. See the ‘nloptr’ package for more details.

beta - The coefficients of the genetic effect on the parameters to be linked to genotypes.

Multiple imputation with spm.impute(…)

The SPM offers longitudinal data imputation with results that are better than from other imputation tools since it preserves data structure, i.e. relation between Y(t) and mu(Y(t),t). Below there are two examples of multiple data imputation with function spm.impute(…).

library(stpm)

#######################################################################
############## One dimensional case (one covariate) ###################
#######################################################################

## Data preparation (short format)#
data <- simdata_discr(N=1000, dt = 2, format="short")

miss.id <- sample(x=dim(data)[1], size=round(dim(data)[1]/4)) # ~25% missing data
incomplete.data <- data
incomplete.data[miss.id,4] <- NA
# End of data preparation #

##### Multiple imputation with SPM #####
imp.data <- spm.impute(x=incomplete.data, id=1, case="xi", t1=3, covariates="y1", minp=1, theta_range=seq(0.075, 0.09, by=0.001))$imputed

##### Look at the incomplete data with missings #####
head(incomplete.data)
##   id xi  t          y1
## 1  1  0 30 83.17665677
## 2  1  0 32 85.35850084
## 3  1  0 34 78.87696269
## 4  1  0 36 77.25905345
## 5  1  0 38 72.27155840
## 6  1  0 40          NA
##### Look at the imputed data #####
head(imp.data)
##   id xi  t          y1
## 1  1  0 30 83.17665677
## 2  1  0 32 85.35850084
## 3  1  0 34 78.87696269
## 4  1  0 36 77.25905345
## 5  1  0 38 72.27155840
## 6  1  0 40 72.63964208
#########################################################
################ Two-dimensional case ###################
#########################################################

## Parameters for data simulation #
a <- matrix(c(-0.05, 0.01, 0.01, -0.05), nrow=2)
f1 <- matrix(c(90, 30), nrow=1, byrow=FALSE)
Q <- matrix(c(1e-7, 1e-8, 1e-8, 1e-7), nrow=2)
f0 <- matrix(c(80, 25), nrow=1, byrow=FALSE)
b <- matrix(c(5, 3), nrow=2, byrow=TRUE)
mu0 <- 1e-04
theta <- 0.07
ystart <- matrix(c(80, 25), nrow=2, byrow=TRUE)

## Data preparation #
data <- simdata_discr(N=1000, a=a, f1=f1, Q=Q, f=f0, b=b, ystart=ystart, mu0 = mu0, theta=theta, dt=2, format="short")

## Delete some observations in order to have approx. 25% missing data
incomplete.data <- data
miss.id <- sample(x=dim(data)[1], size=round(dim(data)[1]/4)) 
incomplete.data <- data
incomplete.data[miss.id,4] <- NA
miss.id <- sample(x=dim(data)[1], size=round(dim(data)[1]/4)) 
incomplete.data[miss.id,5] <- NA
## End of data preparation #

###### Multiple imputation with SPM #####
imp.data <- spm.impute(x=incomplete.data, id=1, case="xi", t1=3, covariates=c("y1", "y2"), minp=1, theta_range=seq(0.060, 0.07, by=0.001))$imputed

###### Look at the incomplete data with missings #####
head(incomplete.data)
##   id xi  t          y1          y2
## 1  1  0 30 72.11169404 27.33882354
## 2  1  0 32 67.63158526 27.77865814
## 3  1  0 34 63.09629646 29.49861518
## 4  1  0 36 66.90660336 33.76277091
## 5  1  0 38 67.11457685 33.68933156
## 6  1  0 40 71.43534015 34.73774375
###### Look at the imputed data #####
head(imp.data)
##   id xi  t          y1          y2
## 1  1  0 30 72.11169404 27.33882354
## 2  1  0 32 67.63158526 27.77865814
## 3  1  0 34 63.09629646 29.49861518
## 4  1  0 36 66.90660336 33.76277091
## 5  1  0 38 67.11457685 33.68933156
## 6  1  0 40 71.43534015 34.73774375

Prediction

We provide a simple function to predict the next value of . Refer to the example below:

#library(stpm)
#data <- simdata_discr(N=100, format="long")
#res <- spm_discrete(data)
#splitted <- split(data, data$id)
#df <- data.frame()
#lapply(1:100, function(i) {df<<-rbind(df,splitted[[i]][dim(splitted[[i]])[1],c("id", "xi", "t1", "y1")])})
#names(df) <- c("id", "xi", "t", "y")
#predicted <- predict(object=res, data=df, dt=3)
#head(predicted)

Hypothesis testing

The package offers following five hypotheses to test for function (Arbeev et al. 2016):

H01: \(Q(t)=0\) (i.e., \(a_Q = 0\) and \(b_Q = 0\),so that there is no quadratic term in the hazard rate and mortality is described by the baseline Gompertz rate \(μ_0(t)\)).

H02: \(Q(t) = a_Q\) (i.e., \(b_Q = 0\)).

H03: \(f_1(t) = 0\) (i.e., \(a_{f1} = 0\) and \(b_{f1} = 0\)).

H04: \(f_1(t) = a_{f_1}\) (i.e., \(b_{f_1} = 0\)).

H05: \(a(t) = a_Y\) (i.e., \(b_Y = 0\)).

To perform hypothesis testing you should put the variable lrtest to TRUE (this is "H01" by default) or to any of the following: "H01", "H02", "H03", "H04", "H05".

Examples of hypothesis testing

library(stpm)
n <- 1000

# Data simulation:
data <- simdata_time_dep(N=n, format="long")
head(data)

# Hypotheses testing

## H01
res <- spm_time_dep(data, verbose=F, 
                    frm = list(at="a", f1t="f1", Qt="Q", ft="f", bt="b", mu0t="mu0"),
                    start=list(a=-0.05, f1=80, Q=1e-8, f=90, b=5, mu0=0.001), 
                    lb=c(a=-1, f1=30, Q=1e-9, f=10, b=1, mu0=1e-6), 
                    ub=c(a=0, f1=120, Q=1e-7, f=150, b=10, mu0=1e-2),
                    opts = list(algorithm = "NLOPT_LN_NELDERMEAD", 
                    maxeval = 200, ftol_rel = 1e-12), lrtest="H01")

res$alternative$lr.test.pval

## H02
res <- spm_time_dep(data, verbose=F, 
                    frm = list(at="a", f1t="f1", Qt="1e-6", ft="f", bt="b", mu0t="mu0"),
                    start=list(a=-0.05, f1=80, Q=1e-8, f=90, b=5, mu0=0.001), 
                    lb=c(a=-1, f1=30, Q=1e-9, f=10, b=1, mu0=1e-6), 
                    ub=c(a=0, f1=120, Q=1e-7, f=150, b=10, mu0=1e-2),
                    opts = list(algorithm = "NLOPT_LN_NELDERMEAD", 
                    maxeval = 200, ftol_rel = 1e-12), lrtest="H02")

res$alternative$lr.test.pval

## H03
res <- spm_time_dep(data, verbose=F, 
                    frm = list(at="a", f1t="f1", Qt="Q", ft="f", bt="b", mu0t="mu0"),
                    start=list(a=-0.05, f1=80, Q=1e-8, f=90, b=5, mu0=0.001), 
                    ub=c(a=0, f1=120, Q=1e-7, f=150, b=10, mu0=1e-2),
                    opts = list(algorithm = "NLOPT_LN_NELDERMEAD", 
                    maxeval = 200, ftol_rel = 1e-12), lrtest="H03")

res$alternative$lr.test.pval

## H04
res <- spm_time_dep(data, verbose=F, 
                    frm = list(at="a", f1t="120", Qt="Q", ft="f", bt="b", mu0t="mu0"),
                    start=list(a=-0.05, f1=80, Q=1e-8, f=90, b=5, mu0=0.001), 
                    lb=list(a=-1, f1=30, Q=1e-9, f=10, b=1, mu0=1e-6), 
                    opts = list(algorithm = "NLOPT_LN_NELDERMEAD", 
                    maxeval = 200, ftol_rel = 1e-12), lrtest="H04")

res$alternative$lr.test.pval

## H05
res <- spm_time_dep(data, verbose=F, 
                    frm = list(at="-0.1", f1t="f1", Qt="Q", ft="f", bt="b", mu0t="mu0"),
                    start=list(a=-0.05, f1=80, Q=1e-8, f=90, b=5, mu0=0.001), 
                    opts = list(algorithm = "NLOPT_LN_NELDERMEAD", 
                    maxeval = 200, ftol_rel = 1e-12), lrtest="H05")

res$alternative$lr.test.pval

References

Akushevich, I., A. Kulminski, and K. G. Manton. 2005. “Life Tables with Covariates: Dynamic Model for Nonlinear Analysis of Longitudinal Data.” Mathematical Population Studies 12 (2): 51–80.
Arbeev, Konstantin G., Igor Akushevich, Alexander M. Kulminski, Liubov S. Arbeeva, Lucy Akushevich, Svetlana V. Ukraintseva, Irina V. Culminskaya, and Anatoli I. Yashin. 2009. “Genetic Model for Longitudinal Studies of Aging, Health, and Longevity and Its Potential Application to Incomplete Data.” Journal of Theoretical Biology 258 (1): 103–11.
Arbeev, Konstantin G., Igor Akushevich, Alexander M. Kulminski, Svetlana V. Ukraintseva, and Anatoliy I. Yashin. 2014. “Joint Analyses of Longitudinal and Time-to-Event Data in Research on Aging: Implications for Predicting Health and Survival.” Frontiers in Public Health 2 (November).
Arbeev, Konstantin G., Alan A. Cohen, Liubov S. Arbeeva, Emmanuel Milot, Eric Stallard, Alexander M. Kulminski, Igor Akushevich, Svetlana V. Ukraintseva, Kaare Christensen, and Anatoliy I. Yashin. 2016. “Optimal Versus Realized Trajectories of Physiological Dysregulation in Aging and Their Relation to Sex-Specific Mortality Risk.” Frontiers in Public Health 4 (January).
He, Liang, Ilya Zhbannikov, Konstantin G. Arbeev, Anatoliy I. Yashin, and Alexander M. Kulminski. 2017. “A Genetic Stochastic Process Model for Genome-Wide Joint Analysis of Biomarker Dynamics and Disease Susceptibility with Longitudinal Data.” Genetic Epidemiology 41 (7): 620–35.
Woodbury, Max A., and Kenneth G. Manton. 1977. “A Random-Walk Model of Human Mortality and Aging.” Theoretical Population Biology 11 (1): 37–48.
Yashin, Anatoli I., Konstantin G. Arbeev, Igor Akushevich, Aliaksandr Kulminski, Lucy Akushevich, and Svetlana V. Ukraintseva. 2007. “Stochastic Model for Analysis of Longitudinal Data on Aging and Mortality.” Mathematical Biosciences 208 (2): 538–51.
Yashin, Anatoli I., Konstantin G. Arbeev, Aliaksandr Kulminski, Igor Akushevich, Lucy Akushevich, and Svetlana V. Ukraintseva. 2007. “Health Decline, Aging and Mortality: How Are They Related?” Biogerontology 8 (3): 291–302.