Title: | Spatio-Temporal Analysis of Point Patterns on Linear Networks |
---|---|
Description: | Statistical analysis of spatio-temporal point processes on linear networks. This packages provides tools to visualise and analyse spatio-temporal point patterns on linear networks using first- and second-order summary statistics. |
Authors: | Mehdi Moradi [aut, cre] , Ottmar Cronie [ctb], Jorge Mateu [ctb] |
Maintainer: | Mehdi Moradi <[email protected]> |
License: | GPL (>= 2) |
Version: | 0.4.0 |
Built: | 2024-11-01 04:39:06 UTC |
Source: | https://github.com/moradii/stlnpp |
This function projects an object of class stlpp
into a linear network.
## S3 method for class 'stlpp' as.lpp(x,...)
## S3 method for class 'stlpp' as.lpp(x,...)
x |
an object of class |
... |
arguments passed to |
This function projects the spatio-temporal point pattern x on the linear network L into L, giving its corresponding spatial point pattern.
An object of class lpp
.
Mehdi Moradi <[email protected]>
data(easynet) x <- runifpointOnLines(40, easynet) t1 <- sample(1:10,40,replace=TRUE) Y <- as.stlpp(x,t=t1,L=easynet) as.lpp.stlpp(Y)
data(easynet) x <- runifpointOnLines(40, easynet) t1 <- sample(1:10,40,replace=TRUE) Y <- as.stlpp(x,t=t1,L=easynet) as.lpp.stlpp(Y)
This function converts data to a spatio-temporal point pattern on a linear network.
as.stlpp(x,y,t,L)
as.stlpp(x,y,t,L)
x , y , t
|
vectors of Cartesian coordinates and time occurrence. Alternatively, x can be of classes |
L |
linear network (object of class |
This function converts data to an object of class stlpp. Data can be of formats:
x is of class class data.frame
with three columns. Then columns are considered as Cartesian coordinates (i.e. x,y,t) and they will be converted to a spatio-temporal point pattern on the linear network L.
x is a planar point pattern (class ppp
). Then x will be converted to a spatio-temporal point pattern on the linear network L and with coresponding time vector t.
x is a linear point pattern (class lpp
). Then x will be converted to a spatio-temporal point pattern on the linear network L and with coresponding time vector t.
x,y,t are vectors of same length where x,y are living on the corresponding network L.
A spatio-temporal point pattern on a linear network. An object of class stlpp
.
Mehdi Moradi <[email protected]>
data(easynet) x <- runifpointOnLines(40, easynet) t1 <- sample(1:10,40,replace=TRUE) Y <- as.stlpp(x,t=t1,L=easynet) Z <- as.lpp.stlpp(Y) t2 <- sample(1:10,40,replace=TRUE) W <- as.stlpp(Z,t=t2)
data(easynet) x <- runifpointOnLines(40, easynet) t1 <- sample(1:10,40,replace=TRUE) Y <- as.stlpp(x,t=t1,L=easynet) Z <- as.lpp.stlpp(Y) t2 <- sample(1:10,40,replace=TRUE) W <- as.stlpp(Z,t=t2)
This function converts an object of class stlpp
to class tpp
.
as.tpp.stlpp(X)
as.tpp.stlpp(X)
X |
an object of class |
This function projects the spatio-temporal point pattern X into its corresponding time domain T.
An object of class tpp.
Mehdi Moradi <[email protected]>
X <- rpoistlpp(10,1,2,easynet) as.tpp.stlpp(X)
X <- rpoistlpp(10,1,2,easynet) as.tpp.stlpp(X)
Kernel density estimation of a spatio-temporal point pattern on a linear network.
## S3 method for class 'stlpp' density(x,lbw,tbw,at=c("points","pixels"),dimt=512,...)
## S3 method for class 'stlpp' density(x,lbw,tbw,at=c("points","pixels"),dimt=512,...)
x |
an object of class |
lbw |
network smoothing bandwidth |
tbw |
time smoothing bandwidth |
at |
string specifying whether to compute the intensity values at a grid of pixel locations and times (at="pixels") or only at the points of x (at="points"). default is to estimate the intensity at pixels |
dimt |
the number of equally spaced points at which the temporal density is to be estimated. see density |
... |
arguments passed to |
Kernel smoothing is applied to the spatio-temporal point pattern x using methods in Moradi et al (2019). The function computes estimated intensities assuming first-order separability. Estimated intensity values of the marginal spatial point pattern on the linear network will be obtained using the fast kernel smoothing technique of Rakshit et al. (2019) and function densityQuick.lpp
, whereas the estimated intensity values of the marginal temporal point pattern will be estimated using the function density
.
If lbw and tbw are not given, then they will be selected using bw.nrd0
and bw.scott.iso
respectively.
If at="points"
: a vector of intensity values at the data points of x.
If at="pixels"
: a list of images on linear network. Each image represents an estimated spatio-temporal intensity at a fixed time.
Check the attributes for more accommodated outputs.
Mehdi Moradi <[email protected]>
Moradi, M., & Mateu, J. (2020). First-and second-order characteristics of spatio-temporal point processes on linear networks. Journal of Computational and Graphical Statistics, 29(3), 432-443.
density
, density.lpp
, bw.nrd0
, bw.scott.iso
X <- rpoistlpp(.2,a=0,b=5,L=easynet) density(X)
X <- rpoistlpp(.2,a=0,b=5,L=easynet) density(X)
Kernel estimation of intensity of one-dimensional point patterns.
## S3 method for class 'tpp' density(x,tbw,at=c("points","pixels"),...)
## S3 method for class 'tpp' density(x,tbw,at=c("points","pixels"),...)
x |
an object of class |
tbw |
time smoothing bandwidth |
at |
string specifying whether to compute the intensity values at a grid of pixel locations (at="pixels") or only at the points of x (at="points"). default is to estimate the intensity at pixels |
... |
arguments passed to density |
A vector of intensity values.
If at="points"
: a vector of intensity values at the data points of x.
If at="pixels"
: a vector of intensity values over a grid.
Mehdi Moradi <[email protected]> and Ottmar Cronie
Mateu, J., Moradi, M., & Cronie, O. (2019). Spatio-temporal point patterns on linear networks: Pseudo-separable intensity estimation. Spatial Statistics, 100400.
X <- tpp(sample(c(1:24),200,replace = TRUE)) plot(density(X))
X <- tpp(sample(c(1:24),200,replace = TRUE)) plot(density(X))
This function performs adaptive intensity estimation for spatio-temporal point patterns on linear networks using Voronoi-Dirichlet tessellation.
## S3 method for class 'stlpp' densityVoronoi(X, f = 1, nrep = 1, separable=FALSE, at=c("points","pixels"), dimt=128,...)
## S3 method for class 'stlpp' densityVoronoi(X, f = 1, nrep = 1, separable=FALSE, at=c("points","pixels"), dimt=128,...)
X |
an object of class |
f |
fraction (between 0 and 1 inclusive) of the data points that will be used to build a tessellation for the intensity estimate |
nrep |
number of independent repetitions of the randomised procedure |
separable |
logical. If FALSE, it then calculates a pseudo-separable estimate |
at |
string specifying whether to compute the intensity values at a grid of pixel locations and time (at="pixels") or only at the points of x (at="points"). default is to estimate the intensity at pixels |
dimt |
the number of equally spaced points at which the temporal density is to be estimated. see density |
... |
arguments passed to |
This function computes intensity estimates for spatio-temporal point patterns on linear networks using Voronoi-Dirichlet tessellation. Both first-order separability and pseudo-separability assumptions are accommodated in the function.
If separable=TRUE, the estimated intensities will be a product of the estimated intensities on the network and those on time. Estimated intensity of the spatial component will be obtained using densityVoronoi.lpp
, whereas estimated intensities of the temporal component will be obtained via densityVoronoi.tpp
. If f=1, the function calculates the estimations based on the original Voronoi intensity estimator.
If separable=FALSE, the estimated intensities will be calculated based on a sub-sampling technique explained in Mateu et al. (2019). nrep sub-samples will be obtained from X based on a given retention probability f, the function densityVoronoi.stlpp
, considering separable=TRUE and f=1, will be applied to each obtained sub-sample, and finally, the estimated intensities will be the sum of all obtained estimated intensities from all sub-samples divided by the (f * nrep).
If at="points"
: a vector of intensity values at the data points of X.
If at="pixels"
: a list of images on a linear network. Each image represents an estimated spatio-temporal intensity at a fixed time.
Mehdi Moradi <[email protected]> and Ottmar Cronie
Mateu, J., Moradi, M., & Cronie, O. (2019). Spatio-temporal point patterns on linear networks: Pseudo-separable intensity estimation. Spatial Statistics, 100400.
densityVoronoi.lpp
, density.stlpp
X <- rpoistlpp(.2,a=0,b=5,L=easynet) densityVoronoi(X)
X <- rpoistlpp(.2,a=0,b=5,L=easynet) densityVoronoi(X)
This function performs adaptive intensity estimation for temporal point patterns using Voronoi-Dirichlet tessellation.
## S3 method for class 'tpp' densityVoronoi(X, f = 1, nrep = 1, at=c("points","pixels"), dimt=128,...)
## S3 method for class 'tpp' densityVoronoi(X, f = 1, nrep = 1, at=c("points","pixels"), dimt=128,...)
X |
an object of class |
f |
fraction (between 0 and 1 inclusive) of the data points that will be used to build a tessellation for the intensity estimate |
nrep |
number of independent repetitions of the randomised procedure |
at |
string specifying whether to compute the intensity values at a grid of pixel locations and time (at="pixels") or only at the points of x (at="points"). default is to estimate the intensity at pixels |
dimt |
the number of equally spaced points at which the temporal density is to be estimated. see density |
... |
arguments passed to |
This function computes intensity estimates for temporal point patterns using Voronoi-Dirichlet tessellation.
IF f<1, then nrep independent sub-samples of X are obtained using the function rthin.stlpp
. Then for each of the obtained sub-samples, we calculate the Voronoi estimate. The final estimation is the sum of all obtained estimated intensities divided by (f*nrep).
If at="points"
: a vector of intensity values at the data points of X.
If at="pixels"
: a vector of intensity values over a grid.
Mehdi Moradi <[email protected]> and Ottmar Cronie
Mateu, J., Moradi, M., & Cronie, O. (2019). Spatio-temporal point patterns on linear networks: Pseudo-separable intensity estimation. Spatial Statistics, 100400.
densityVoronoi.lpp
, density.stlpp
X <- rpoistlpp(0.2,a=0,b=5,L=easynet) Y <- as.tpp.stlpp(X) densityVoronoi(Y)
X <- rpoistlpp(0.2,a=0,b=5,L=easynet) Y <- as.tpp.stlpp(X) densityVoronoi(Y)
This dataset represents the spatio-temporal locations of traffic accidents in the down-town of Eastbourne (UK) in the period of 2005-2010. The network was provided by “OS OpenData" at www.ordnancesurvey.co.uk and is usable under the terms of the OS OpenData license. The traffic locations were collected by the UK Department for Transport at www.data.gov.uk and obtained through kaggle at www.kaggle.com.
The dataset Eastbourne
is an object of class stlpp
.
data(Eastbourne)
data(Eastbourne)
Mehdi Moradi <[email protected]>
Usability: The network of Eastbourne was provided by OS OpenData and contains OS data © Crown copyright and database right (2018). The traffic accident locations in Eastbourne were collected by the UK Department for Transport and were provided by kaggle.
This data is a part of enitre data which is selected and converted to this format by Mehdi Moradi.
Moradi, M., & Mateu, J. (2020). First-and second-order characteristics of spatio-temporal point processes on linear networks. Journal of Computational and Graphical Statistics, 29(3), 432-443.
data(Eastbourne) plot(Eastbourne)
data(Eastbourne) plot(Eastbourne)
A simple and not real network.
data(easynet)
data(easynet)
Created by Mehdi Moradi
This dataset represents the spatio-temporal locations of traffic accidents in an area near the pontifical bolivarian university in Medellin (Colombia) during 2016. The entire data were published in the OpenData portal of Medellin Town Hall at https://www.medellin.gov.co/geomedellin/index.hyg.
The dataset Medellin
is an object of class stlpp
.
data(Medellin)
data(Medellin)
Mehdi Moradi <[email protected]>
This data is a part of enitre data which is selected and converted to this format by Mehdi Moradi.
Moradi, M., & Mateu, J. (2020). First-and second-order characteristics of spatio-temporal point processes on linear networks. Journal of Computational and Graphical Statistics, 29(3), 432-443.
data(Medellin) plot(Medellin)
data(Medellin) plot(Medellin)
Methods for space-time point patterns on a linear network.
## S3 method for class 'stlpp' plot(x,xlab = xlab,...) ## S3 method for class 'stlppint' plot(x,style=style,xlab=xlab,xlim=xlim,ylim=ylim,bar=TRUE,...) ## S3 method for class 'sumstlpp' plot(x,style=c("level","contour","perspective"), theta = 35, phi = 10, facets = FALSE, ticktype = "detailed", resfac = 5,xlab="r = distance",ylab="t = time",...) ## S3 method for class 'stlpp' print(x,...) ## S3 method for class 'stlppint' print(x,...) ## S3 method for class 'sumstlpp' print(x,...) ## S3 method for class 'stlpp' x[i] ## S3 method for class 'stlppint' x[i] ## S3 method for class 'stlppint' as.linim(X,...) ## S3 method for class 'stlppint' as.tppint(x) ## S3 method for class 'sumstlpp' as.data.frame(x,...)
## S3 method for class 'stlpp' plot(x,xlab = xlab,...) ## S3 method for class 'stlppint' plot(x,style=style,xlab=xlab,xlim=xlim,ylim=ylim,bar=TRUE,...) ## S3 method for class 'sumstlpp' plot(x,style=c("level","contour","perspective"), theta = 35, phi = 10, facets = FALSE, ticktype = "detailed", resfac = 5,xlab="r = distance",ylab="t = time",...) ## S3 method for class 'stlpp' print(x,...) ## S3 method for class 'stlppint' print(x,...) ## S3 method for class 'sumstlpp' print(x,...) ## S3 method for class 'stlpp' x[i] ## S3 method for class 'stlppint' x[i] ## S3 method for class 'stlppint' as.linim(X,...) ## S3 method for class 'stlppint' as.tppint(x) ## S3 method for class 'sumstlpp' as.data.frame(x,...)
x , X
|
an object of classes |
style |
style of plot |
theta , phi
|
see persp3D |
facets , ticktype
|
see persp3D |
resfac |
see persp3D |
xlab , ylab
|
the x,y label of the plot |
xlim |
giving the x limits for the plot |
ylim |
giving the y limits for the plot |
bar |
if TRUE, bar plot of rounded time occurances will be added to the density plot |
i |
numeric, logical, or an object of class |
... |
either ignore for |
Mehdi Moradi <[email protected]>
Methods for one-dimensional point patterns.
## S3 method for class 'tpp' plot(x,xlab="time",ylab="",main = "cumulative number",...) ## S3 method for class 'tppint' plot(x,xlab=xlab,xlim=xlim,line=2.5,main="NULL",...) ## S3 method for class 'tpp' print(x,...) ## S3 method for class 'tppint' print(x,...) ## S3 method for class 'tpp' x[i] ## S3 method for class 'tppint' x[i]
## S3 method for class 'tpp' plot(x,xlab="time",ylab="",main = "cumulative number",...) ## S3 method for class 'tppint' plot(x,xlab=xlab,xlim=xlim,line=2.5,main="NULL",...) ## S3 method for class 'tpp' print(x,...) ## S3 method for class 'tppint' print(x,...) ## S3 method for class 'tpp' x[i] ## S3 method for class 'tppint' x[i]
x |
an object of class tpp or tppint. |
xlab , ylab
|
the x,y label of the plot. |
main |
overall title for the plot. |
xlim |
giving the x limits for the plot. |
line |
specifying a value for line overrides the default placement of y label, and places it this many lines outwards from the plot edge. |
i |
numeric, logical, or an object of class |
... |
graphics parameters passed to plot/print function. |
Mehdi Moradi <[email protected]>
X <- tpp(sample(c(1:24),200,replace = TRUE)) plot(X) plot(density(X))
X <- tpp(sample(c(1:24),200,replace = TRUE)) plot(X) plot(density(X))
This function simulates realisations of a spatio-temporal Poisson point process on a linear network.
rpoistlpp(lambda,a,b,L,check=FALSE,lmax=NULL,nsim=1)
rpoistlpp(lambda,a,b,L,check=FALSE,lmax=NULL,nsim=1)
lambda |
intensity of the point process. it can be either a number, function of location and time, or an abject of class |
a |
lower bound of time period |
b |
upper bound of time period |
L |
a linear network |
check |
logical value indicating whether to check that all the (x,y) points lie inside the specified window. see |
lmax |
upper bound for the values of |
nsim |
number of simulated patterns to generate |
This function generates realisations of a spatio-temporal poisson point process on a linear network based on an intensity function lambda and lower/upper bounds a and b.
an object of class stlpp
if nsim=1, otherwise a list of objects of class stlpp
.
Mehdi Moradi <[email protected]>
Moradi, M., & Mateu, J. (2020). First-and second-order characteristics of spatio-temporal point processes on linear networks. Journal of Computational and Graphical Statistics, 29(3), 432-443.
X <- rpoistlpp(0.2,a=0,b=5,L=easynet) X
X <- rpoistlpp(0.2,a=0,b=5,L=easynet) X
This function simulates realisations of an one-dimensional Poisson point process.
rpoistpp(lambda,a,b,check=FALSE,lmax=NULL,nsim=1)
rpoistpp(lambda,a,b,check=FALSE,lmax=NULL,nsim=1)
lambda |
intensity of the point process. it can be either a number, a function of location and time, or an object of class |
a |
lower bound of time period |
b |
upper bound of time period |
check |
logical value indicating whether to check that all the (x,y) points lie inside the specified time period. |
lmax |
upper bound for the values of |
nsim |
number of simulated patterns to generate |
This function generates realisations of a temporal poisson point process based on a given intensity function lambda and lower/upper bounds a and b.
an object of class tpp
if nsim=1, otherwise a list of objects of class tpp
.
Mehdi Moradi <[email protected]>
Moradi, M., & Mateu, J. (2020). First-and second-order characteristics of spatio-temporal point processes on linear networks. Journal of Computational and Graphical Statistics, 29(3), 432-443.
f <- function(t){0.1*exp(t)} X <- rpoistpp(f,a=1,b=10)
f <- function(t){0.1*exp(t)} X <- rpoistpp(f,a=1,b=10)
This function applies independent random thinning to a spatio-temporal point pattern on a linear network.
## S3 method for class 'stlpp' rthin(X, P = P, nsim = 1)
## S3 method for class 'stlpp' rthin(X, P = P, nsim = 1)
X |
a spatio-temporal point pattern of class |
P |
retention probability |
nsim |
number of simulated realisations to be generated |
See rthin
.
An object of the same kind as X if nsim=1, or a list of such objects if nsim > 1.
Mehdi Moradi <[email protected]>
Moradi, M., & Mateu, J. (2020). First-and second-order characteristics of spatio-temporal point processes on linear networks. Journal of Computational and Graphical Statistics, 29(3), 432-443.
data(Medellin) rthin(Medellin,P=.5)
data(Medellin) rthin(Medellin,P=.5)
This function computes the pair correlation function for spatio-temporal point patterns on linear networks.
STLg(X,r=NULL,t=NULL,nxy=10)
STLg(X,r=NULL,t=NULL,nxy=10)
X |
a spatio-temporal point pattern of class |
r |
values of argument r where pair correlation function will be evaluated. optional |
t |
values of argument t where pair correlation function will be evaluated. optional |
nxy |
pixel array dimensions. optional |
This function calculates the pair correlation function for a homogeneous spatio-temporal point patterns on a linear network.
An object of class sumstlpp
.
Mehdi Moradi <[email protected]>
Moradi, M., & Mateu, J. (2020). First-and second-order characteristics of spatio-temporal point processes on linear networks. Journal of Computational and Graphical Statistics, 29(3), 432-443.
X <- rpoistlpp(.2,a=0,b=5,L=easynet) g <- STLg(X) plot(g)
X <- rpoistlpp(.2,a=0,b=5,L=easynet) g <- STLg(X) plot(g)
This function computes the inhomogeneous pair correlation function for spatio-temporal point patterns on linear networks.
STLginhom(X,lambda,normalize=FALSE,r=NULL,t=NULL,nxy=10)
STLginhom(X,lambda,normalize=FALSE,r=NULL,t=NULL,nxy=10)
X |
a spatio-temporal point pattern of class |
lambda |
values of estimated intensity at data points |
normalize |
normalization factor to be considered |
r |
values of argument r where pair correlation function will be evaluated. optional |
t |
values of argument t where pair correlation function will be evaluated. optional |
nxy |
pixel array dimensions. optional |
This function calculates the inhomogeneous pair correlation function for a spatio-temporal point patterns on a linear network.
An object of class sumstlpp
.
Mehdi Moradi <[email protected]>
Moradi, M., & Mateu, J. (2020). First-and second-order characteristics of spatio-temporal point processes on linear networks. Journal of Computational and Graphical Statistics, 29(3), 432-443.
X <- rpoistlpp(.2,a=0,b=5,L=easynet) d <- density(X,at="points") g <- STLginhom(X,lambda=d,normalize=TRUE) plot(g)
X <- rpoistlpp(.2,a=0,b=5,L=easynet) d <- density(X,at="points") g <- STLginhom(X,lambda=d,normalize=TRUE) plot(g)
This function computes the K-function for spatio-temporal point patterns on linear networks.
STLK(X,r=NULL,t=NULL,nxy=10)
STLK(X,r=NULL,t=NULL,nxy=10)
X |
a spatio-temporal point pattern of class |
r |
values of argument r where pair correlation function will be evaluated. optional |
t |
values of argument t where pair correlation function will be evaluated. optional |
nxy |
pixel array dimensions. optional |
This function calculates the K-function for a homogeneous spatio-temporal point patterns on a linear network.
An object of class sumstlpp
.
Mehdi Moradi <[email protected]>
Moradi, M., & Mateu, J. (2020). First-and second-order characteristics of spatio-temporal point processes on linear networks. Journal of Computational and Graphical Statistics, 29(3), 432-443.
X <- rpoistlpp(.2,a=0,b=5,L=easynet) k <- STLK(X) plot(k)
X <- rpoistlpp(.2,a=0,b=5,L=easynet) k <- STLK(X) plot(k)
This function computes the inhomogeneous K-function for spatio-temporal point patterns on linear networks.
STLKinhom(X,lambda=lambda,normalize=FALSE,r=NULL,t=NULL,nxy=10)
STLKinhom(X,lambda=lambda,normalize=FALSE,r=NULL,t=NULL,nxy=10)
X |
a spatio-temporal point pattern of class |
lambda |
values of estimated intensity at data points |
normalize |
normalization factor to be considered |
r |
values of argument r where pair correlation function will be evaluated. optional |
t |
values of argument t where pair correlation function will be evaluated. optional |
nxy |
pixel array dimensions. optional |
This function calculates the inhomogeneous K-function for a spatio-temporal point patterns on a linear network.
An object of class sumstlpp
.
Mehdi Moradi <[email protected]>
Moradi, M., & Mateu, J. (2020). First-and second-order characteristics of spatio-temporal point processes on linear networks. Journal of Computational and Graphical Statistics, 29(3), 432-443.
X <- rpoistlpp(.2,a=0,b=5,L=easynet) lambda <- density(X,at="points") k <- STLKinhom(X,lambda=lambda,normalize=TRUE) plot(k)
X <- rpoistlpp(.2,a=0,b=5,L=easynet) lambda <- density(X,at="points") k <- STLKinhom(X,lambda=lambda,normalize=TRUE) plot(k)
Create an object of class stlpp
representing a spatio-temporal point pattern on a linear network.
stlpp(X, L, T, ...)
stlpp(X, L, T, ...)
X |
Locations of the points. a matrix or data frame of coordinates, or a point pattern object (of class "ppp") or other data acceptable to |
L |
linear network (object of class |
T |
time occurrence of the points |
... |
ignored |
This function creates an object of class stlpp
. For details about X see lpp
. T
represents the time occurrences of data points.
An object of class stlpp
.
Mehdi Moradi <[email protected]>
data(easynet) X <- rpoislpp(1,easynet) t <- runif(npoints(X)) stlpp(X,T=t,L=easynet)
data(easynet) X <- rpoislpp(1,easynet) t <- runif(npoints(X)) stlpp(X,T=t,L=easynet)
Create an object of class tpp
representing a one-dimensional point pattern.
tpp(X,a,b)
tpp(X,a,b)
X |
|
a |
lower band of the time domain. if not given by the user, it will be the minimum of X |
b |
upper bound of the time domain. if not given by the user, it will be the maximum of X |
Create a one-dimensional point pattern.
An object of class tpp.
Mehdi Moradi <[email protected]>
tpp(runif(10))
tpp(runif(10))
This function extracts unique points from a spatio-temporal point pattern on a linear network.
## S3 method for class 'stlpp' unique(x,...)
## S3 method for class 'stlpp' unique(x,...)
x |
a spatio-temporal point pattern of class |
... |
arguments for |
This function extracts unique points from a spatio-temporal point pattern on a linear network.
A spatio-temporal point pattern on a linear network with no duplicated point.
Mehdi Moradi <[email protected]>
Moradi, M., & Mateu, J. (2020). First-and second-order characteristics of spatio-temporal point processes on linear networks. Journal of Computational and Graphical Statistics, 29(3), 432-443.
X <- rpoistlpp(0.1,0,5,L=easynet) df <- as.data.frame(X) df_dup <- df[sample(nrow(df), 20,replace = TRUE), ] Y <- as.stlpp(df_dup,L=easynet) npoints(Y) npoints(unique(Y))
X <- rpoistlpp(0.1,0,5,L=easynet) df <- as.data.frame(X) df_dup <- df[sample(nrow(df), 20,replace = TRUE), ] Y <- as.stlpp(df_dup,L=easynet) npoints(Y) npoints(unique(Y))