PSM.simulate: Create simulation data for multiple individuals

Description Usage Arguments Details Value Note Author(s) References See Also Examples

Description

Simulates data for multiple individuals in a mixed effects model based on stochastic differential equations using an euler scheme.

Usage

1
PSM.simulate(Model, Data, THETA, deltaTime, longX=TRUE) 

Arguments

Model

A list containing the model components either Linear or Non-Linear Model list.*

Data

List with elements described below. No Data$Y is needed as it is generated through the simulation. The number of individuals simulated is equal to length(Data).

Time

Time vector

U

Input list for the Model

covar

Covariates list

THETA

Vector of population parameters

deltaTime

Time Step in the Euler scheme

longX

Boolean. Toggles output of the entire simulated outcome of the states

* See description in PSM.estimate.

Details

The eta is drawn from the multivariate normal distribution N(0,OMEGA). The simulation is an euler based method but for every time interval dt the model is predicted and the states affected by system noise (SIG).

The measurements are added an normal error term belonging to N(0,S).

The function mvrnorm from the MASS pacakge is used to to generate random numbers fra multivariate normal distributions.

Value

The simulated outcome of the model is returned in a list, where each element is the data for an individual.

X

Simulated states sampled at time points for measurements

Y

Simulated measurements

Time

Time points for measurements

U

Input vector used in the simulation

eta

The random effects used in the simulation

Dose

The dose list used in the simulation

longX

Entire outcome of simulated states

longTime

Time points for longX.

Note

For further details please also read the package vignette pdf-document by writing vignette("PSM") in R.

Author(s)

Stig B. Mortensen and Soeren Klim

References

Please visit http://www.imm.dtu.dk/psm or refer to the help page for PSM.

See Also

PSM, PSM.estimate, PSM.smooth, PSM.plot, PSM.template

Examples

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
#specify pharmacokinetic model
#2 state equations, 1 observation equation, 1 random effect

mod = vector(mode="list")
mod$Matrices = function(phi) {
  list(
    matA=matrix(c(-phi$ka, 0, phi$ka, -phi$ke), nrow=2, ncol=2, byrow=TRUE),
    matC=matrix(c(0, 1), nrow=1, ncol=2)
  )
}
mod$h = function(eta, theta, covar) {
  phi = theta
  phi$dose = theta$dose * exp(eta[1])
  phi
}
mod$S = function(phi) {
  matrix(c(phi$sigma), nrow=1, ncol=1)
}
mod$SIG = function(phi) {
  matrix(c(0, 0, 0, phi$omega), nrow=2, ncol=2, byrow=TRUE)
}
mod$X0 = function(Time, phi, U) {
  matrix(c(phi$dose, 0), nrow=2, ncol=1)
}
mod$ModelPar = function(THETA) {
  list(theta=list(dose = THETA["dose"], ka = THETA["ka"], ke = THETA["ke"],
                  omega = THETA["omega"], sigma = THETA["sigma"]),
       OMEGA=matrix(c(THETA["BSV_dose"]), nrow=1, ncol=1)
  )
}


#specify sampling scheme and RNG

TheophPSM <- list()
TheophPSM[[1]] <- list(Time = seq(0,25,5))
set.seed(12345)

#simulate and visualize ODE model (no volatility)

parM <- c(ka = 1.58, ke = 0.08, dose = 9.54, omega = 0, sigma = 1.05, BSV_dose = 0)
TheophSim <- PSM.simulate(mod, TheophPSM, THETA = parM, deltaTime = 0.1)
plot(TheophSim[[1]]$longTime, TheophSim[[1]]$longX[2,],
     type="l", ylab="concentration", xlab="time")

#contrast it to SDE model

parM <- c(ka = 1.58, ke = 0.08, dose = 9.54, omega = 0.34, sigma = 1.05, BSV_dose = 0)
TheophSim <- PSM.simulate(mod, TheophPSM, THETA = parM, deltaTime = 0.1)
lines(TheophSim[[1]]$longTime, TheophSim[[1]]$longX[2,],
     ylab="concentration", xlab="time")

PSM documentation built on May 2, 2019, 9:47 a.m.