Title: | Statistics for Long-Memory Processes (Book Jan Beran), and Related Functionality |
---|---|
Description: | Datasets and Functionality from 'Jan Beran' (1994). Statistics for Long-Memory Processes; Chapman & Hall. Estimation of Hurst (and more) parameters for fractional Gaussian noise, 'fARIMA' and 'FEXP' models. |
Authors: | Jan Beran [aut] (original S functions and scripts), Martin Maechler [cre, aut] (Toplevel R functions and much more, <https://orcid.org/0000-0002-8685-9910>), Brandon Whitcher [ctb] (Datasets) |
Maintainer: | Martin Maechler <[email protected]> |
License: | GPL (>= 2) |
Version: | 1.1-3 |
Built: | 2024-11-08 05:28:30 UTC |
Source: | https://github.com/cran/longmemo |
Compute the covariance matrix of for a fractional
ARIMA process.
CetaARIMA(eta, p, q, m = 10000, delta = 1e-9)
CetaARIMA(eta, p, q, m = 10000, delta = 1e-9)
eta |
parameter vector |
p , q
|
integer scalars giving the AR and MA order respectively. |
m |
integer specifying the length of the Riemann sum, with step
size |
delta |
step size for numerical derivative computation. |
builds on calling specARIMA(eta,p,q,m)
the (square) matrix containg covariances up to ...
Jan Beran (principal) and Martin Maechler (fine tuning)
Beran(1984), listing on p.224–225.
(C.7 <- CetaARIMA(0.7, m = 256, p = 0, q = 0)) (C.5 <- CetaARIMA(eta = c(H = 0.5, phi=c(-.06, 0.42, -0.36), psi=0.776), m = 256, p = 3, q = 1))
(C.7 <- CetaARIMA(0.7, m = 256, p = 0, q = 0)) (C.5 <- CetaARIMA(eta = c(H = 0.5, phi=c(-.06, 0.42, -0.36), psi=0.776), m = 256, p = 3, q = 1))
Covariance matrix of for fractional Gaussian
noise (fGn).
CetaFGN(eta, m = 10000, delta = 1e-9)
CetaFGN(eta, m = 10000, delta = 1e-9)
eta |
parameter vector |
m |
integer specifying the length of the Riemann sum, with step
size |
delta |
step size for numerical derivative computation. |
Currently, the step size for numerical derivative is the same in all
coordinate directions of eta
. In principle, this can be far
from optimal.
Variance-covariance matrix of the estimated parameter vector
.
Jan Beran (principal) and Martin Maechler (speedup, fine tuning)
(C.7 <- CetaFGN(0.7, m = 256)) (C.5 <- CetaFGN(eta = c(H = 0.5), m = 256)) (C.5. <- CetaFGN(eta = c(H = 0.5), m = 1024))
(C.7 <- CetaFGN(0.7, m = 256)) (C.5 <- CetaFGN(eta = c(H = 0.5), m = 256)) (C.5. <- CetaFGN(eta = c(H = 0.5), m = 1024))
Compute the Autocovariances of a fractional ARIMA(0,d,0) process (d = H - 1/2).
ckARMA0(n, H)
ckARMA0(n, H)
n |
sample size (length of time series). |
H |
self-similarity (‘Hurst’) parameter. |
The theoretical formula,
where ,
leads to over-/underflow for larger lags
;
hence use the asymptotical formula there.
numeric vector of length n
of covariances
.
Jan Beran (principal) and Martin Maechler (speedup, fine tuning)
Jan Beran (1994), p.63, (2.35) and (2.39).
ckFGN0
which does the same for fractional
Gaussian noise.
str(C.8 <- ckARMA0(50, H = 0.8)) yl <- c(0,max(C.8)) plot(0:49, C.8, type = "h", ylim = yl) plot(0:49, C.8, type = "h", log = "xy", main = "Log-Log ACF for ARIMA(0,d,0)")
str(C.8 <- ckARMA0(50, H = 0.8)) yl <- c(0,max(C.8)) plot(0:49, C.8, type = "h", ylim = yl) plot(0:49, C.8, type = "h", log = "xy", main = "Log-Log ACF for ARIMA(0,d,0)")
Compute the Autocovariances of a fractional Gaussian process
ckFGN0(n, H)
ckFGN0(n, H)
n |
sample size (length of time series). |
H |
self-similarity (‘Hurst’) parameter. |
numeric vector of covariances upto lag n-1.
Jan Beran (principal) and Martin Maechler (fine tuning)
ckARMA0
which does the same for a fractional
ARIMA process.
str(C.8 <- ckFGN0(50, H = 0.8)) plot(0:49, C.8, type = "h", ylim = 0:1) plot(0:49, C.8, type = "h", log = "xy", main = "Log-Log ACF for frac.GaussNoise(H = 0.8)")
str(C.8 <- ckFGN0(50, H = 0.8)) plot(0:49, C.8, type = "h", ylim = 0:1) plot(0:49, C.8, type = "h", log = "xy", main = "Log-Log ACF for frac.GaussNoise(H = 0.8)")
Ethernet traffic data from a LAN at Bellcore, Morristown (Leland et al. 1993, Leland and Wilson 1991). The data are listed in chronological sequence by row.
data(ethernetTraffic)
data(ethernetTraffic)
A times series of length 4000.
Jan Beran and Brandon Whitcher by E-mail in fall 1995.
data(ethernetTraffic) str(ethernetTraffic) plot(ethernetTraffic)## definitely special
data(ethernetTraffic) str(ethernetTraffic) plot(ethernetTraffic)## definitely special
FEXPest(x, *)
computes Beran's Fractional EXP or ‘FEXP’
model estimator.
.ffreq(n)
returns the Fourier frequencies
(of a time series of length
n
).
FEXPest(x, order.poly, pvalmax, verbose = FALSE) ## S3 method for class 'FEXP' print(x, digits = getOption("digits"), ...) .ffreq(n, full = FALSE)
FEXPest(x, order.poly, pvalmax, verbose = FALSE) ## S3 method for class 'FEXP' print(x, digits = getOption("digits"), ...) .ffreq(n, full = FALSE)
x |
numeric vector representing a time series. |
order.poly |
integer specifying the maximal polynomial order that
should be taken into account. |
pvalmax |
maximal P-value – the other iteration stopping
criterion and “model selection tuning parameter”.
Setting this to |
verbose |
logical indicating if iteration output should be printed. |
digits , ...
|
optional arguments for |
n |
a positive integer, typically the length of a time series. |
full |
logical indicating if |
FEXPest(x,..)
returns an object of class FEXP
, basically a list with components
call |
the function |
n |
time series length |
H |
the “Hurst” parameter which is simply |
coefficients |
numeric 4-column matrix as returned from
|
order.poly |
the effective polynomial order used. |
max.order.poly |
the original |
early.stop |
logical indicating if |
spec |
the spectral estimate |
yper |
raw periodogram of (centered and scaled |
There currently are methods for print()
,
plot
and lines
(see
plot.FEXP
) for objects of class "FEXP"
.
Martin Maechler, using Beran's “main program” in Beran(1994), p.234 ff
Beran, Jan (1993) Fitting long-memory models by generalized linear regression. Biometrika 80, 817–822.
Beran, Jan (1994). Statistics for Long-Memory Processes; Chapman & Hall.
WhittleEst
;
the plot method, plot.FEXP
.
data(videoVBR) (fE <- FEXPest(videoVBR, order = 3, pvalmax = .5)) (fE3 <- FEXPest(videoVBR, order = 3, pvalmax = 1 )) (fE7 <- FEXPest(videoVBR, order = 3, pvalmax = 0.10)) ##--> this also choses order 2, as "FE" : all.equal(fE $coef, fE7$coef) # -> TRUE confint(fE) confint(fE7, level = 0.99) .ffreq(8) .ffreq(8, TRUE) stopifnot(all.equal((1:3)/4, .ffreq(8) / pi))
data(videoVBR) (fE <- FEXPest(videoVBR, order = 3, pvalmax = .5)) (fE3 <- FEXPest(videoVBR, order = 3, pvalmax = 1 )) (fE7 <- FEXPest(videoVBR, order = 3, pvalmax = 0.10)) ##--> this also choses order 2, as "FE" : all.equal(fE $coef, fE7$coef) # -> TRUE confint(fE) confint(fE7, level = 0.99) .ffreq(8) .ffreq(8, TRUE) stopifnot(all.equal((1:3)/4, .ffreq(8) / pi))
Log-Log and “Log-X” plot of spectrum. Very simple utilities, kept here mainly for back compatibility, as they appear in the book scripts.
llplot(yper, spec) lxplot(yper, spec)
llplot(yper, spec) lxplot(yper, spec)
yper |
periodogram values |
spec |
spectrum values |
Jan Beran (principal) and Martin Maechler (speedup, fine tuning)
spectrum()
from standard R (package
stats).
NBS weight measurements - deviation from 1 kg in micrograms, see the references. The data are listed in chronological sequence by row.
data(NBSdiff1kg)
data(NBSdiff1kg)
A time series of length 289.
Jan Beran and Brandon Whitcher by E-mail in fall 1995.
H.P. Graf, F.R. Hampel, and J.Tacier (1984). The problem of unsuspected serial correlations. In J. Franke, W. Härdle, and R.D. Martin, editors, Robust and Nonlinear Time Series Analysis, Lecture Notes in Statistics 26, 127–145; Springer.
Pollak, M., Croakin, C., and Hagwood, C. (1993). Surveillance schemes with applications to mass calibration. NIST report 5158; Gaithersburg, MD.
data(NBSdiff1kg) plot(NBSdiff1kg)
data(NBSdiff1kg) plot(NBSdiff1kg)
Monthly temperature for the northern hemisphere for the years 1854-1989, from the data base held at the Climate Research Unit of the University of East Anglia, Norwich, England. The numbers consist of the temperature (degrees C) difference from the monthly average over the period 1950-1979.
data(NhemiTemp)
data(NhemiTemp)
Time-Series (ts
) of length 1632, frequency 12,
starting 1854, ending 1990.
Jan Beran and Brandon Whitcher by E-mail in fall 1995.
Jones, P.D. and Briffa, K.R. (1992) Global surface air temperature variations during the twentieth century, part 1. The Holocene 2, 165–179.
Jan Beran (1994). Dataset no. 5, p.29–31.
data(NhemiTemp) plot(NhemiTemp) mean(window(NhemiTemp, 1950,1979))# (about) 0 ``by definition''
data(NhemiTemp) plot(NhemiTemp) mean(window(NhemiTemp, 1950,1979))# (about) 0 ``by definition''
Yearly minimal water levels of the Nile river for the years 622 to 1281, measured at the Roda gauge near Cairo, (Tousson, p. 366–385).
data(NileMin)
data(NileMin)
Time-Series (ts
) of length 663.
The original Nile river data supplied by Beran only contained only 500 observations (622 to 1121). However, the book claimed to have 660 observations (622 to 1281). First added the remaining observations from the book by hand, and still came up short with only 653 observations (622 to 1264). Finally have 663 observations : years 622–1284 (as in orig. source)
Tousson, O. (1925). Mémoire sur l'Histoire du Nil; Mémoire de l'Institut d'Egypte.
Jan Beran (1994). Dataset no.1, p.20–22.
data(NileMin) plot(NileMin, main = "Nile River Minima 622 - 1284")
data(NileMin) plot(NileMin, main = "Nile River Minima 622 - 1284")
Simply estimate the periodogram via the Fast Fourier Transform.
per(z)
per(z)
z |
numeric vector with the series to compute the periodogram from. |
This is basically the same asspec.pgram(z, fast = FALSE, detrend = FALSE, taper = 0) $ spec
,
and not really recommended to use — exactly for the reason that
spec.pgram
has the defaults differently,
fast = TRUE, detrend = TRUE, taper = 0.1
, see that help page.
numeric vector of length where
n = length(z)
.
Jan Beran (principal) and Martin Maechler (fine tuning)
a more versatile periodogram estimate by
spec.pgram
.
data(NileMin) plot(10*log10(per(NileMin)), type='l')
data(NileMin) plot(10*log10(per(NileMin)), type='l')
(S3) methods for the generic functions plot
(and
lines
) applied to fractional EXP (FEXP) and
"WhittleEst"
(FEXPest
,
models.
plot()
plots the data periodogram and the ‘FEXP'’ model
estimated spectrum, where lines()
and does the latter.
## S3 method for class 'FEXP' plot(x, log = "xy", type = "l", col.spec = 4, lwd.spec = 2, xlab = NULL, ylab = expression(hat(f)(nu)), main = paste(deparse(x$call)[1]), sub = NULL, ...) ## (With identical argument list:) ## S3 method for class 'WhittleEst' plot(x, log = "xy", type = "l", col.spec = 4, lwd.spec = 2, xlab = NULL, ylab = expression(hat(f)(nu)), main = paste(deparse(x$call)[1]), sub = NULL, ...) ## S3 method for class 'FEXP' lines(x, type = "l", col = 4, lwd = 2, ...) ## S3 method for class 'WhittleEst' lines(x, type = "l", col = 4, lwd = 2, ...)
## S3 method for class 'FEXP' plot(x, log = "xy", type = "l", col.spec = 4, lwd.spec = 2, xlab = NULL, ylab = expression(hat(f)(nu)), main = paste(deparse(x$call)[1]), sub = NULL, ...) ## (With identical argument list:) ## S3 method for class 'WhittleEst' plot(x, log = "xy", type = "l", col.spec = 4, lwd.spec = 2, xlab = NULL, ylab = expression(hat(f)(nu)), main = paste(deparse(x$call)[1]), sub = NULL, ...) ## S3 method for class 'FEXP' lines(x, type = "l", col = 4, lwd = 2, ...) ## S3 method for class 'WhittleEst' lines(x, type = "l", col = 4, lwd = 2, ...)
x |
an R object of |
log |
character specifying log scale should be used, see
|
type |
plot type for the periodogram, see |
col.spec , lwd.spec , col , lwd
|
graphical parameters used for drawing the
estimated spectrum, see |
xlab , ylab , main , sub
|
labels for annotating the plot, see
|
... |
further arguments passed to |
Martin Maechler
FEXPest
, WhittleEst
;
plot.default
and spectrum
.
data(videoVBR) fE <- FEXPest(videoVBR, order = 3, pvalmax = .5) plot(fE) fE3 <- FEXPest(videoVBR, order = 3, pvalmax = 1)#-> order 3 lines(fE3, col = "red3", lty=2) f.GN <- WhittleEst(videoVBR) f.am21 <- WhittleEst(videoVBR, model = "fARIMA", start= list(H= .5, AR = c(.5,0), MA= .5)) lines(f.GN, col = "blue4") lines(f.am21, col = "goldenrod") ##--- Using a tapered periodogram --------- spvVBR <- spec.pgram(videoVBR, fast=FALSE, plot=FALSE) fam21 <- WhittleEst(periodogr.x = head(spvVBR$spec, -1), n = length(videoVBR), model = "fARIMA", start= list(H= .5, AR = c(.5,0), MA= .5)) fam21 f.am21 # similar but slightly different plot(fam21) ## Now, comparing to traditional ("log-X", not "log-log") spectral plot: plot(fam21, log="y") ## compared to the standard R spectral plot : s if(dev.interactive(TRUE)) getOption("device")()# a new graphics window plot(spvVBR, log = "yes", col="gray50") all.equal(.ffreq(fE$n) / (2*pi) -> ffr., head(spvVBR$freq, -1))# TRUE lines(ffr., fam21$spec, col=4, lwd=2) ## need to adjust for different 'theta1': lines(ffr., f.am21$spec * fam21$theta1 / f.am21$theta1, col = adjustcolor("tomato", 0.6), lwd=2)
data(videoVBR) fE <- FEXPest(videoVBR, order = 3, pvalmax = .5) plot(fE) fE3 <- FEXPest(videoVBR, order = 3, pvalmax = 1)#-> order 3 lines(fE3, col = "red3", lty=2) f.GN <- WhittleEst(videoVBR) f.am21 <- WhittleEst(videoVBR, model = "fARIMA", start= list(H= .5, AR = c(.5,0), MA= .5)) lines(f.GN, col = "blue4") lines(f.am21, col = "goldenrod") ##--- Using a tapered periodogram --------- spvVBR <- spec.pgram(videoVBR, fast=FALSE, plot=FALSE) fam21 <- WhittleEst(periodogr.x = head(spvVBR$spec, -1), n = length(videoVBR), model = "fARIMA", start= list(H= .5, AR = c(.5,0), MA= .5)) fam21 f.am21 # similar but slightly different plot(fam21) ## Now, comparing to traditional ("log-X", not "log-log") spectral plot: plot(fam21, log="y") ## compared to the standard R spectral plot : s if(dev.interactive(TRUE)) getOption("device")()# a new graphics window plot(spvVBR, log = "yes", col="gray50") all.equal(.ffreq(fE$n) / (2*pi) -> ffr., head(spvVBR$freq, -1))# TRUE lines(ffr., fam21$spec, col=4, lwd=2) ## need to adjust for different 'theta1': lines(ffr., f.am21$spec * fam21$theta1 / f.am21$theta1, col = adjustcolor("tomato", 0.6), lwd=2)
Qeta()
( of Beran(1994), p.117)
is up to scaling the negative log likelihood function of the specified
model, i.e., fractional Gaussian noise or fractional ARIMA.
Qeta(eta, model = c("fGn","fARIMA"), n, yper, pq.ARIMA, give.B.only = FALSE)
Qeta(eta, model = c("fGn","fARIMA"), n, yper, pq.ARIMA, give.B.only = FALSE)
eta |
parameter vector = (H, phi[1:p], psi[1:q]). |
model |
character specifying the kind model class. |
n |
data length |
yper |
numeric vector of length |
pq.ARIMA |
integer, = c(p,q) specifying models orders of AR and
MA parts — only used when |
give.B.only |
logical, indicating if only the |
Calculation of and
where
,
and the sum is taken over all Fourier frequencies
, (
).
is the spectral density of fractional Gaussian noise or
fractional ARIMA(p,d,q) with self-similarity parameter
(and
AR and
MA parameters in the latter case), and is
computed either by
specFGN
or specARIMA
.
a list with components
n |
= input |
H |
(input) Hurst parameter, = |
eta |
= input |
A , B
|
defined as above. |
Tn |
the goodness of fit test statistic
|
z |
the standardized test statistic |
pval |
the corresponding p-value P(W > z) |
theta1 |
the scale parameter
such that |
spec |
scaled spectral density |
yper[1] must be the periodogram at
the frequency
, i.e., not the frequency zero !
Jan Beran (principal) and Martin Maechler (fine tuning)
Jan Beran (1992). A Goodness-of-Fit Test for Time Series with Long Range Dependence. JRSS B 54, 749–760.
Beran, Jan (1994). Statistics for Long-Memory Processes; Chapman & Hall. (Section 6.1, p.116–119; 12.1.3, p.223 ff)
WhittleEst
computes an approximate MLE for fractional
Gaussian noise / fractional ARIMA, by minimizing Qeta
.
data(NileMin) y <- NileMin n <- length(y) yper <- per(scale(y))[2:(1+ (n-1) %/% 2)] eta <- c(H = 0.3) q.res <- Qeta(eta, n=n, yper=yper) str(q.res)
data(NileMin) y <- NileMin n <- length(y) yper <- per(scale(y))[2:(1+ (n-1) %/% 2)] eta <- c(H = 0.3) q.res <- Qeta(eta, n=n, yper=yper) str(q.res)
Simulation of a Gaussian series . Whereas
simGauss
works from autocovariances, where simFGN0
and
simARMA0
call it,
for simulating a fractional ARIMA(0,d,0) process (),
or fractional Gaussian noise, respectively.
simARMA0 (n, H) simFGN0 (n, H) simFGN.fft(n, H, ...) simGauss(autocov)
simARMA0 (n, H) simFGN0 (n, H) simFGN.fft(n, H, ...) simGauss(autocov)
n |
length of time series |
H |
self-similarity parameter |
... |
optional arguments passed to |
autocov |
numeric vector of auto covariances
|
simGauss
implements the method by Davies and Harte which is
relatively fast using the FFT (fft
) twice.
To simulate ARIMA(p, d, q), (for d in (-1/2, 1,2), you can use
arima.sim(n, model = list(ar= .., ma = ..),
innov= simARMA0(n,H=d+1/2) , n.start = 0)
.
simFGN.fft()
is about twice as fast as simFGN0()
and
uses Paxson's proposal, by default via
B.specFGN(*, k.approx=3, adjust=TRUE)
.
The simulated series , an R object of class
"ts"
, constructed from ts()
.
Jan Beran (original) and Martin Maechler (simGauss
,
speedup, simplication).
simFGN.fft
: Vern Paxson.
Beran (1994), 11.3.3, p.216 f, referring to
Davis, R.B. and Harte, D.S. (1987). Tests for Hurst effect, Biometrika 74, 95–102.
Vern Paxson (1997). Fast, Approximate Synthesis of Fractional Gaussian Noise for Generating Self-Similar Network Traffic; Computer Communications Review 27 5, 5–18.
ckARMA0
on which simARMA0
relies, and
ckFGN0
on which simFGN0
relies.
x1 <- simFGN0(100, 0.7) x2 <- simARMA0(100, 0.7) plot(simFGN0(1000, 0.8)) #- time series plot
x1 <- simFGN0(100, 0.7) x2 <- simARMA0(100, 0.7) plot(simFGN0(1000, 0.8)) #- time series plot
Calculate the spectral density of a fractional ARMA process with standard normal innovations and self-similarity parameter H.
specARIMA(eta, p, q, m)
specARIMA(eta, p, q, m)
eta |
parameter vector |
p , q
|
integers giving AR and MA order respectively. |
m |
sample size determining Fourier frequencies. |
at the Fourier frequencies , (
),
cov(X(t),X(t+k)) = (sigma/(2*pi))*integral(exp(iuk)g(u)du).
— or rather – FIXME –
1. cov(X(t),X(t+k)) = integral[ exp(iuk)f(u)du ]
2. f() = theta1 * f*() ; spec = f*(), and integral[log(f*())] = 0
an object of class "spec"
(see also spectrum
)
with components
freq |
the Fourier frequencies (in |
spec |
the scaled values spectral density |
theta1 |
the scale factor |
pq |
a vector of length two, |
eta |
a named vector |
method |
a character indicating the kind of model used. |
Jan Beran (principal) and Martin Maechler (fine tuning)
Beran (1994) and more, see ....
The spectral estimate for fractional Gaussian noise,
specFGN
.
In general, spectrum
and spec.ar
.
str(r.7 <- specARIMA(0.7, m = 256, p = 0, q = 0)) str(r.5 <- specARIMA(eta = c(H = 0.5, phi=c(-.06, 0.42, -0.36), psi=0.776), m = 256, p = 3, q = 1)) plot(r.7) plot(r.5)
str(r.7 <- specARIMA(0.7, m = 256, p = 0, q = 0)) str(r.5 <- specARIMA(eta = c(H = 0.5, phi=c(-.06, 0.42, -0.36), psi=0.776), m = 256, p = 3, q = 1)) plot(r.7) plot(r.5)
Calculation of the spectral density of
normalized fractional Gaussian noise with self-similarity parameter
at the Fourier frequencies 2*pi*j/m (j=1,...,(m-1)).
B.specFGN
computes (approximations of) the
component of the spectrum
.
specFGN(eta, m, ...) B.specFGN(lambd, H, k.approx=3, adjust = (k.approx == 3), nsum = 200)
specFGN(eta, m, ...) B.specFGN(lambd, H, k.approx=3, adjust = (k.approx == 3), nsum = 200)
eta |
parameter vector |
m |
sample size determining Fourier frequencies. |
... |
optional arguments for |
lambd |
numeric vector of frequencies in [0, pi] |
H |
Hurst parameter in |
k.approx |
either integer (the order of the Paxson approximation), or
|
adjust |
logical indicating (only for |
nsum |
if the slow sum is used (e.g. for k.approx = NA), the number of terms. |
Note that
cov(X(t),X(t+k)) = integral[ exp(iuk)f(u)du ]
f=theta1*spec and integral[log(spec)]=0.
Since longmemo version 1.1-0, a fast approximation is available
(and default), using k.approx
terms and an adjustment
(adjust=TRUE
in the default case of k.approx=3
),
which is due to the analysis and S code from Paxson (1997).
specFGN()
returns an object of class "spec"
(see also
spectrum
) with components
freq |
the Fourier frequencies |
spec |
the scaled values spectral density |
theta1 |
the scale factor |
H |
the self-similarity parameter from input. |
method |
a character indicating the kind of model used. |
B.specFGN()
returns a vector of (approximate) values
.
Jan Beran originally (using the slow sum); Martin Maechler, based on Vern Paxson (1997)'s code.
Jan Beran (1994). Statistics for Long-Memory Processes; Chapman & Hall, NY.
Vern Paxson (1997). Fast, Approximate Synthesis of Fractional Gaussian Noise for Generating Self-Similar Network Traffic; Computer Communications Review 27 5, 5–18.
The spectral estimate for fractional ARIMA,
specARIMA
; more generally, spectrum
.
str(rg.7 <- specFGN(0.7, m = 100)) str(rg.5 <- specFGN(0.5, m = 100))# { H = 0.5 <--> white noise ! } plot(rg.7) ## work around plot.spec() `bug' in R < 1.6.0 plot(rg.5, add = TRUE, col = "blue") text(2, mean(rg.5$spec), "H = 0.5 [white noise]", col = "blue", adj = c(0,-1/4)) ## This was the original method in longmemo, upto version 1.0-0 (incl): rg.7.o <- specFGN(0.7, m = 100, k.approx=NA, nsum = 200) ## quite accurate (but slightly slower): rg.7f <- specFGN(0.7, m = 100, k.approx=NA, nsum = 10000) ## comparing old and new default : all.equal(rg.7, rg.7.o)# different in about 5th digit all.equal(rg.7, rg.7f )# ==> new default is *more* accurate: 1.42 e-6 ## takes about 7 sec {in 2011}: rg.7ff <- specFGN(0.7, m = 100, k.approx=NA, nsum = 500000) all.equal(rg.7f, rg.7ff)# ~ 10 ^ -7 all.equal(rg.7 $spec, rg.7ff$spec)# ~ 1.33e-6 -- Paxson is accurate! all.equal(rg.7.o$spec, rg.7ff$spec)# ~ 2.40e-5 -- old default is less so
str(rg.7 <- specFGN(0.7, m = 100)) str(rg.5 <- specFGN(0.5, m = 100))# { H = 0.5 <--> white noise ! } plot(rg.7) ## work around plot.spec() `bug' in R < 1.6.0 plot(rg.5, add = TRUE, col = "blue") text(2, mean(rg.5$spec), "H = 0.5 [white noise]", col = "blue", adj = c(0,-1/4)) ## This was the original method in longmemo, upto version 1.0-0 (incl): rg.7.o <- specFGN(0.7, m = 100, k.approx=NA, nsum = 200) ## quite accurate (but slightly slower): rg.7f <- specFGN(0.7, m = 100, k.approx=NA, nsum = 10000) ## comparing old and new default : all.equal(rg.7, rg.7.o)# different in about 5th digit all.equal(rg.7, rg.7f )# ==> new default is *more* accurate: 1.42 e-6 ## takes about 7 sec {in 2011}: rg.7ff <- specFGN(0.7, m = 100, k.approx=NA, nsum = 500000) all.equal(rg.7f, rg.7ff)# ~ 10 ^ -7 all.equal(rg.7 $spec, rg.7ff$spec)# ~ 1.33e-6 -- Paxson is accurate! all.equal(rg.7.o$spec, rg.7ff$spec)# ~ 2.40e-5 -- old default is less so
Amount of coded information (variable bit rate) per frame for a certain video sequence. There were about 25 frames per second.
data(videoVBR)
data(videoVBR)
a time-series of length 1000.
Heeke, H. (1991) Statistical multiplexing gain for variable bit rate codecs in ATM networks. Int. J. Digit. Analog. Commun. Syst. 4, 261–268.
Heyman, D., Tabatabai, A., and Lakshman, T.V. (1991) Statistical analysis and simulation of video teleconferencing in ATM networks. IEEE Trans. Circuits. Syst. Video Technol. 2, 49–59.
Jan Beran (1994). Dataset no. 2, p.22–23.
data(videoVBR) plot(log(videoVBR), main="VBR Data (log)")
data(videoVBR) plot(log(videoVBR), main="VBR Data (log)")
Computes Whittle's approximate MLE for fractional Gaussian noise or fractional ARIMA (=: fARIMA) models, according to Beran's prescript.
Relies on minmizing Qeta()
(,
which itself is based on the “true” spectrum of the
corresponding process; for the spectrum, either
specFGN
or specARIMA
is used.
WhittleEst(x, periodogr.x = per(if(scale) x/sd(x) else x)[2:((n+1) %/% 2)], n = length(x), scale = FALSE, model = c("fGn", "fARIMA"), p, q, start = list(H= 0.5, AR= numeric(), MA=numeric()), verbose = getOption("verbose")) ## S3 method for class 'WhittleEst' print(x, digits = getOption("digits"), ...)
WhittleEst(x, periodogr.x = per(if(scale) x/sd(x) else x)[2:((n+1) %/% 2)], n = length(x), scale = FALSE, model = c("fGn", "fARIMA"), p, q, start = list(H= 0.5, AR= numeric(), MA=numeric()), verbose = getOption("verbose")) ## S3 method for class 'WhittleEst' print(x, digits = getOption("digits"), ...)
x |
numeric vector representing a time series. Maybe omitted if
|
periodogr.x |
the (raw) periodogram of |
n |
length of the time series, |
scale |
logical indicating if |
model |
numeric vector representing a time series. |
p , q
|
optional integers specifying the AR and MA orders of the
fARIMA model, i.e., only applicable when |
start |
list of starting values; currently necessary for |
verbose |
logical indicating if iteration output should be printed. |
digits , ...
|
optional arguments for |
An object of class WhittleEst
which is basically a list with components
call |
the function |
model |
= input |
n |
time series length |
p , q
|
for "fARIMA": order of AR and MA parts, respectively. |
coefficients |
numeric 4-column matrix of coefficients
with estimate of the full parameter vector |
theta1 |
the scale parameter |
vcov |
the variance-covariance matrix for |
periodogr.x |
= input (with default). |
spec |
the spectral estimate |
There is a print
method, and coef
,
confint
or vcov
methods work as well for
objects of class "WhittleEst"
.
Martin Maechler, based on Beran's “main program” in Beran(1994).
Beran, Jan (1994). Statistics for Long-Memory Processes; Chapman & Hall. (Section 6.1, p.116–119; 12.1.3, p.223 ff)
Qeta
is the function minimized by these Whittle
estimators.
FEXPest
for an alternative model with Hurst parameter,
also estimated by a “Whittle” approximate MLE, i.e., a
Whittle's estimator in the more general sense.
The plot method, plot.WhittleEst
.
data(NileMin) (f.Gn.N <- WhittleEst(NileMin)) # H = 0.837 (f.A00.N <- WhittleEst(NileMin, model = "fARIMA", p=0, q=0)) # H = 0.899 confint(f.Gn.N) confint(f.A00.N) data(videoVBR) (f.GN <- WhittleEst(videoVBR)) ## similar {but faster !}: (f.am00 <- WhittleEst(videoVBR, model = "fARIMA", p=0, q=0)) rbind(f.am00$coef, f.GN $coef)# really similar f.am11 <- WhittleEst(videoVBR, model = "fARIMA", start= list(H= .5, AR = .5, MA= .5)) f.am11 vcov(f.am11) op <- if(require("sfsmisc")) mult.fig(3, main = "Whittle Estimators for videoVBR data")$old.par else par(mar = c(3,1), mgp = c(1.5, 0.6, 0), mar = c(4,4,2,1)+.1) plot(f.GN) plot(f.am00) plot(f.am11) et <- as.list(coef(f.am11)) et$AR <- c(et$AR, 0, 0) # two more AR coefficients .. f.am31 <- WhittleEst(videoVBR, model = "fARIMA", start = et) ## ... warning non nonconvergence, but "kind of okay": lines(f.am31, col = "red3") ## drawing on top of ARMA(1,1) above - *small* diff f.am31 # not all three are "significant" round(cov2cor(vcov(f.am31)), 3) # and they are highly correlated et <- as.list(coef(f.am31)) et$AR <- unname(unlist(et[c("AR1", "AR2")])) f.am21 <- WhittleEst(videoVBR, model = "fARIMA", start = c(et[c("H","AR", "MA")])) f.am21 lines(f.am21, col = adjustcolor("gold", .4), lwd=4) par(op)## (reset graphic layout) ##--> ?plot.WhittleEst for an example using 'periodogr.x'
data(NileMin) (f.Gn.N <- WhittleEst(NileMin)) # H = 0.837 (f.A00.N <- WhittleEst(NileMin, model = "fARIMA", p=0, q=0)) # H = 0.899 confint(f.Gn.N) confint(f.A00.N) data(videoVBR) (f.GN <- WhittleEst(videoVBR)) ## similar {but faster !}: (f.am00 <- WhittleEst(videoVBR, model = "fARIMA", p=0, q=0)) rbind(f.am00$coef, f.GN $coef)# really similar f.am11 <- WhittleEst(videoVBR, model = "fARIMA", start= list(H= .5, AR = .5, MA= .5)) f.am11 vcov(f.am11) op <- if(require("sfsmisc")) mult.fig(3, main = "Whittle Estimators for videoVBR data")$old.par else par(mar = c(3,1), mgp = c(1.5, 0.6, 0), mar = c(4,4,2,1)+.1) plot(f.GN) plot(f.am00) plot(f.am11) et <- as.list(coef(f.am11)) et$AR <- c(et$AR, 0, 0) # two more AR coefficients .. f.am31 <- WhittleEst(videoVBR, model = "fARIMA", start = et) ## ... warning non nonconvergence, but "kind of okay": lines(f.am31, col = "red3") ## drawing on top of ARMA(1,1) above - *small* diff f.am31 # not all three are "significant" round(cov2cor(vcov(f.am31)), 3) # and they are highly correlated et <- as.list(coef(f.am31)) et$AR <- unname(unlist(et[c("AR1", "AR2")])) f.am21 <- WhittleEst(videoVBR, model = "fARIMA", start = c(et[c("H","AR", "MA")])) f.am21 lines(f.am21, col = adjustcolor("gold", .4), lwd=4) par(op)## (reset graphic layout) ##--> ?plot.WhittleEst for an example using 'periodogr.x'