Title: | Methods for Fitting Network Time Series Models |
---|---|
Description: | Simulation of, and fitting models for, Generalised Network Autoregressive (GNAR) time series models which take account of network structure, potentially with exogenous variables. Such models are described in Knight et al. (2020) <doi:10.18637/jss.v096.i05> and Nason and Wei (2021) <doi:10.1111/rssa.12875>. Diagnostic tools for GNAR(X) models can be found in Nason et al. (2023) <doi:10.48550/arXiv.2312.00530>. |
Authors: | Kathryn Leeming [aut], Guy Nason [aut], Matt Nunes [aut, cre], Marina Knight [ctb], James Wei [aut], Daniel Salnikov [aut], Mario Cortina Borja [ctb] |
Maintainer: | Matt Nunes <[email protected]> |
License: | GPL-2 |
Version: | 1.1.4 |
Built: | 2024-12-02 06:53:28 UTC |
Source: | CRAN |
Produces an active node matrix heat-map, which compares the local impact each node has on all the other ones (i.e., regressing on
) once a model order has been chosen. The local relevance indes is
which is closer to one the more relevant
is when forecasting
.
active_node_plot(vts, network, max_lag, r_stages)
active_node_plot(vts, network, max_lag, r_stages)
vts |
Vector time series under study. |
network |
GNAR network object, which is the underlying network for the time series under study. |
max_lag |
Maximum lag of the fitted GNAR model - i.e., |
r_stages |
Neighbourhood regression oreder of the fitted GNAR model - i.e., |
Produces the local influence matrix heat-map for a specific model order. Does not return any values.
Daniel Salnikov and Guy Nason
Nason, G.P., Salnikov, D. and Cortina-Borja, M. (2023) New tools for network time series with an application to COVID-19 hospitalisations. https://arxiv.org/abs/2312.00530
# # Produces an active node heat-map matrix from a stationary GNAR(2, [2, 1]) simulation. # gnar_simulation <- GNARsim(n = 100, net=fiveNet, alphaParams = list(rep(0.25, 5), rep(0.12, 5)), betaParams = list(c(0.25, 0.13), c(0.20)), sigma=1) # # Active node plot # active_node_plot(gnar_simulation, fiveNet, 2, c(2, 1))
# # Produces an active node heat-map matrix from a stationary GNAR(2, [2, 1]) simulation. # gnar_simulation <- GNARsim(n = 100, net=fiveNet, alphaParams = list(rep(0.25, 5), rep(0.12, 5)), betaParams = list(c(0.25, 0.13), c(0.20)), sigma=1) # # Active node plot # active_node_plot(gnar_simulation, fiveNet, 2, c(2, 1))
Function calculating AIC for GNARfit
models.
## S3 method for class 'GNARfit' AIC(object, ..., k=2)
## S3 method for class 'GNARfit' AIC(object, ..., k=2)
object |
a |
... |
additional arguments, not used here. |
k |
the penalty for the criterion, the default |
Smaller AIC values correspond to better fit.
A numeric value corresponding to the AIC
(or other criterion if k
is set to something other than 2
). Note that the value returned is the “time-normalised” AIC for the GNAR model, and also removes any proportionality constants in the calculation.
Occasionally it has been observed that when the forecast horizon has been set too high, this can result in NA
AIC values. This can be resolved by reducing the forecast horizon.
If package users find this problem persists or they experience any other unexpected behaviour, please contact the package maintainer.
#AIC for two different GNAR fits for fiveNet data #GNAR(2,[1,1]) AIC(GNARfit()) #GNAR(2,[1,0]) AIC(GNARfit(betaOrder=c(1,0)))
#AIC for two different GNAR fits for fiveNet data #GNAR(2,[1,1]) AIC(GNARfit()) #GNAR(2,[1,0]) AIC(GNARfit(betaOrder=c(1,0)))
Takes an input GNARnet and neighbour stage and outputs the corresponding adjacency matrix.
## S3 method for class 'GNARnet' as.matrix(x, stage=1, normalise=FALSE, ...)
## S3 method for class 'GNARnet' as.matrix(x, stage=1, normalise=FALSE, ...)
x |
the network |
stage |
the neighbour set that the adjacency matrix is created for. |
normalise |
whether to normalise each to non-zero row to have sum one. |
... |
additional arguments, unused here. |
S3 method for class "GNARnet".
With normalisation this is a non-invertible transform. See NofNeighbours for neighbour set definition.
as.matrix
performed on a GNARnet
returns a square matrix with the number of rows and columns equal to the length of the $edges
list. Entry i,j
of the matrix will be non-zero if node j
is in the stage
neighbour set of i
.
#fiveNet as an adjacency matrix as.matrix(fiveNet)
#fiveNet as an adjacency matrix as.matrix(fiveNet)
Function calculating BIC for GNARfit
models.
## S3 method for class 'GNARfit' BIC(object, ...)
## S3 method for class 'GNARfit' BIC(object, ...)
object |
a |
... |
additional arguments, not used here. |
Smaller BIC values correspond to better fit.
A numeric value corresponding to the BIC
. Note that this is the “time-normalised” value of the AIC for the GNAR model, and also removes any proportionality constants in the calculation.
#BIC for two different GNAR fits for fiveNet data #GNAR(2,[1,1]) BIC(GNARfit()) #GNAR(2,[1,0]) BIC(GNARfit(betaOrder=c(1,0)))
#BIC for two different GNAR fits for fiveNet data #GNAR(2,[1,1]) BIC(GNARfit()) #GNAR(2,[1,0]) BIC(GNARfit(betaOrder=c(1,0)))
coef.GNARfit
returns the vector of coefficients from a GNARfit object.
## S3 method for class 'GNARfit' coef(object,...)
## S3 method for class 'GNARfit' coef(object,...)
object |
the output of a GNARfit call |
... |
additional arguments, unused here. |
S3 method for class "GNARfit".
coef.GNARfit
returns a vector of coefficient values.
#get the coefficients of the fiveNode data GNAR fit coef(GNARfit())
#get the coefficients of the fiveNode data GNAR fit coef(GNARfit())
Plots the GNAR network autcorrelation funciton for a choice of maximum lag and maximum r-stage depth in the network. Using the nacf
function for network autocorrelation and pnacf
for partial network autocorrelation.
corbit_plot(vts, net, max_lag, max_stage, weight_matrix, viridis_color_option="viridis", size_option="absolute_val", partial="no", r_corbit="no", line_segment="no", rectangular_plot="no")
corbit_plot(vts, net, max_lag, max_stage, weight_matrix, viridis_color_option="viridis", size_option="absolute_val", partial="no", r_corbit="no", line_segment="no", rectangular_plot="no")
vts |
Vector time series observations for which one wishes to plot the network autocorrelation or partial network autocorrelation. |
net |
GNAR network object linked to the time series under study. |
max_lag |
Maximum lag the Corbit plot produces (i.e., number of time-steps considered for the network autocorrelaiton.) |
max_stage |
Maximum r-stage depth considered for the Corbit plot (i.e., the number of rings in the plot). Corresponds to the length of paths in the underlying network. |
weight_matrix |
A matrix which entries correspond to the weights between nodes. If this term is NULL, then this argument is equal weights between r-stage neighbours. |
viridis_color_option |
Colour scale for the Corbit plot. The default option is |
size_option |
Point size scale for the Corbit plot. The default is the absolute value of the network autocorrelation function (i.e., |
partial |
Option for selecting between computing the network autocorrelation function or the partial network autocorrelation funciton. Default choice is network autocorrelation (i.e., partial="no"), change argument to "yes" for computing the partial network autocorrelation funciton. |
r_corbit |
Choice for distinguishing between Corbit and R-Corbit plots, default is set to Corbit (inner function call). For producing R-Corbit plots one should use |
line_segment |
Default is set to no (i.e., |
rectangular_plot |
Option for producing alternative rectangular plots. Default is set to no (i.e., |
This function computes network autocorrelation or partial network autocorrelation function values for a specific choice of maximum lag and r-stage depth, and produces the corresponing Corbit plot. Each point in a Corbit plot corresponds to the (P)NACF value, i.e., (p)nacf(h, r)
, at a h-lag and r-stage pair. The ring number starting from the inside corresponds to r-stage depth (path length), and the numbers on the outside ring indicate time lag. The colour scale is based on the overall network autocorrelation values (i.e., the colour is set to highlight strong correlations). The dot at the centre has (p)nacf(h, r)=0
and the smallest size for aiding comparison. Please see the RMarkdown document mentioned in GNAR
references for plot examples.
Produces the specified, i.e, NACF or PNACF for a choice of lag and -stage depth,
, Corbit plot. Does not print (P)NACF values, these are stored as an invisble data frame (matrix), and can be accessed by printing or calling the object produced by the
corbit_plot
call.
Daniel Salnikov and Guy Nason.
Nason, G.P., Salnikov, D. and Cortina-Borja, M. (2023) New tools for network time series with an application to COVID-19 hospitalisations. https://arxiv.org/abs/2312.00530
# # Simulate 100 observations from a stationary GNAR(2, [2, 1]), where # fiveNet is the underlying network. # gnar_simulation <- GNARsim(n = 100, net=fiveNet, alphaParams = list(rep(0.25, 5), rep(0.12, 5)), betaParams = list(c(0.25, 0.13), c(0.20)), sigma=1) # We produce the corresponding Corbit plots. corbit_plot(gnar_simulation, fiveNet, 20, 3) corbit_plot(gnar_simulation, fiveNet, 20, 3, partial = "yes") # If the network object comes with its own weights, then these can be added by including # the option weight_matrix in the corbit call. # corbit_plot(vts, net, max_lag, max_stage, weight_matrix = object_weights_matrix)
# # Simulate 100 observations from a stationary GNAR(2, [2, 1]), where # fiveNet is the underlying network. # gnar_simulation <- GNARsim(n = 100, net=fiveNet, alphaParams = list(rep(0.25, 5), rep(0.12, 5)), betaParams = list(c(0.25, 0.13), c(0.20)), sigma=1) # We produce the corresponding Corbit plots. corbit_plot(gnar_simulation, fiveNet, 20, 3) corbit_plot(gnar_simulation, fiveNet, 20, 3, partial = "yes") # If the network object comes with its own weights, then these can be added by including # the option weight_matrix in the corbit call. # corbit_plot(vts, net, max_lag, max_stage, weight_matrix = object_weights_matrix)
fitted.GNARfit
returns the fitted values of a GNARfit object as a matrix.
## S3 method for class 'GNARfit' fitted(object,...)
## S3 method for class 'GNARfit' fitted(object,...)
object |
the output of a GNARfit call |
... |
additional arguments, unused here. |
S3 method for class "GNARfit".
fitted.GNARfit
returns a ts object of fitted values, with t-alphaOrder
rows and nnodes
columns.
#get the fitted values of the fiveNode GNAR fit fitted(GNARfit())
#get the fitted values of the fiveNode GNAR fit fitted(GNARfit())
A multivariate time series fiveVTS
and corresponding network fiveNet
.
data("fiveNode")
data("fiveNode")
This dataset contains two R objects:fiveVTS
is a ts object with a matrix of 200 rows (t=200) and 5 columns (n=5)
fiveNet
is a GNARnet object containing $edges
and $dist
. edges
is a list of length five, with edges[[i]]
containing the vertices that node i is connected to. dist
is a list of length five, with dist[[i]]
containing the length of the vertices that node i is connected to.
plot(fiveNet) image(fiveVTS)
plot(fiveNet) image(fiveVTS)
This dataset is from the OECD (OECD (2018), Quarterly GDP (indicator). <doi:10.1787/b86d1fc8-en> (Accessed on 29 January 2018)) and is differenced annual growth rate for 35 countries for 1962-2013.
gdpVTS
gdpVTS
gdpVTS
is a ts object with a matrix of 52 rows (t=52) and 35 columns (n=35)
#Plot using 'ts' S3 function, can only plot up to 10 columns at once plot(gdpVTS[,1:5]) #Plot as heatmap image(gdpVTS)
#Plot using 'ts' S3 function, can only plot up to 10 columns at once plot(gdpVTS[,1:5]) #Plot as heatmap image(gdpVTS)
Computes a list of r-stage adjacency matrices, each matrix in the list inidicates whether or not nodes and
are r-stage neighbours in the underlying network. Essentially
if and only if
get_k_stages_adjacency_tensor(St_1, r)
get_k_stages_adjacency_tensor(St_1, r)
St_1 |
One-stage adjacency matrix (i.e., the adjacency matrix of the underlying network). |
r |
Maximum r-stage for which one wishes to compute the r-stage adjacency matrix. |
List containing the adjacency matrices in ascending order.
\eqn{ \{ \mathbf{S}_q \}_{q=1}^{q = r} } |
Each entry is the r-stage adjacency matrix at depth |
Daniel Salnikov and Guy Nason
Nason, G.P., Salnikov, D. and Cortina-Borja, M. (2023) New tools for network time series with an application to COVID-19 hospitalisations. https://arxiv.org/abs/2312.00530
# # Produce the r-stage adjacency tensors for the fiveNet network. # get_k_stages_adjacency_tensor(as.matrix(GNARtoigraph(fiveNet)), 3) #
# # Produce the r-stage adjacency tensors for the fiveNet network. # get_k_stages_adjacency_tensor(as.matrix(GNARtoigraph(fiveNet)), 3) #
A package to fit, predict, and simulate time series using the Generalised Network AutoRegressive (GNAR) model, potentially with exogenous variables. The main functions are GNARfit and GNARXfit, which fits the model to a time series and network(s), S3 method predict.GNARfit which predicts from a fitted GNAR model, and GNARsim which simulates from a GNAR model with specified parameters. For details of the model, see GNARfit. The package also contains an example network time series in data file fiveNode, with network fiveNet, and simulated time series fiveVTS.
Knight, M.I., Nunes, M.A. and Nason, G.P. (2015) Modelling, detrending and decorrelation of network time series. arXiv preprint.
Knight, M.I., Leeming, K., Nason, G.P. and Nunes, M. A. (2020) Generalised Network Autoregressive Processes and the GNAR package. Journal of Statistical Software, 96 (5), 1–36.
Nason G.P. and Wei J. (2022) Quantifying the economic response to COVID-19 mitigations and death rates via forecasting Purchasing Managers’ Indices using Generalised Network Autoregressive models with exogenous variables. Journal of the Royal Statistical Society Series A, 185 (4), 1778–1792.
Creates the design matrix necessary for fitting the GNAR model.
GNARdesign(vts = GNAR::fiveVTS, net = GNAR::fiveNet, alphaOrder = 2, betaOrder = c(1,1), fact.var = NULL, globalalpha=TRUE, tvnets=NULL, netsstart=NULL)
GNARdesign(vts = GNAR::fiveVTS, net = GNAR::fiveNet, alphaOrder = 2, betaOrder = c(1,1), fact.var = NULL, globalalpha=TRUE, tvnets=NULL, netsstart=NULL)
vts |
a matrix or ts object containing the multivariate time series to be modelled. The |
net |
the (first) network associated with the time series, containing a list with entries |
alphaOrder |
a non-negative integer specifying the maximum time-lag to model. |
betaOrder |
a vector of length |
fact.var |
a vector of factors indicating which nodes belong to each set with different parameters to be fitted. |
globalalpha |
a TRUE/FALSE value indivating whether to use global alpha parameters. |
tvnets |
a list of additional networks. Currently only NULL (the static network case) is supported. |
netsstart |
a vector of times corresponding to the first time points for each network of |
GNARdesign |
returns a matrix containing |
#Design matrix to fit GNAR(2,[1,1]) to the fiveVTS data GNARdesign()
#Design matrix to fit GNAR(2,[1,1]) to the fiveVTS data GNARdesign()
Fits the GNAR model to the given inputs using GNARdesign
and lm
.
GNARfit(vts=GNAR::fiveVTS, net=GNAR::fiveNet, alphaOrder=2, betaOrder=c(1,1), fact.var=NULL, globalalpha=TRUE, tvnets=NULL, netsstart=NULL, ErrorIfNoNei=TRUE)
GNARfit(vts=GNAR::fiveVTS, net=GNAR::fiveNet, alphaOrder=2, betaOrder=c(1,1), fact.var=NULL, globalalpha=TRUE, tvnets=NULL, netsstart=NULL, ErrorIfNoNei=TRUE)
vts |
a matrix containing the multivariate time series to be modelled. The |
net |
the (first) network associated with the time series, containing a list with entries |
alphaOrder |
a non-negative integer specifying the maximum time-lag to model. |
betaOrder |
a vector of length |
fact.var |
a vector of factors indicating which nodes belong to different sets with different parameters to be fitted. |
globalalpha |
a TRUE/FALSE value indivating whether to use global alpha parameters. |
tvnets |
a list of additional networks. Currently only NULL (the static network case) is supported. |
netsstart |
a vector of times corresponding to the first time points for each network of |
ErrorIfNoNei |
a TRUE/FALSE value indicating whether to stop the function call with an error if betaOrder specifies more neighbour sets than exist in the network. If FALSE the function will continue and some parameters will be NA. |
The GNAR model of order is defined as
where is the maximum time lag,
and
is the maximum stage of neighbour dependence for time lag
,
is the
th stage neighbour set of node
at time
,
is the connection weight between node
and node
at time
if the path corresponds
to covariate
. Here, we consider a sum from one to zero to be zero and
, are assumed to be independent and identically distributed at each node
, with mean zero and variance
.
Currently, only a single network GNAR model can be fitted.
The connection weight,
, is the inverse of the distance between nodes
i
and q
, normalised so that they sum to 1 for each i,t
.
See is.GNARnet for GNARnet
object information and example construction.
mod |
the |
y |
the original response values, with NAs left in. |
dd |
the output of |
frbic |
inputs to other |
Knight, M.I., Nunes, M.A. and Nason, G.P. Modelling, detrending and decorrelation of network time series.
arXiv preprint.
Knight, M.I., Leeming, K., Nason, G.P. and Nunes, M. A. (2020) Generalised Network Autoregressive Processes and the GNAR package. Journal of Statistical Software, 96 (5), 1–36.
#Fit the GNAR(2,[1,1]) model to the fiveVTS data GNARfit() #Convert the residuals to matrix form residToMat(GNARfit())$resid
#Fit the GNAR(2,[1,1]) model to the fiveVTS data GNARfit() #Convert the residuals to matrix form residToMat(GNARfit())$resid
Simulates a GNAR process with Normally distributed innovations.
GNARsim(n=200, net=GNAR::fiveNet, alphaParams=list(c(rep(0.2,5))), betaParams=list(c(0.5)), sigma=1, tvnets=NULL, netsstart=NULL)
GNARsim(n=200, net=GNAR::fiveNet, alphaParams=list(c(rep(0.2,5))), betaParams=list(c(0.5)), sigma=1, tvnets=NULL, netsstart=NULL)
n |
time length of simulation. |
net |
network used for the GNAR simulation. |
alphaParams |
a list containing vectors of auto-regression parameters for each time-lag. |
betaParams |
a list of equal length as |
sigma |
the standard deviation for the innovations. |
tvnets |
Only NULL is currently supported. |
netsstart |
Only NULL is currently supported. |
Parameter lists should not be NULL, set unused parameters to be zero. See GNARfit for model description.
GNARsim
returns the multivariate time series as a ts object, with n
rows and a column for each of the nodes in the network.
Knight, M.I., Nunes, M.A. and Nason, G.P. Modelling, detrending and decorrelation of network time series.
arXiv preprint.
Knight, M.I., Leeming, K., Nason, G.P. and Nunes, M. A. (2020) Generalised Network Autoregressive Processes and the GNAR package. Journal of Statistical Software, 96 (5), 1–36.
#Simulate a GNAR(1,[1]) process with the fiveNet network GNARsim()
#Simulate a GNAR(1,[1]) process with the fiveNet network GNARsim()
Takes an input network and neighbour stage and returns it in igraph form.
GNARtoigraph(net=GNAR::fiveNet, stage=1, normalise=FALSE)
GNARtoigraph(net=GNAR::fiveNet, stage=1, normalise=FALSE)
net |
a GNARnet object containing |
stage |
the neighbour set that the adjacency matrix is created for. |
normalise |
whether to normalise each to non-zero row to have sum one. |
With normalisation this is a non-invertible transform. See NofNeighbours for neighbour set definition. See is.GNARnet for GNARnet
object information and example construction.
GNARtoigraph
returns an 'igraph' object with weights as the inverse distances of the input network.
#fiveNet as an igraph object GNARtoigraph()
#fiveNet as an igraph object GNARtoigraph()
Creates the design matrix necessary for fitting the GNAR model.
GNARXdesign(vts = GNAR::fiveVTS, net = GNAR::fiveNet, alphaOrder = 2, betaOrder = c(1,1), fact.var = NULL, globalalpha=TRUE, tvnets=NULL, netsstart=NULL, lambdaOrder=NULL, xvts=NULL)
GNARXdesign(vts = GNAR::fiveVTS, net = GNAR::fiveNet, alphaOrder = 2, betaOrder = c(1,1), fact.var = NULL, globalalpha=TRUE, tvnets=NULL, netsstart=NULL, lambdaOrder=NULL, xvts=NULL)
vts |
a matrix or ts object containing the multivariate time series to be modelled. The |
net |
the (first) network associated with the time series, containing a list with entries |
alphaOrder |
a non-negative integer specifying the maximum time-lag to model. |
betaOrder |
a vector of length |
fact.var |
a vector of factors indicating which nodes belong to each set with different parameters to be fitted. |
globalalpha |
a TRUE/FALSE value indivating whether to use global alpha parameters. |
tvnets |
a list of additional networks. Currently only NULL (the static network case) is supported. |
netsstart |
a vector of times corresponding to the first time points for each network of |
lambdaOrder |
a vector of the same length as |
xvts |
a list of matrices containing values of the exogenous regressors for each vertex/node. The |
GNARdesign |
returns a matrix containing |
Nason G.P. and Wei J. (2022) Quantifying the economic response to COVID-19 mitigations and death rates via forecasting Purchasing Managers’ Indices using Generalised Network Autoregressive models with exogenous variables. Journal of the Royal Statistical Society Series A, 185, 1778–1792.
set.seed(1) n = 1000 xvts=list() xvts[[1]] = matrix(rnorm(5*n, mean=0, sd=2), nrow=n, ncol=5) xvts[[2]] = matrix(rnorm(5*n, mean=0, sd=2), nrow=n, ncol=5) lambdaParams=list() lambdaParams[[1]] = c(0.5, -0.5) lambdaParams[[2]] = c(0.3, 0.1) # Simulate the GNARX using the exogenous variables xvts with associated parameters lambdaParams Y_data <- GNARXsim(n=n, net=GNAR::fiveNet, alphaParams=list(c(rep(0.2,5))), betaParams=list(c(0.5)), sigma=1, xvts=xvts, lambdaParams=lambdaParams) #Design matrix to fit GNARX(2,[1,1]) to the fiveVTS data Xdesign <- GNARXdesign(vts = Y_data, net = GNAR::fiveNet, globalalpha = TRUE, alphaOrder = 1, betaOrder = 1, xvts = xvts, lambdaOrder = c(1,1)) Xdesign
set.seed(1) n = 1000 xvts=list() xvts[[1]] = matrix(rnorm(5*n, mean=0, sd=2), nrow=n, ncol=5) xvts[[2]] = matrix(rnorm(5*n, mean=0, sd=2), nrow=n, ncol=5) lambdaParams=list() lambdaParams[[1]] = c(0.5, -0.5) lambdaParams[[2]] = c(0.3, 0.1) # Simulate the GNARX using the exogenous variables xvts with associated parameters lambdaParams Y_data <- GNARXsim(n=n, net=GNAR::fiveNet, alphaParams=list(c(rep(0.2,5))), betaParams=list(c(0.5)), sigma=1, xvts=xvts, lambdaParams=lambdaParams) #Design matrix to fit GNARX(2,[1,1]) to the fiveVTS data Xdesign <- GNARXdesign(vts = Y_data, net = GNAR::fiveNet, globalalpha = TRUE, alphaOrder = 1, betaOrder = 1, xvts = xvts, lambdaOrder = c(1,1)) Xdesign
Fits the GNARX model to the given inputs using GNARdesignX
and lm
.
GNARXfit(vts=GNAR::fiveVTS, net=GNAR::fiveNet, alphaOrder=2, betaOrder=c(1,1), fact.var=NULL, globalalpha=TRUE, tvnets=NULL, netsstart=NULL, ErrorIfNoNei=TRUE, xvts=NULL, lambdaOrder=NULL)
GNARXfit(vts=GNAR::fiveVTS, net=GNAR::fiveNet, alphaOrder=2, betaOrder=c(1,1), fact.var=NULL, globalalpha=TRUE, tvnets=NULL, netsstart=NULL, ErrorIfNoNei=TRUE, xvts=NULL, lambdaOrder=NULL)
vts |
a matrix containing the multivariate time series to be modelled. The |
net |
the (first) network associated with the time series, containing a list with entries |
alphaOrder |
a non-negative integer specifying the maximum time-lag to model. |
betaOrder |
a vector of length |
fact.var |
a vector of factors indicating which nodes belong to different sets with different parameters to be fitted. |
globalalpha |
a TRUE/FALSE value indivating whether to use global alpha parameters. |
tvnets |
a list of additional networks. Currently only NULL (the static network case) is supported. |
netsstart |
a vector of times corresponding to the first time points for each network of |
ErrorIfNoNei |
a TRUE/FALSE value indicating whether to stop the function call with an error if betaOrder specifies more neighbour sets than exist in the network. If FALSE the function will continue and some parameters will be NA. |
xvts |
a list of matrices containing values of the exogenous regressors for each vertex/node. The |
lambdaOrder |
a vector of the same length as |
The GNAR model of order is defined as
where is the maximum time lag,
and
is the maximum stage of neighbour dependence for time lag
,
is the
th stage neighbour set of node
at time
,
is the connection weight between node
and node
at time
if the path corresponds
to covariate
. Here, we consider a sum from one to zero to be zero and
, are assumed to be independent and identically distributed at each node
, with mean zero and variance
.
Currently, only a single network GNAR model can be fitted.
The connection weight,
, is the inverse of the distance between nodes
i
and q
, normalised so that they sum to 1 for each i,t
.
See is.GNARnet for GNARnet
object information and example construction.
A GNARX process of order , neighbourhood order vector
of length
and
node-specific time series exogenous variables with order vector
, denoted
GNARX
, is given by
where are assumed to be set of mutually uncorrelated random variables with
mean zero and variance of
.
mod |
the |
y |
the original response values, with NAs left in. |
dd |
the output of |
frbic |
inputs to other |
Knight, M.I., Nunes, M.A. and Nason, G.P. Modelling, detrending and decorrelation of network time series.
arXiv preprint.
Knight, M.I., Leeming, K., Nason, G.P. and Nunes, M. A. (2020) Generalised Network Autoregressive Processes and the GNAR package. Journal of Statistical Software, 96 (5), 1–36.
Nason G.P. and Wei J. (2022) Quantifying the economic response to COVID-19 mitigations and death rates via forecasting Purchasing Managers’ Indices using Generalised Network Autoregressive models with
exogenous variables. Journal of the Royal Statistical Society Series A, 185, 1778–1792.
set.seed(1) n = 1000 xvts=list() xvts[[1]] = matrix(rnorm(5*n, mean=0, sd=2), nrow=n, ncol=5) xvts[[2]] = matrix(rnorm(5*n, mean=0, sd=2), nrow=n, ncol=5) lambdaParams=list() lambdaParams[[1]] = c(0.5, -0.5) lambdaParams[[2]] = c(0.3, 0.1) # Simulate the GNARX using the exogenous variables xvts with associated parameters lambdaParams Y_data <- GNARXsim(n=n, net=GNAR::fiveNet, alphaParams=list(c(rep(0.2,5))), betaParams=list(c(0.5)), sigma=1, xvts=xvts, lambdaParams=lambdaParams) # fit a GNARX to the model model <- GNARXfit(vts = Y_data, net = GNAR::fiveNet,globalalpha = TRUE, alphaOrder = 1, betaOrder = 1, xvts = xvts, lambdaOrder = c(1,1)) # look at the residuals residToMat(model)$resid
set.seed(1) n = 1000 xvts=list() xvts[[1]] = matrix(rnorm(5*n, mean=0, sd=2), nrow=n, ncol=5) xvts[[2]] = matrix(rnorm(5*n, mean=0, sd=2), nrow=n, ncol=5) lambdaParams=list() lambdaParams[[1]] = c(0.5, -0.5) lambdaParams[[2]] = c(0.3, 0.1) # Simulate the GNARX using the exogenous variables xvts with associated parameters lambdaParams Y_data <- GNARXsim(n=n, net=GNAR::fiveNet, alphaParams=list(c(rep(0.2,5))), betaParams=list(c(0.5)), sigma=1, xvts=xvts, lambdaParams=lambdaParams) # fit a GNARX to the model model <- GNARXfit(vts = Y_data, net = GNAR::fiveNet,globalalpha = TRUE, alphaOrder = 1, betaOrder = 1, xvts = xvts, lambdaOrder = c(1,1)) # look at the residuals residToMat(model)$resid
Simulates a GNAR process with Normally distributed innovations.
GNARXsim(n=200, net=GNAR::fiveNet, alphaParams=list(c(rep(0.2,5))), betaParams=list(c(0.5)), sigma=1, tvnets=NULL, netsstart=NULL, xvts=NULL, lambdaParams=NULL)
GNARXsim(n=200, net=GNAR::fiveNet, alphaParams=list(c(rep(0.2,5))), betaParams=list(c(0.5)), sigma=1, tvnets=NULL, netsstart=NULL, xvts=NULL, lambdaParams=NULL)
n |
time length of simulation. |
net |
network used for the GNAR simulation. |
alphaParams |
a list containing vectors of auto-regression parameters for each time-lag. |
betaParams |
a list of equal length as |
sigma |
the standard deviation for the innovations. |
tvnets |
a list of additional networks. Currently only NULL (the static network case) is supported. |
netsstart |
a vector of times corresponding to the first time points for each network of |
xvts |
a list of matrices containing values of the exogenous regressors for each vertex/node.
The |
lambdaParams |
a list containing vectors of parameters associated to effect of the exogenous regressor variables for each time-lag. |
Parameter lists should not be NULL, set unused parameters to be zero. See GNARXfit for model description.
GNARXsim
returns the multivariate time series as a ts object, with n
rows and a column for each of the nodes in the network.
Knight, M.I., Nunes, M.A. and Nason, G.P. Modelling, detrending and decorrelation of network time series.
arXiv preprint.
Knight, M.I., Leeming, K., Nason, G.P. and Nunes, M. A. (2020) Generalised Network Autoregressive Processes and the GNAR package. Journal of Statistical Software, 96 (5), 1–36.
Nason G.P. and Wei J. (2022) Quantifying the economic response to COVID-19 mitigations and death rates via forecasting Purchasing Managers’ Indices using Generalised Network Autoregressive models with
exogenous variables. Journal of the Royal Statistical Society Series A, 185, 1778–1792.
#Simulate a GNARX process with the fiveNet network set.seed(1) n = 1000 xvts=list() xvts[[1]] = matrix(rnorm(5*n, mean=0, sd=2), nrow=n, ncol=5) xvts[[2]] = matrix(rnorm(5*n, mean=0, sd=2), nrow=n, ncol=5) lambdaParams=list() lambdaParams[[1]] = c(0.5, -0.5) lambdaParams[[2]] = c(0.3, 0.1) # Simulate the GNARX using the exogenous variables xvts with associated parameters lambdaParams Y_data <- GNARXsim(n=n, net=GNAR::fiveNet, alphaParams=list(c(rep(0.2,5))), betaParams=list(c(0.5)), sigma=1, xvts=xvts, lambdaParams=lambdaParams) # now try to refit the model model <- GNARXfit(vts = Y_data, net = GNAR::fiveNet,globalalpha = TRUE, alphaOrder = 1, betaOrder = 1, xvts = xvts, lambdaOrder = c(1,1)) model
#Simulate a GNARX process with the fiveNet network set.seed(1) n = 1000 xvts=list() xvts[[1]] = matrix(rnorm(5*n, mean=0, sd=2), nrow=n, ncol=5) xvts[[2]] = matrix(rnorm(5*n, mean=0, sd=2), nrow=n, ncol=5) lambdaParams=list() lambdaParams[[1]] = c(0.5, -0.5) lambdaParams[[2]] = c(0.3, 0.1) # Simulate the GNARX using the exogenous variables xvts with associated parameters lambdaParams Y_data <- GNARXsim(n=n, net=GNAR::fiveNet, alphaParams=list(c(rep(0.2,5))), betaParams=list(c(0.5)), sigma=1, xvts=xvts, lambdaParams=lambdaParams) # now try to refit the model model <- GNARXfit(vts = Y_data, net = GNAR::fiveNet,globalalpha = TRUE, alphaOrder = 1, betaOrder = 1, xvts = xvts, lambdaOrder = c(1,1)) model
Converts an 'igraph' to the GNARnet
form for use as an input to GNAR functions.
igraphtoGNAR(ig)
igraphtoGNAR(ig)
ig |
an 'igraph' object. |
The values in the $dist
list are the reciprocal of the values from the weighted adjacency matrix.
igraphtoGNAR
returns a GNARnet
: a list with elements $edges
and $dist
.
#Convert fiveNet to igraph and back again igraphtoGNAR(GNARtoigraph(fiveNet))
#Convert fiveNet to igraph and back again igraphtoGNAR(GNARtoigraph(fiveNet))
is.GNARfit
returns either TRUE or FALSE according to a series of GNARfit checks.
is.GNARfit(x)
is.GNARfit(x)
x |
the object to be tested |
The is.GNARfit
function checks whether the object passes a series of tests that correspond to it being the output of GNARfit:
Is it a list containing $mod
and $frbic
Does it contain either $y
and $dd
or $ys
and $ds
Is $mod
a lm object
Does $frbic
have the components to calculate the BIC with BIC.GNARfit
is.GNARfit
returns TRUE
or FALSE
corresponding to passing the above tests.
#check that the example fit meets the criteria above is.GNARfit(GNARfit())
#check that the example fit meets the criteria above is.GNARfit(GNARfit())
is.GNARnet
returns either TRUE or FALSE according to a series of GNARnet checks. as.GNARnet
returns a GNARnet object from an input weights matrix, 'igraph' object, or a GNARnet without assigned class.
is.GNARnet(x) as.GNARnet(x)
is.GNARnet(x) as.GNARnet(x)
x |
the network to be tested or object to be converted |
The is.GNARnet
function checks whether the network passes a series of tests:
Is it a list containing $edges and $dist
Are the $edges and $dist lists of the same length
Are each of the elements of $edges the same length as the corresponding $dist element
Do the edges only contain valid entries, 1,...,nnodes (or NULL
)
Is it labelled as GNARnet
class
Are no duplicate edges present
Are all distances positive
Are there no self-loops in the network
The as.GNARnet
function converts igraph objects to GNARnet form, other possible inputs are adjacency matrices, and lists with $edges
and $dist
entries of the correct form.
is.GNARnet
returns TRUE
or FALSE
corresponding to passing the above tests.
as.GNARnet
returns a GNARnet
object.
#check that the example network meets the criteria above is.GNARnet(fiveNet) #convert to igraph and back again as.GNARnet(GNARtoigraph(fiveNet)) #generate a new network with three nodes #edges 1->2, 2->1, 2->3 #dist 1, 2, 1 #note 1->2 and 2->1 are of different lengths threeEdge <- list(c(2), c(1,3), NULL) threeDist <- list(c(1), c(2,1), NULL) threeNet <- list(edges=threeEdge, dist=threeDist) #check if this is a GNARnet is.GNARnet(threeNet) #use as.GNARnet to change the class threeNet <- as.GNARnet(threeNet) #check if this is a GNARnet now is.GNARnet(threeNet)
#check that the example network meets the criteria above is.GNARnet(fiveNet) #convert to igraph and back again as.GNARnet(GNARtoigraph(fiveNet)) #generate a new network with three nodes #edges 1->2, 2->1, 2->3 #dist 1, 2, 1 #note 1->2 and 2->1 are of different lengths threeEdge <- list(c(2), c(1,3), NULL) threeDist <- list(c(1), c(2,1), NULL) threeNet <- list(edges=threeEdge, dist=threeDist) #check if this is a GNARnet is.GNARnet(threeNet) #use as.GNARnet to change the class threeNet <- as.GNARnet(threeNet) #check if this is a GNARnet now is.GNARnet(threeNet)
Produces a local neighbourhood relevance plot based on the distances in the underlying network. The heat-map matrix should reflect clusters if a GNAR model is valid. The size of the clusters depends on the maximum r-stage depth for neighbourhood regression, as gets larger, the clusters grow or intersect and cover more nodes. The relative strength of conditionally correlated nodes is
.
local_relevance_plot(network, r_star) cross_correlation_plot(h, vts)
local_relevance_plot(network, r_star) cross_correlation_plot(h, vts)
network |
GNAR network object, which is the underlying network for the time series under study. |
r_star |
Maximum active r-stage depth for neighbourhood regression. |
h |
The lag in the cross correlation plot. |
vts |
The vector time series to compute the cross correlation plot on. |
Produces the local relevance plot. Does not return any values.
Daniel Salnikov and Guy Nason
Nason, G.P., Salnikov, D. and Cortina-Borja, M. (2023) New tools for network time series with an application to COVID-19 hospitalisations. https://arxiv.org/abs/2312.00530
# # Produces a local relevance plot, which is a heat-map matrix from a stationary # GNAR(1, [1]) simulation. # gnar_simulation <- GNARsim(n = 100, net=fiveNet, alphaParams = list(rep(0.35, 5)), betaParams = list(c(0.25)), sigma=1) # Active node plot local_relevance_plot(fiveNet, 1) # Compare to the cross-correlation plot at one-lag cross_correlation_plot(1, gnar_simulation)
# # Produces a local relevance plot, which is a heat-map matrix from a stationary # GNAR(1, [1]) simulation. # gnar_simulation <- GNARsim(n = 100, net=fiveNet, alphaParams = list(rep(0.35, 5)), betaParams = list(c(0.25)), sigma=1) # Active node plot local_relevance_plot(fiveNet, 1) # Compare to the cross-correlation plot at one-lag cross_correlation_plot(1, gnar_simulation)
Returns the log-likelihood for a GNARfit
object.
## S3 method for class 'GNARfit' logLik(object,...)
## S3 method for class 'GNARfit' logLik(object,...)
object |
A |
... |
Optional additional arguments, not used here. |
S3 method for the GNARfit class. The function returns the value of
where is the time length of the observations,
is the number of nodes,
is the estimated covariance matrix and
is the matrix of residuals.
A logLik
object.
#calculate log-likelihood for fiveNode data #global alphas logLik(GNARfit()) #individual alphas logLik(GNARfit(globalalpha=FALSE))
#calculate log-likelihood for fiveNode data #global alphas logLik(GNARfit()) #individual alphas logLik(GNARfit(globalalpha=FALSE))
This matrix contains a multivariate/vector time series that counts the number of daily admissions to mechanical ventilation beds in one of 140 NHS Trusts from 2nd April 2020 to 27th June 2021.
The dimension of the matrix is 452x140. The dates and Trust ID codes are stored in the dimnames (first and second) respectively.
corbit_plot
,NHSTrustMVCAug120.net
Daniel Salnikov and Guy Nason
UK Coronavirus website https://ukhsa-dashboard.data.gov.uk
Nason, G.P., Salnikov, D. and Cortina-Borja, M. (2023) New tools for network time series with an application to COVID-19 hospitalisations. https://arxiv.org/abs/2312.00530
## Not run: data(logMVbedMVC.vts) data(NHSTrustMVCAug120.net) # # Do a corbit plot with this data, with only three lags and one stage # # Note, normally max_lag and max_stage would be higher than this, the # values are artificially small here, as otherwise the run-time restrictions # for CRAN packaging might be exceeded. corbit_plot(vts=logMVbedMVC.vts, net=NHSTrustMVCAug120.net, max_lag=3, max_stage=1) ## End(Not run)
## Not run: data(logMVbedMVC.vts) data(NHSTrustMVCAug120.net) # # Do a corbit plot with this data, with only three lags and one stage # # Note, normally max_lag and max_stage would be higher than this, the # values are artificially small here, as otherwise the run-time restrictions # for CRAN packaging might be exceeded. corbit_plot(vts=logMVbedMVC.vts, net=NHSTrustMVCAug120.net, max_lag=3, max_stage=1) ## End(Not run)
Converts an adjacency matrix to the GNARnet
form for use as an input to GNAR functions.
matrixtoGNAR(input.mat)
matrixtoGNAR(input.mat)
input.mat |
an adjacency matrix whose dimension equals the number of nodes in the resulting network. |
The values in the $dist
list are the reciprocal of the values from the weighted adjacency matrix. Any self-loops (diagonal entries) and negatively weighted edges are removed.
matrixtoGNAR
returns a GNARnet
list with elements $edges
and $dist
.
#Convert fiveNet to an adjacency matrix and back again matrixtoGNAR(as.matrix(fiveNet))
#Convert fiveNet to an adjacency matrix and back again matrixtoGNAR(as.matrix(fiveNet))
Returns a vector with elements TRUE/FALSE identifying which rows contain NA elements.
na.row(mat)
na.row(mat)
mat |
a matrix object. |
This function is used in the unstacking of residuals into a residual matrix and replacing NAs where they were previously present.
na.row
returns a vector of length equal to the number of rows in mat
. Each element is either TRUE or FALSE.
#Check if there are and NAs in fiveVTS na.row(fiveVTS)
#Check if there are and NAs in fiveVTS na.row(fiveVTS)
Computes the NACF for a choice of lag and r-stage depth
, the NACF is given by
where
is the autocovariance bound of the GNAR process.
nacf(h, s, weight_matrix, stages_tensor, nts_data)
nacf(h, s, weight_matrix, stages_tensor, nts_data)
h |
Lag (i.e., time-steps behind) at which the NACF is computed. |
s |
r-stage depth at which the NACF is computed (i.e., shortest distance between nodes). |
weight_matrix |
Weight matrix of the GNAR process, each entry corresponds to the weight between two nodes; see |
stages_tensor |
List of r-stage adjacency matrices |
nts_data |
Network time series observations, the number of rows is equal to the number of time steps, and the number of columns is equal to the number of series (variables). |
If the network time series contains missing values, then the weights matrix and are adjusted, so that missing values do not contribute to the network autocorrelation. This is done by setting to zero the weights which correspond to a missing value and computing the new weight matrix and
value.
Returns a number between and
, which measures the linear correlation between lagged network observations - i.e.,
.
Daniel Salnikov and Guy Nason.
Nason, G.P., Salnikov, D. and Cortina-Borja, M. (2023) New tools for network time series with an application to COVID-19 hospitalisations. https://arxiv.org/abs/2312.00530
# # Compute the NACF with respect to a stationary GNAR simulation # gnar_simulation <- GNARsim(n = 100, net=fiveNet, alphaParams = list(rep(0.35, 5)), betaParams = list(c(0.25)), sigma=1) W = weights_matrix(fiveNet) stages_list = get_k_stages_adjacency_tensor(as.matrix(GNARtoigraph(fiveNet)), 3) # NACF nacf(3, 1, W, stages_list, gnar_simulation)
# # Compute the NACF with respect to a stationary GNAR simulation # gnar_simulation <- GNARsim(n = 100, net=fiveNet, alphaParams = list(rep(0.35, 5)), betaParams = list(c(0.25)), sigma=1) W = weights_matrix(fiveNet) stages_list = get_k_stages_adjacency_tensor(as.matrix(GNARtoigraph(fiveNet)), 3) # NACF nacf(3, 1, W, stages_list, gnar_simulation)
This matrix contains a multivariate/vector time series that counts the number of daily admissions to mechanical ventilation beds in one of 140 NHS Trusts from 2nd April 2020 to 27th June 2021.
An object of class GNARnet from the GNAR package
Daniel Salnikov and Guy Nason
UK Coronavirus website https://ukhsa-dashboard.data.gov.uk
Nason, G.P., Salnikov, D. and Cortina-Borja, M. (2023) New tools for network time series with an application to COVID-19 hospitalisations. https://arxiv.org/abs/2312.00530
corbit_plot
,NHSTrustMVCAug120.net
## Not run: data(logMVbedMVC.vts) data(NHSTrustMVCAug120.net) # # Plot the network and see what it is like # plot(NHSTrustMVCAug120.net) # # Do a corbit plot with this data, with only three lags and one stage # # Note, normally max_lag and max_stage would be higher than this, the # values are artificially small here, as otherwise the run-time restrictions # for CRAN packaging might be exceeded. corbit_plot(vts=logMVbedMVC.vts, net=NHSTrustMVCAug120.net, max_lag=3, max_stage=1) ## End(Not run)
## Not run: data(logMVbedMVC.vts) data(NHSTrustMVCAug120.net) # # Plot the network and see what it is like # plot(NHSTrustMVCAug120.net) # # Do a corbit plot with this data, with only three lags and one stage # # Note, normally max_lag and max_stage would be higher than this, the # values are artificially small here, as otherwise the run-time restrictions # for CRAN packaging might be exceeded. corbit_plot(vts=logMVbedMVC.vts, net=NHSTrustMVCAug120.net, max_lag=3, max_stage=1) ## End(Not run)
nobs
returns the number of obervations (T) of the input multivariate time series in the GNARfit
function.
## S3 method for class 'GNARfit' nobs(object,...)
## S3 method for class 'GNARfit' nobs(object,...)
object |
the output of a GNARfit or GNARpredict call |
... |
additional arguments, unused here. |
S3 method for class "GNARfit".
An integer specifying the number of rows in the input vts
to the GNARfit function.
#observations of example fiveVTS nobs(GNARfit()) #check this is the same as number of rows in fiveVTS all.equal(nobs(GNARfit()), nrow(fiveVTS))
#observations of example fiveVTS nobs(GNARfit()) #check this is the same as number of rows in fiveVTS all.equal(nobs(GNARfit()), nrow(fiveVTS))
Produces a node relevance plot based on the node relevance index which computes the ratio between nodes
column sum for nodes in neighbourhood regressions. Nodes are ordered according to the relative contribution eahc has to the autocovariance. The nodes are ordered in ascending order.
node_relevance_plot(network, r_star, node_names, node_label_size = 2)
node_relevance_plot(network, r_star, node_names, node_label_size = 2)
network |
GNAR network object, which is the underlying network for the time series under study. |
r_star |
Maximum active r-stage depth for neighbourhood regression. |
node_names |
Names corresponding to each, this makes identifying nodes in the plot easier. If this argument is NULL, then the plot links to each node a number. |
node_label_size |
Text size when producing the plot. Default is 2, however, depending on the number of nodes it might be necessary to adjust the size. |
Data Frame consisting of two variable, the node name and the node relevance value.
Daniel Salnikov and Guy Nason.
Nason, G.P., Salnikov, D. and Cortina-Borja, M. (2023) New tools for network time series with an application to COVID-19 hospitalisations. https://arxiv.org/abs/2312.00530
# # Produces a node relevance plot with respect to a stationary GNAR process # with underlying network fiveNet # # GNAR simulation gnar_simulation <- GNARsim(n = 100, net=fiveNet, alphaParams = list(rep(0.25, 5), rep(0.12, 5)), betaParams = list(c(0.25, 0.13), c(0.20)), sigma=1) # Node relevance plot without names node_relevance_plot(network = fiveNet, r_star = 2, node_label_size = 10) # # Node relevance plot with names # node_relevance_plot(network = fiveNet, r_star = 2, node_names = c("A", "B", "C", "D", "E"), node_label_size = 10)
# # Produces a node relevance plot with respect to a stationary GNAR process # with underlying network fiveNet # # GNAR simulation gnar_simulation <- GNARsim(n = 100, net=fiveNet, alphaParams = list(rep(0.25, 5), rep(0.12, 5)), betaParams = list(c(0.25, 0.13), c(0.20)), sigma=1) # Node relevance plot without names node_relevance_plot(network = fiveNet, r_star = 2, node_label_size = 10) # # Node relevance plot with names # node_relevance_plot(network = fiveNet, r_star = 2, node_names = c("A", "B", "C", "D", "E"), node_label_size = 10)
Calculates neighbour sets of a particular node in the network and their distances.
NofNeighbours(node=1, stage=2, net=GNAR::fiveNet)
NofNeighbours(node=1, stage=2, net=GNAR::fiveNet)
node |
is an integer specifying which node to calculate the neighbours of. |
stage |
is an integer specifying the maximum neighbour-stage to calculate to. |
net |
a |
Note that the distances are calculated as the sum along the shortest path; do not use this with a weights (rather than distance) list. Stage-r
neighbours of node i
are denoted , and are nodes that are
r
edges (but no fewer) away from i
. Hence stage-1 neighbours are the immediate neighbours, stage-2 neighbours are the neighbours of neighbours and so on.
edges |
is a list of length |
dist |
is a list of length |
#First and second stage neighbours of node 1 in fiveNet NofNeighbours()
#First and second stage neighbours of node 1 in fiveNet NofNeighbours()
Plots a GNAR network using the 'igraph' package.
## S3 method for class 'GNARnet' plot(x, ...)
## S3 method for class 'GNARnet' plot(x, ...)
x |
the network |
... |
additional arguments for the |
S3 method for class "GNARnet".
#Plot fiveNet plot(fiveNet)
#Plot fiveNet plot(fiveNet)
Computes the PNACF for a choice of lag and r-stage depth
, the PNACF is given by
where
,
, and
are the empirical residuals corresponding to GNAR(h -1, [r-1, ..., r - 1]) fits,
is the same as for the NACF; see
nacf
, and is the mean of the fitted residuals
pnacf(h, s, weight_matrix, stages_tensor, nts_data)
pnacf(h, s, weight_matrix, stages_tensor, nts_data)
h |
Lag (i.e., time-steps behind) at which the NACF is computed. |
s |
r-stage depth at which the NACF is computed (i.e., shortest distance between nodes). |
weight_matrix |
Weight matrix of the GNAR process, each entry corresponds to the weight between two nodes; see |
stages_tensor |
List of r-stage adjacency matrices |
nts_data |
Network time series observations, the number of rows is equal to the number of time steps, and the number of columns is equal to the number of series (variables). |
If the network time series contains missing values, then the weights matrix and are adjusted, so that missing values do not contribute to the partial network autocorrelation. This is done by setting to zero the weights which correspond to a missing value and computing the new weight matrix and
value.
Daniel Salnikov and Guy Nason
Nason, G.P., Salnikov, D. and Cortina-Borja, M. (2023) New tools for network time series with an application to COVID-19 hospitalisations. https://arxiv.org/abs/2312.00530
# # Compute the PNACF with respect to a stationary GNAR simulation # gnar_simulation <- GNARsim(n = 100, net=fiveNet, alphaParams = list(rep(0.35, 5)), betaParams = list(c(0.25)), sigma=1) W = weights_matrix(fiveNet) stages_list = get_k_stages_adjacency_tensor(as.matrix(GNARtoigraph(fiveNet)), 3) # PNACF pnacf(3, 1, W, stages_list, gnar_simulation)
# # Compute the PNACF with respect to a stationary GNAR simulation # gnar_simulation <- GNARsim(n = 100, net=fiveNet, alphaParams = list(rep(0.35, 5)), betaParams = list(c(0.25)), sigma=1) W = weights_matrix(fiveNet) stages_list = get_k_stages_adjacency_tensor(as.matrix(GNARtoigraph(fiveNet)), 3) # PNACF pnacf(3, 1, W, stages_list, gnar_simulation)
Predicts future observations from a GNARfit object, based on the fitted GNAR model.
## S3 method for class 'GNARfit' predict(object, n.ahead=1, ...)
## S3 method for class 'GNARfit' predict(object, n.ahead=1, ...)
object |
the output of a GNARfit call |
n.ahead |
the time length of the predictions |
... |
further arguments passed to the simulate.GNARfit function, such as |
S3 method for class "GNARfit". This function calls simulate.GNARfit.
A multivariate time series of dimension n.ahead
x nnodes.
#simulate 5 future observations from fiveVTS predict(GNARfit(), n.ahead=5)
#simulate 5 future observations from fiveVTS predict(GNARfit(), n.ahead=5)
print.GNARfit
prints model, call, and coefficients of a GNARfit object.
## S3 method for class 'GNARfit' print(x,...)
## S3 method for class 'GNARfit' print(x,...)
x |
the output of a GNARfit call |
... |
additional arguments, unused here. |
S3 method for class "GNARfit".
#print the information of the fiveNode GNAR fit print(GNARfit())
#print the information of the fiveNode GNAR fit print(GNARfit())
Prints information about a GNAR network.
## S3 method for class 'GNARnet' print(x, ...)
## S3 method for class 'GNARnet' print(x, ...)
x |
the network |
... |
additional arguments, unused here. |
S3 method for class "GNARnet".
#print fiveNet information print(fiveNet)
#print fiveNet information print(fiveNet)
Produces a R-Corbit plot for comparing the network autocorrelation (NACF) and partial network autocorrelation function (PNACF) values for a choice of maximum lag and maximum -stage depth. Starting from the first and continuing to the outermost ring, each ring corresponds to said choice of
-stage depth. The numbers on the outermost ring are time-lags, and each dot corresponds to a specific time-slice or covariate-level.
r_corbit_plot(vts_frames, network_list, max_lag, max_stage, weight_matrices, frame_names, same_net="no", viridis_color_option="viridis", size_option="absolute_val", partial="no", r_corbit="yes")
r_corbit_plot(vts_frames, network_list, max_lag, max_stage, weight_matrices, frame_names, same_net="no", viridis_color_option="viridis", size_option="absolute_val", partial="no", r_corbit="yes")
vts_frames |
List containing the vector time series linked to each of the covariate-levels and/or time-slices, which the R-Corbit plot compares. |
network_list |
List of network objects for which the R-Corbit plot compares network autocorrelation or partial network autocorrelation. |
max_lag |
Maximum lag for the R-Corbit plot. |
max_stage |
Maximum |
weight_matrices |
List of weigth matrices, each weight matrix corresponds to a particular choice of time-slice or covariate-level. If all the time-slices have the same weight matrix, then the argument is a list, where all the entries are equal to the unique weight matrix. |
frame_names |
Indicates the name of each time-slice or covariate-level time series. Order should be the same as in the weight matrices and vector time series lists. |
same_net |
Indicates whether or not all time-slices or covariate-levels share the same weight matrix. Default choice is no, if the time-slices or covariate-levels share the same weight matrix, then this argument should be set to "yes" (i.e., same_net = "yes"). |
viridis_color_option |
Colour scale for the R-Corbit plot. The default option is |
size_option |
Point size scale for the R-Corbit plot. Default is the absolute value of the network autocorrelation function (i.e., |
partial |
Option for selecting between computing the network autocorrelation function or the partial network autocorrelation function. Default choice is network autocorrelation (i.e., partial="no"). Change argument to "yes" for computing the partial network autocorrelation function (PNACF). |
r_corbit |
Choice for distinguishing between Corbit and R-Corbit plots, default is set to Corbit (inner function call). For producing R-Corbit plots one should use |
R-Corbit plots compare the network autocorrelation function (NACF) and partial network autocorrelation function (PNACF) values for a choice of different time-slices and/or covariate-levels. R-Corbit plots are read in the same manner as Corbit plots corbit_plot
, and include a legend on the right-hand side for distinguishing between covariate-levels and/or time-slices. The point at the centre is the mean value of the NACF or PNACF values arising from the time-slices and/or covariate-levels data splits. Essentially, if , where
is the number of covariate-levels or time-slices, then the value at the centre is
where
is the (P)NACF value corresponding to the covariate-level/time-slice
. The number of covariate-levels and time-slices
must be equal to the length of the lists used for producing the R-Corbit plot.
Produces the specified, i.e., NACF or PNACF, values for a choice of lag and -stage depth,
, R-Corbit plot. Does not print (P)NACF values, these are stored as invisble data frames (matrices), and can be accessed by printing or calling the object produced by the
r_corbit_plot
call. The invisible object is a list of matrices, one matrix for each covariate-level/time-slice.
Guy Nason and Daniel Salnikov
Nason, G.P., Salnikov, D. and Cortina-Borja, M. (2023) New tools for network time series with an application to COVID-19 hospitalisations. https://arxiv.org/abs/2312.00530
## Not run: # # Produces a R-Corbit plot, which compares three stationary GNAR simulations, where # the underlying network is fiveNet. # # Compute the weight matrix W = weights_matrix(fiveNet) # # Simulate three stationary GNAR processe sim1 <- GNARsim(n = 100, net=fiveNet, alphaParams = list(c(0.1, 0.12, 0.16, 0.075, 0.21), c(0.12, 0.14, 0.15, 0.6, 0.22)), betaParams = list(c(0.1, 0.16), c(0.11, 0.14))) sim2 <- GNARsim(n = 100, net=fiveNet, alphaParams = list(rep(.25, 5)), betaParams = list(c(0.1, 0.16))) sim3 <- GNARsim(n = 100, net=fiveNet, alphaParams = list(rep(.25, 5), rep(0.13, 5)), betaParams = list(c(0.1, 0.16), c(0.11))) # Produce NACF R-Corbit plot with the same network and weights matrix r_corbit_plot(list(sim1, sim2, sim3), list(fiveNet), 10, 3, list(W), c("sim1", "sim2", "sim3"), same_net = "yes") # # Produce PNACF R-Corbit with different networks and weight matrices print(r_corbit_plot(list(sim1, sim2, sim3), list(fiveNet, fiveNet, fiveNet), 10, 3, list(W, W, W), c("sim1", "sim2", "sim3"), same_net = "no", partial = "yes")) ## End(Not run)
## Not run: # # Produces a R-Corbit plot, which compares three stationary GNAR simulations, where # the underlying network is fiveNet. # # Compute the weight matrix W = weights_matrix(fiveNet) # # Simulate three stationary GNAR processe sim1 <- GNARsim(n = 100, net=fiveNet, alphaParams = list(c(0.1, 0.12, 0.16, 0.075, 0.21), c(0.12, 0.14, 0.15, 0.6, 0.22)), betaParams = list(c(0.1, 0.16), c(0.11, 0.14))) sim2 <- GNARsim(n = 100, net=fiveNet, alphaParams = list(rep(.25, 5)), betaParams = list(c(0.1, 0.16))) sim3 <- GNARsim(n = 100, net=fiveNet, alphaParams = list(rep(.25, 5), rep(0.13, 5)), betaParams = list(c(0.1, 0.16), c(0.11))) # Produce NACF R-Corbit plot with the same network and weights matrix r_corbit_plot(list(sim1, sim2, sim3), list(fiveNet), 10, 3, list(W), c("sim1", "sim2", "sim3"), same_net = "yes") # # Produce PNACF R-Corbit with different networks and weight matrices print(r_corbit_plot(list(sim1, sim2, sim3), list(fiveNet, fiveNet, fiveNet), 10, 3, list(W, W, W), c("sim1", "sim2", "sim3"), same_net = "no", partial = "yes")) ## End(Not run)
Unstacks the entries of the GNARfit fitted and residual values to return matrices of a similar form to the multivariate time series input.
residToMat(GNARobj=GNARfit(), nnodes=5)
residToMat(GNARobj=GNARfit(), nnodes=5)
GNARobj |
the output from the GNARfit function |
nnodes |
the number of nodes in the original network time series |
This function also replaces the NAs that were removed in fitting.
resid |
is the matrix of residual values, with |
fit |
is the matrix of fitted values, with |
#Get residual and fitted matrices from GNARpredict fit of fiveVTS residToMat()
#Get residual and fitted matrices from GNARpredict fit of fiveVTS residToMat()
residuals.GNARfit
returns the residuals of a GNARfit object as a matrix.
## S3 method for class 'GNARfit' residuals(object,...)
## S3 method for class 'GNARfit' residuals(object,...)
object |
the output of a GNARfit call |
... |
additional arguments, unused here. |
The function first checks if the object is of GNARfit class, then uses residToMat to return the residuals.
residuals.GNARfit
returns a 'ts' object of residuals, with t-alphaOrder
rows and nnodes
columns.
#get the residuals of the fiveNode GNAR fit residuals(GNARfit())
#get the residuals of the fiveNode GNAR fit residuals(GNARfit())
Seed numbers for reproducible random graphs.
seed.nos
seed.nos
seed.nos
is a vector of length 10,000 containing integers.
g <- seedToNet(seed.nos[1], nnodes=35, graph.prob=0.15) plot(g, vertex.label=colnames(gdpVTS), vertex.size=0)
g <- seedToNet(seed.nos[1], nnodes=35, graph.prob=0.15) plot(g, vertex.label=colnames(gdpVTS), vertex.size=0)
Produces a reproducible undirected Erdos-Reyni random network using a particular seed value.
seedToNet(seed.no, nnodes=34, graph.prob=0.5)
seedToNet(seed.no, nnodes=34, graph.prob=0.5)
seed.no |
a valid number to set the seed to. |
nnodes |
the number of nodes in the produced network. |
graph.prob |
the probability that each pair of nodes is connected. |
graph.prob
effectively controls the sparsity of the network. All distances are set to 1.
A GNARnet object.
#Generate the random graph from seed 10, with 5 nodes and connection prob 0.5 seedToNet(10,nnodes=5,graph.prob=0.5)
#Generate the random graph from seed 10, with 5 nodes and connection prob 0.5 seedToNet(10,nnodes=5,graph.prob=0.5)
Simulates from a GNARfit object, either creating a new series or future observations of the original series based upon the fitted GNAR model.
## S3 method for class 'GNARfit' simulate(object, nsim=object$frbic$time.in, seed=NULL, future=TRUE, set.noise=NULL, allcoefs=FALSE, ...)
## S3 method for class 'GNARfit' simulate(object, nsim=object$frbic$time.in, seed=NULL, future=TRUE, set.noise=NULL, allcoefs=FALSE, ...)
object |
the output of a GNARfit call |
nsim |
the time length of the simulations |
seed |
either |
future |
whether the simulations follow on from the original time series ( |
set.noise |
a value to set the standard deviation of the noise to, or if |
allcoefs |
if |
... |
additional arguments, unused here. |
S3 method for class "GNARfit".
A multivariate time series of dimension nsim
x nnodes.
#simulate 5 future observations from fiveVTS simulate(GNARfit(), nsim=5)
#simulate 5 future observations from fiveVTS simulate(GNARfit(), nsim=5)
Returns the summary of a GNARfit object, including BIC.
## S3 method for class 'GNARfit' summary(object, ...)
## S3 method for class 'GNARfit' summary(object, ...)
object |
output of a GNARfit call. |
... |
additional arguments, unused here. |
The output is the summary of the fit using summary.lm, and BIC calculated using BIC.GNARfit.
summary.GNARfit
prints the model summary and the value of the BIC.
#summary for the GNAR(2,[1,1]) model using GNARfit on fiveVTS summary(GNARfit())
#summary for the GNAR(2,[1,1]) model using GNARfit on fiveVTS summary(GNARfit())
Prints brief information about a GNAR network.
## S3 method for class 'GNARnet' summary(object, ...)
## S3 method for class 'GNARnet' summary(object, ...)
object |
the network |
... |
additional arguments, unused here. |
S3 method for class "GNARnet".
#print fiveNet summary information summary(fiveNet)
#print fiveNet summary information summary(fiveNet)
Returns the varaince-covariance matrix of the parameters of a GNARfit object.
## S3 method for class 'GNARfit' vcov(object,...)
## S3 method for class 'GNARfit' vcov(object,...)
object |
a GNARfit object, the output from a GNARfit call. |
... |
further arguments passed to vcov. |
S3 method for class "GNARfit".
A matrix of estimated covariances between the parameter estimates, this is calculated using vcov for lm objects.
#covariance matrix of fiveNode fit vcov(GNARfit())
#covariance matrix of fiveNode fit vcov(GNARfit())
A suite of data objects concerning wind speed analysis. The dataset contains a multivariate time series of wind speeds, two network descriptions, a vector of names for weather stations, and the coordinates of the weather stations.
data("vswind")
data("vswind")
This dataset contains six R objects:vswindts
is a ts object with a matrix of 721 rows (t=721) and 102 columns (n=102). This corresponds to 721 observations made through time at
102 weather stations.
vswindnetD
is a GNARnet
object containing $edges
and $dist
. edges
is a list of length 102, with edges[[i]]
containing the vertices that node i is connected to. dist
is a list of length 102, with dist[[i]]
containing the length of the vertices that node i is connected to.
vswindnet
is the same as vswindnetD
except that all the distances
are replaced by 1.
vswindnames
is a character vector of length 102 containing the wind speed
site names and
vswindcoords
is a matrix with 102 rows (one for each wind station) and
two columns providing the x and y coordinates of the weather stations.
The base data were obtained from the http://wow.metoffice.gov.uk UK Met Office WeatherObservationsWebsite distributed under the UK Open Government License https://www.nationalarchives.gov.uk/doc/open-government-licence/version/1/open-government-licence.htm Contains public sector information licensed under the Open Goverment Licence v1.0.
# # The name entry for Bristol # vswindnames[77] #[1] "BRIST" # # plot the distance network # ## Not run: windnetplot()
# # The name entry for Bristol # vswindnames[77] #[1] "BRIST" # # plot the distance network # ## Not run: windnetplot()
Computes the weights matrix with normalised weights (i.e., add up to one) for the network time series with underlying network provided by the user. If the network is unweighted, then each r-stage neighbour is considered equally relevant, i.e., =
, where
is the indicator function and the distance is the shortest path in the underlying network.
weights_matrix(network, max_r_stage)
weights_matrix(network, max_r_stage)
network |
Network linked to the vector time series under study, must be a GNARnet object. |
max_r_stage |
Longest shortest path for which weights are non-zero. If not specified, then its set equal to the upper bound, which is the longest shortest path in the underlying network. |
Weight matrix , each entry is the weight
between a pair of nodes. The matrix is not symmetric, and each row adds up to one when considering r-stage neighbours for a particular r.
Daniel Salnikov and Guy Nason.
Nason, G.P., Salnikov, D. and Cortina-Borja, M. (2023) New tools for network time series with an application to COVID-19 hospitalisations. https://arxiv.org/abs/2312.00530
# # Weights matrix linked to the mechanical ventilation beds time series. # This network has a longest shortest path equal to six. # #data(fiveNet) W_norm = weights_matrix(fiveNet, 6)
# # Weights matrix linked to the mechanical ventilation beds time series. # This network has a longest shortest path equal to six. # #data(fiveNet) W_norm = weights_matrix(fiveNet, 6)
Plots the wind speed data network with distance information.
windnetplot()
windnetplot()
None.
The wind speed data is to be found in the vswind
data set. This function contains commands, using
functionality from the wordcloud
package, to plot the
network, with node names and edges. Distances between nodes are
plotted next to the edges.
## Not run: windplotnet()
## Not run: windplotnet()