Stem.Kriging: Dynamical spatial mapping

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

Description

This functions performs spatial prediction in a set of new S spatial locations for a fixed time point.

Usage

1
2
Stem.Kriging(StemModel, coord.newlocations, covariates.newlocations,
K.newlocations, time.point, cov.spat = Sigmastar.exp)

Arguments

StemModel

an object of class “Stem.Model” given as output by the Stem.Estimation function.

coord.newlocations

a matrix or a data frame of dimension S by 2.

covariates.newlocations

a matrix of dimension S by r, where r is the number of covariates as in the StemModel object.
It has the same structure of the StemModel$data$covariates object in the sense that each station data set is stacked under the others; see the DETAILS of the Stem.Model function.

K.newlocations

a loading matrix of dimension S by p.

time.point

the time point between 1 and n for which the spatial prediction is performed.

cov.spat

type of spatial covariance function. For the moment only the exponential function is implemented.

Details

Given the observation matrix and using the multivariate Normal distribution standard theory, the predictor in the new generic spatial location s_0 at time t is an univariate Gaussian distribution with mean z(s_0,t) and variance σ^2(s_0) given by:

z(s_0,t) = X(s_0,t)*β+K(s_0)*y_t+Ω '*Σ_e^{-1}*(z_t-X_t*β+K*y_t)

σ^2(s_0)=σ^2{ω}-Ω '*Σ_e^{-1}*Ω

where Ω is the d*1 constant in time covariance vector, whose i-th generic element (i=1,...,d) is Cov(z(s_i,t),z(s_0,t)). Moreover, X(s_0,t) is the 1*r vector of covariates for the new site s_0 and K(s_0) is the 1*p loading vector. Note that all the parameters in the previous formula are ML estimates and the latent process y_t is the output of the Kalman filtering procedures for each time point t.

Value

The function returns a list which is given by:

data.newlocations

a list of five objects related to the new spatial locations: the coordinates (coordinates), the covariates (covariates), the K matrix, the predictions z and the prediction standard errors (se.pred).

time.point

the time point for which the spatial prediction is performed.

Author(s)

Michela Cameletti michela.cameletti@unibg.it

References

Amisigo, B.A., Van De Giesen, N.C. (2005) Using a spatio-temporal dynamic state-space model with the EM algorithm to patch gaps in daily riverflow series. Hydrology and Earth System Sciences 9, 209–224.

Fasso', A., Cameletti, M., Nicolis, O. (2007) Air quality monitoring using heterogeneous networks. Environmetrics 18, 245–264.

Fasso', A., Cameletti, M. (2007) A general spatio-temporal model for environmental data. Tech.rep. n.27 Graspa - The Italian Group of Environmental Statistics - http://www.graspa.org.

Fasso', A., Cameletti, M. (2009) A unified statistical approach for simulation, modelling, analysis and mapping of environmental data. Accepted for publication by Simulation: transaction of the Society for Modeling and Simulation International.

Mc Lachlan, G.J., Krishnan, T. (1997) The EM Algorithm and Extensions. Wiley, New York.

Shumway, R.H., Stoffer, D.S. (2006) Time Series Analysis and Its Applications: with R Examples. Springer, New York.

Xu, K., Wikle, C.K. (2007) Estimation of parameterized spatio-temporal dynamic models. Journal of Statistical Inference and Planning 137, 567–588.

See Also

See Also pm10, Stem.Model and Stem.Estimation

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
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
#load the data
data(pm10)

#extract the data
coordinates <- pm10$coords
covariates <- pm10$covariates
z <- pm10$z

#build the parameter list 
#(the phi list is used for the algorithm starting values)
phi <- list(beta=matrix(c(3.65,0.046,-0.904),3,1),
			sigma2eps=0.1,
			sigma2omega=0.2,
			theta=0.01,
			G=matrix(0.77,1,1),
			Sigmaeta=matrix(0.3,1,1),
			m0=as.matrix(0),
			C0=as.matrix(1))

K<-matrix(1,ncol(z),1)

mod1 <- Stem.Model(z=z,covariates=covariates,coordinates=coordinates,phi=phi,K=K)
class(mod1)
is.Stem.Model(mod1)

