fit_stelfi {stelfi}R Documentation

Modelling spatiotemporal self-excitement

Description

Fits spatiotemporal Hawkes models. The self-excitement is Gaussian in space and exponentially decaying in time.

Usage

fit_stelfi(
  times,
  locs,
  sf,
  smesh,
  parameters,
  covariates,
  GMRF = FALSE,
  time_independent = TRUE,
  tmb_silent = TRUE,
  nlminb_silent = TRUE,
  ...
)

Arguments

times

A vector of numeric observed time points.

locs

A data.frame of x and y locations, 2xn.

sf

An sf of type POLYGON specifying the spatial region of the domain.

smesh

A Delaunay triangulation of the spatial domain returned by fmesher::fm_mesh_2d().

parameters

A list of named parameters:

  • coefs, logged base rate of the Hawkes process and coefficients of covariates

  • alpha, intensity jump after an event occurrence

  • beta, rate of exponential decay of intensity after event occurrence

  • tau, \tau parameter for the GMRF (supplied only if GMRF = TRUE)

  • kappa, \kappa parameter for the GMRF (supplied only if GMRF = TRUE)

  • xsigma, standard deviation on x-axis of self-exciting kernel (if time_independent = FALSE, it is the s.d. after 1 time unit)

  • ysigma, standard deviation on y-axis of self-exciting kernel (if time_independent = FALSE, it is the s.d. after 1 time unit)

  • rho, correlation between x and y for the self-exciting kernel (the off-diagonal elements in the kernel's covariate matrix are xsigma * ysigma * rho)

covariates

Optional, a matrix of covariates at each smesh node.

GMRF

Logical, default FALSE. If TRUE, a Gaussian Markov Random Field is included as a latent spatial effect.

time_independent

Logical, default TRUE. If FALSE, Gaussian kernels have a covariate matrix that is proportional to time since the event. Warning, this is very memory intensive.

tmb_silent

Logical, if TRUE (default) then TMB inner optimisation tracing information will be printed.

nlminb_silent

Logical, if TRUE (default) then for each iteration nlminb() output will be printed.

...

Additional arguments to pass to optim()

Details

Temporal self-excitement follows an exponential decay function. The self-excitement over space follows a Gaussian distribution centered at the triggering event. There are two formulations of this model. The default is that the Gaussian function has a fixed spatial covariance matrix, independent of time. Alternatively, covariance can be directly proportional to time, meaning that the self-excitement radiates out from the center over time. This can be appropriate when the mechanism causing self-excitement travels at a finite speed, but is very memory-intensive. The spatiotemporal intensity function used by stelfi is \lambda(s,t) = \mu + \alpha \Sigma_{i:\tau_i<t}(\textrm{exp}(-\beta * (t-\tau_i)) G_i(s-x_i, t - \tau_i)) where

G_i(.,.) can take two forms:

Value

A list containing components of the fitted model, see TMB::MakeADFun. Includes

See Also

fit_hawkes and fit_lgcp

Examples


## No GMRF
if(requireNamespace("fmesher")){
data(xyt, package = "stelfi")
N <- 50
locs <- data.frame(x = xyt$x[1:N], y = xyt$y[1:N])
times <- xyt$t[1:N]
domain <- sf::st_as_sf(xyt$window)
bnd <- fmesher::fm_as_segm(as.matrix(sf::st_coordinates(domain)[, 1:2]))
smesh <- fmesher::fm_mesh_2d(boundary = bnd, max.edge = 0.75, cutoff = 0.3) 
param <- list( mu = 3, alpha = 1, beta = 3, xsigma = 0.2, ysigma = 0.2, rho = 0.8)
fit <- fit_stelfi(times = times, locs = locs, sf = domain, smesh = smesh, parameters = param) 
get_coefs(fit)
## GMRF
param <- list( mu = 5, alpha = 1, beta = 3, kappa = 0.9, tau = 1, xsigma = 0.2,
ysigma = 0.2, rho = 0.8)
fit <- fit_stelfi(times = times, locs = locs, sf = domain, smesh = smesh,
parameters = param, GMRF = TRUE)
get_coefs(fit)
}


[Package stelfi version 1.0.1 Index]