tsmooth {TSSS}R Documentation

Prediction and Interpolation of Time Series

Description

Predict and interpolate time series based on state space model by Kalman filter.

Usage

tsmooth(y, f, g, h, q, r, x0 = NULL, v0 = NULL, filter.end = NULL,
        predict.end = NULL, minmax = c(-1.0e+30, 1.0e+30), missed = NULL,
        np = NULL, plot = TRUE, ...)

Arguments

y

a univariate time series yny_n.

f

state transition matrix FnF_n.

g

matrix GnG_n.

h

matrix HnH_n.

q

system noise variance QnQ_n.

r

observational noise variance RR.

x0

initial state vector X(00)X(0\mid0).

v0

initial state covariance matrix V(00)V(0\mid0).

filter.end

end point of filtering.

predict.end

end point of prediction.

minmax

lower and upper limits of observations.

missed

start position of missed intervals.

np

number of missed observations.

plot

logical. If TRUE (default), mean vectors of the smoother and estimation error are plotted.

...

graphical arguments passed to plot.smooth.

Details

The linear Gaussian state space model is

xn=Fnxn1+Gnvn,x_n = F_n x_{n-1} + G_n v_n,

yn=Hnxn+wn,y_n = H_n x_n + w_n,

where yny_n is a univariate time series, xnx_n is an mm-dimensional state vector.

FnF_n, GnG_n and HnH_n are m×mm \times m, m×km \times k matrices and a vector of length mm , respectively. QnQ_n is k×kk \times k matrix and RnR_n is a scalar. vnv_n is system noise and wnw_n is observation noise, where we assume that E(vn,wn)=0E(v_n, w_n) = 0, vnN(0,Qn)v_n \sim N(0, Q_n) and wnN(0,Rn)w_n \sim N(0, R_n). User should give all the matrices of a state space model and its parameters. In current version, FnF_n, GnG_n, HnH_n, QnQ_n, RnR_n should be time invariant.

Value

An object of class "smooth", which is a list with the following components:

mean.smooth

mean vectors of the smoother.

cov.smooth

variance of the smoother.

esterr

estimation error.

llkhood

log-likelihood.

aic

AIC.

References

Kitagawa, G. (2020) Introduction to Time Series Modeling with Applications in R. Chapman & Hall/CRC.

Kitagawa, G. and Gersch, W. (1996) Smoothness Priors Analysis of Time Series. Lecture Notes in Statistics, No.116, Springer-Verlag.

Examples

## Example of prediction (AR model)
data(BLSALLFOOD)
BLS120 <- BLSALLFOOD[1:120]
z1 <- arfit(BLS120, plot = FALSE)
tau2 <- z1$sigma2

# m = maice.order, k=1
m1 <- z1$maice.order
arcoef <- z1$arcoef[[m1]]
f <- matrix(0.0e0, m1, m1)
f[1, ] <- arcoef
if (m1 != 1)
  for (i in 2:m1) f[i, i-1] <- 1
g <- c(1, rep(0.0e0, m1-1))
h <- c(1, rep(0.0e0, m1-1))
q <- tau2[m1+1]
r <- 0.0e0
x0 <- rep(0.0e0, m1)
v0 <- NULL

s1 <- tsmooth(BLS120, f, g, h, q, r, x0, v0, filter.end = 120, predict.end = 156)
s1

plot(s1, BLSALLFOOD)

## Example of interpolation of missing values (AR model)
z2 <- arfit(BLSALLFOOD, plot = FALSE)
tau2 <- z2$sigma2

# m = maice.order, k=1
m2 <- z2$maice.order
arcoef <- z2$arcoef[[m2]]
f <- matrix(0.0e0, m2, m2)
f[1, ] <- arcoef
if (m2 != 1)
  for (i in 2:m2) f[i, i-1] <- 1
g <- c(1, rep(0.0e0, m2-1))
h <- c(1, rep(0.0e0, m2-1))
q <- tau2[m2+1]
r <- 0.0e0
x0 <- rep(0.0e0, m2)
v0 <- NULL

tsmooth(BLSALLFOOD, f, g, h, q, r, x0, v0, missed = c(41, 101), np = c(30, 20))

[Package TSSS version 1.3.4-5 Index]