Hessian Matrix via the Harvey (1989) Recursion
MARSSharveyobsFI.Rd
Calculates the observed Fisher Information analytically via the recursion by Harvey (1989) as adapted by Holmes (2017) for MARSS models with linear constraints. This is the same as the Hessian of the negative log-likelihood function at the MLEs. This is a utility function in the MARSS-package
and is not exported. Use MARSShessian()
to access.
Value
The observed Fisher Information matrix computed via equation 3.4.69 in Harvey (1989). The differentials in the equation are computed in the recursion in equations 3.4.73a to 3.4.74b. See Holmes (2016c) for a discussion of the Harvey (1989) algorithm and Holmes (2017) for the specific implementation of the algorithm for MARSS models with linear constraints.
Harvey (1989) discusses missing observations in section 3.4.7. However, the MARSSharveyobsFI()
function implements the approach of Shumway and Stoffer (2006) in section 6.4 for the missing values. See Holmes (2012) for a full discussion of the missing values modifications.
References
R. H. Shumway and D. S. Stoffer (2006). Section 6.4 (Missing Data Modifications) in Time series analysis and its applications. Springer-Verlag, New York.
Harvey, A. C. (1989) Section 3.4.5 (Information matrix) in Forecasting, structural time series models and the Kalman filter. Cambridge University Press, Cambridge, UK.
See also J. E. Cavanaugh and R. H. Shumway (1996) On computing the expected Fisher information matrix for state-space model parameters. Statistics & Probability Letters 26: 347-355. This paper discusses the Harvey (1989) recursion (and proposes an alternative).
Holmes, E. E. (2012). Derivation of the EM algorithm for constrained and unconstrained multivariate autoregressive state-space (MARSS) models. Technical Report. arXiv:1302.3919 [stat.ME]
Holmes, E. E. 2016c. Notes on computing the Fisher Information matrix for MARSS models. Part III Overview of Harvey 1989. https://eeholmes.github.io/posts/2016-6-16-FI-recursion-3/
Holmes, E. E. 2017. Notes on computing the Fisher Information matrix for MARSS models. Part IV Implementing the Recursion in Harvey 1989. https://eeholmes.github.io/posts/2017-5-31-FI-recursion-4/
Examples
dat <- t(harborSeal)
dat <- dat[c(2, 11), ]
fit <- MARSS(dat)
#> Success! abstol and log-log tests passed at 26 iterations.
#> Alert: conv.test.slope.tol is 0.5.
#> Test with smaller values (<0.1) to ensure convergence.
#>
#> MARSS fit is
#> Estimation method: kem
#> Convergence test: conv.test.slope.tol = 0.5, abstol = 0.001
#> Estimation converged in 26 iterations.
#> Log-likelihood: 11.74016
#> AIC: -9.480311 AICc: -6.3692
#>
#> Estimate
#> R.diag 0.0115
#> U.X.CoastalEstuaries 0.0613
#> U.X.OR.NorthCoast 0.0510
#> Q.(X.CoastalEstuaries,X.CoastalEstuaries) 0.0147
#> Q.(X.OR.NorthCoast,X.OR.NorthCoast) 0.0122
#> x0.X.CoastalEstuaries 7.3823
#> x0.X.OR.NorthCoast 6.2707
#> Initial states (x0) defined at t=0
#>
#> Standard errors have not been calculated.
#> Use MARSSparamCIs to compute CIs and bias estimates.
#>
MARSS:::MARSSharveyobsFI(fit)
#> R.diag U.X.CoastalEstuaries
#> R.diag 52914.2357125 -154.19136
#> U.X.CoastalEstuaries -154.1913603 1662.00142
#> U.X.OR.NorthCoast 169.2513722 0.00000
#> Q.(X.CoastalEstuaries,X.CoastalEstuaries) 12164.4113628 120.60426
#> Q.(X.OR.NorthCoast,X.OR.NorthCoast) 13737.8120129 0.00000
#> x0.X.CoastalEstuaries 1.9507559 68.10553
#> x0.X.OR.NorthCoast -0.7477657 0.00000
#> U.X.OR.NorthCoast
#> R.diag 169.2514
#> U.X.CoastalEstuaries 0.0000
#> U.X.OR.NorthCoast 2332.1767
#> Q.(X.CoastalEstuaries,X.CoastalEstuaries) 0.0000
#> Q.(X.OR.NorthCoast,X.OR.NorthCoast) -159.5791
#> x0.X.CoastalEstuaries 0.0000
#> x0.X.OR.NorthCoast 82.0963
#> Q.(X.CoastalEstuaries,X.CoastalEstuaries)
#> R.diag 12164.411363
#> U.X.CoastalEstuaries 120.604263
#> U.X.OR.NorthCoast 0.000000
#> Q.(X.CoastalEstuaries,X.CoastalEstuaries) 16821.992564
#> Q.(X.OR.NorthCoast,X.OR.NorthCoast) 0.000000
#> x0.X.CoastalEstuaries -1.525828
#> x0.X.OR.NorthCoast 0.000000
#> Q.(X.OR.NorthCoast,X.OR.NorthCoast)
#> R.diag 13737.8120129
#> U.X.CoastalEstuaries 0.0000000
#> U.X.OR.NorthCoast -159.5790963
#> Q.(X.CoastalEstuaries,X.CoastalEstuaries) 0.0000000
#> Q.(X.OR.NorthCoast,X.OR.NorthCoast) 23236.9701088
#> x0.X.CoastalEstuaries 0.0000000
#> x0.X.OR.NorthCoast 0.7050328
#> x0.X.CoastalEstuaries
#> R.diag 1.950756
#> U.X.CoastalEstuaries 68.105533
#> U.X.OR.NorthCoast 0.000000
#> Q.(X.CoastalEstuaries,X.CoastalEstuaries) -1.525828
#> Q.(X.OR.NorthCoast,X.OR.NorthCoast) 0.000000
#> x0.X.CoastalEstuaries 44.922134
#> x0.X.OR.NorthCoast 0.000000
#> x0.X.OR.NorthCoast
#> R.diag -0.7477657
#> U.X.CoastalEstuaries 0.0000000
#> U.X.OR.NorthCoast 82.0962990
#> Q.(X.CoastalEstuaries,X.CoastalEstuaries) 0.0000000
#> Q.(X.OR.NorthCoast,X.OR.NorthCoast) 0.7050328
#> x0.X.CoastalEstuaries 0.0000000
#> x0.X.OR.NorthCoast 21.5655793