A pomp
object implements a partially observed Markov process (POMP) model.
Basic operations on this model (with shorthand terms) include:
simulation of the state process given parameters (rprocess
)
evaluation of the likelihood of a given state trajectory given parameters (dprocess
)
simulation of the observation process given the states and parameters (rmeasure
)
evaluation of the likelihood of a set of observations given the states and parameters (dmeasure
)
simulation from the prior probability distribution (rprior
)
evaluation of the prior probability density (dprior
)
simulation from the distribution of initial states, given parameters (init.state
)
evaluation of the deterministic skeleton at a point in state space, given parameters (skeleton
)
computation of a trajectory of the deterministic skeleton given parameters (trajectory
)
pomp provides S4 methods that implement each of these basic operations.
These operations can be combined to implement statistical inference methods that depend only on a model's POMP structure.
For convenience, parameter transformations may also be enclosed in a pomp
object.
This page documents these elements.
# S4 method for pomp
rprocess(object, xstart, times, params, offset = 0, …)
# S4 method for pomp
dprocess(object, x, times, params, log = FALSE, …)
# S4 method for pomp
rmeasure(object, x, times, params, …)
# S4 method for pomp
dmeasure(object, y, x, times, params, log = FALSE, …)
# S4 method for pomp
dprior(object, params, log = FALSE, …)
# S4 method for pomp
rprior(object, params, …)
# S4 method for pomp
init.state(object, params, t0, nsim, …)
# S4 method for pomp
skeleton(object, x, t, params, …)
# S4 method for pomp
trajectory(object, params, times, t0, as.data.frame = FALSE, …,
verbose = getOption("verbose", FALSE))
# S4 method for pomp
pompLoad(object, …)
# S4 method for pomp
pompUnload(object, …)
an object of class pomp
.
an nvar
x nrep
matrix containing the starting state of the system.
Columns of xstart
correspond to states; rows to components of the state vector.
One independent simulation will be performed for each column.
Note that in this case, params
must also have nrep
columns.
a rank-3 array containing states of the unobserved process.
The dimensions of x
are nvars
x nrep
x ntimes
, where nvars
is the number of state variables, nrep
is the number of replicates, and ntimes
is the length of times
.
a matrix containing observations.
The dimensions of y
are nobs
x ntimes
, where nobs
is the number of observables and ntimes
is the length of times
.
a numeric vector (length ntimes
) containing times.
These must be in non-decreasing order.
a npar
x nrep
matrix of parameters.
Each column is an independent parameter set and is paired with the corresponding column of x
or xstart
.
In the case of init.state
, params
is a named vector of parameters.
integer;
the first offset
times in times
will not be returned.
the initial time at which initial states are requested.
optional integer; the number of initial states to simulate.
By default, this is equal to the number of columns of params
.
if TRUE, log probabilities are returned.
logical; if TRUE
, return the result as a data-frame.
In trajectory
, additional arguments are passed to the ODE integrator (if the skeleton is a vectorfield) and ignored if it is a map.
See ode
for a description of the additional arguments accepted.
In all other cases, additional arguments are ignored.
logical; if TRUE
, more information will be displayed.
rprocess
simulates the process-model portion of partially-observed Markov process.
When rprocess
is called, the first entry of times
is taken to be the initial time
(i.e., that corresponding to xstart
).
Subsequent times are the additional times at which the state of the simulated processes are required.
rprocess
returns a rank-3 array with rownames.
Suppose x
is the array returned.
Then
dim(x)=c(nvars,nrep,ntimes-offset),
where nvars
is the number of state variables (=nrow(xstart)
), nrep
is the number of independent realizations simulated (=ncol(xstart)
), and ntimes
is the length of the vector times
.
x[,j,k]
is the value of the state process in the j
-th realization at time times[k+offset]
.
The rownames of x
must correspond to those of xstart
.
dprocess
evaluates the probability density of a sequence of consecutive state transitions.
dprocess
returns a matrix of dimensions nrep
x ntimes-1
.
If d
is the returned matrix, d[j,k]
is the likelihood of the transition from state x[,j,k-1]
at time times[k-1]
to state x[,j,k]
at time times[k]
.
rmeasure
simulate the measurement model given states and parameters.
rmeasure
returns a rank-3 array of dimensions nobs
x nrep
x ntimes
, where nobs
is the number of observed variables.
dmeasure
evaluates the probability density of observations given states.
dmeasure
returns a matrix of dimensions nreps
x ntimes
.
If d
is the returned matrix, d[j,k]
is the likelihood of the observation y[,k]
at time times[k]
given the state x[,j,k]
.
dprior
evaluates the prior probability density and rprior
simulates from the prior.
init.state
returns an nvar
x nsim
matrix of state-process initial conditions when given an npar
x nsim
matrix of parameters, params
, and an initial time t0
.
By default, t0
is the initial time defined when the pomp
object ws constructed.
If nsim
is not specified, then nsim=ncol(params)
.
The method skeleton
evaluates the deterministic skeleton at a point or points in state space, given parameters.
In the case of a discrete-time system, the skeleton is a map.
In the case of a continuous-time system, the skeleton is a vectorfield.
NB: skeleton
just evaluates the deterministic skeleton;
it does not iterate or integrate.
skeleton
returns an array of dimensions nvar
x nrep
x ntimes
.
If f
is the returned matrix, f[i,j,k]
is the i-th component of the deterministic skeleton at time times[k]
given the state x[,j,k]
and parameters params[,j]
.
trajectory
computes a trajectory of the deterministic skeleton of a Markov process.
In the case of a discrete-time system, the deterministic skeleton is a map and a trajectory is obtained by iterating the map.
In the case of a continuous-time system, the deterministic skeleton is a vector-field; trajectory
uses the numerical solvers in deSolve to integrate the vectorfield.
trajectory
returns an array of dimensions nvar
x nrep
x ntimes
.
If x
is the returned matrix, x[i,j,k]
is the i-th component of the state vector at time times[k]
given parameters params[,j]
.
When the skeleton is a vectorfield, trajectory
integrates it using ode
.
When the skeleton is a map, trajectory
iterates it.
By default, time is advanced 1 unit per iteration.
The user can change this behavior by specifying the desired timestep using the argument skelmap.delta.t
in the construction of the pomp
object.
User-defined parameter transformations enclosed in the pomp
object can be accessed via partrans
.
pompLoad
and pompUnload
cause compiled codes associated with object
to be dynamically linked or unlinked, respectively.
When C snippets are used in the construction of a pomp
object, the resulting shared-object library is dynamically loaded (linked) before each use, and unloaded afterward.
These functions are provided because in some instances, greater control may be desired.
These functions have no effect on shared-object libraries linked by the user.
# NOT RUN {
pompExample(ricker)
p <- parmat(c(r=42,c=1,phi=10,sigma=0.3,N_0=7,e.0=0),10)
t <- c(1:10,20,30)
t0 <- 0
x0 <- init.state(ricker,params=p,t0=t0)
x <- rprocess(ricker,xstart=x0,times=c(t0,t),params=p,offset=1)
y <- rmeasure(ricker,params=p,x=x,times=t)
ll <- dmeasure(ricker,y=y[,3,,drop=FALSE],x=x,times=t,params=p,log=TRUE)
apply(ll,1,sum)
f <- skeleton(ricker,x=x,t=t,params=p)
z <- trajectory(ricker,params=p,times=t,t0=t0)
## short arguments are recycled:
p <- c(r=42,phi=10,c=1,sigma=0.3,N_0=7,e.0=0)
t <- c(1:10,20,30)
t0 <- 0
x0 <- init.state(ricker,params=p,t0=t0)
x <- rprocess(ricker,xstart=x0,times=c(t0,t),params=p,offset=1)
y <- rmeasure(ricker,params=p,x=x,times=t)
ll <- dmeasure(ricker,y=y,x=x,times=t,params=p,log=TRUE)
f <- skeleton(ricker,x=x,t=t,params=p)
z <- trajectory(ricker,params=p,times=t,t0=t0)
# }
Run the code above in your browser using DataLab