#mod1 is given as output by the Stem.Model function
mod1.est <- Stem.Estimation(mod1)

#coordinates of the 25 new points displaced in a regular grid (S=25)
xxx  <- seq(400,470,length=5)
yyy <- seq(5000,5070,length=5)
coord.new <- expand.grid(x=xxx,y=yyy)

#plot of the spatial locations
plot(pm10$coords[,1],pm10$coords[,2],xlab=colnames(pm10$coords)[1],
ylab=colnames(pm10$coords)[2])
points(coord.new[,1],coord.new[,2],col=2,pch=19)
legend("topleft",col=c(1,2),lty=c(0,0), pch=c(21,19),
legend=c("Original spatial locations","New spatial locations"))

#the covariates matrix for the new 25 spatial locations for the 10th time point
covariates.new <- cbind(rep(1,25),  
c(37.98348, 18.14824, 15.32287, 11.00458, 6.67696,
29.120820, 10.487590, 2.401088, 26.112971, 1.683525,
19.211907, 31.363448, 3.629172, 10.352472, 48.289624,
7.199692, 3.524810, 25.546621, 19.598600, 10.521586,
0.004736363, 0.365510044, 0.975484255, 25.523642458, 4.671496566),
c(0.227688, 0.173037, 0.139985, 0.116392, 0.102476,
0.278325, 0.256422, 0.168136, 0.129460, 0.121040,
0.722656, 0.238780, 0.202586, 0.166547, 0.154638,
0.733208, 1.467990, 0.380001, 0.251896, 0.240350,
 2.292299 ,2.275844 ,1.382322, 0.300729, 0.208798))

K.new<-matrix(1,25,1)

#dynamical spatial prediction (10th day)
mod1.pred <-Stem.Kriging(StemModel=mod1.est,coord.newlocations=coord.new,
			covariates.newlocations=covariates.new,
K.newlocations<-K.new,time.point=10)

#post-processing: build an image map
image(x=xxx,y=yyy,z=matrix(mod1.pred$data.newlocations$z,
length(xxx),length(yyy)),
xlab=colnames(pm10$coords)[1],ylab=colnames(pm10$coords)[2],
xlim=range(mod1.est$data$coordinates[,1])+5,
ylim=range(mod1.est$data$coordinates[,2]+5))
points(pm10$coords[,1],pm10$coords[,2])
points(coord.new[,1],coord.new[,2],col=2,pch=19)

byline <- min((range(xxx)[2]-range(xxx)[1])/4,(range(yyy)[2]-range(yyy)[1])/4)
abline(v=seq(range(xxx)[1],range(xxx)[2],by=byline),col="grey",lty=2)
abline(h=seq(range(yyy)[1],range(yyy)[2],by=byline),col="grey",lty=2)

Example output

Loading required package: mvtnorm
Loading required package: MASS
[1] "Stem.Model"
[1] TRUE
****************EM Algorithm - iteration n. 1 
***NR Algorithm - iteration n. 1 
***NR Algorithm - iteration n. 2 
***NR Algorithm - iteration n. 3 
***NR Algorithm - iteration n. 4 
***NR Algorithm - iteration n. 5 
****************EM Algorithm - iteration n. 2 
***NR Algorithm - iteration n. 1 
***NR Algorithm - iteration n. 2 
***NR Algorithm - iteration n. 3 
***NR Algorithm - iteration n. 4 
****************EM Algorithm - iteration n. 3 
***NR Algorithm - iteration n. 1 
***NR Algorithm - iteration n. 2 
***NR Algorithm - iteration n. 3 
***NR Algorithm - iteration n. 4 
***NR Algorithm - iteration n. 5 
****************EM Algorithm - iteration n. 4 
***NR Algorithm - iteration n. 1 
***NR Algorithm - iteration n. 2 
****************EM Algorithm - iteration n. 5 
***NR Algorithm - iteration n. 1 
***NR Algorithm - iteration n. 2 
****************EM Algorithm - iteration n. 6 
***NR Algorithm - iteration n. 1 
***NR Algorithm - iteration n. 2 
****************EM Algorithm - iteration n. 7 

Stem documentation built on May 2, 2019, 8:56 a.m.