Title: | Generalization of the Convex Hull of a Sample of Points in the Plane |
---|---|
Description: | Computation of the alpha-shape and alpha-convex hull of a given sample of points in the plane. The concepts of alpha-shape and alpha-convex hull generalize the definition of the convex hull of a finite set of points. The programming is based on the duality between the Voronoi diagram and Delaunay triangulation. The package also includes a function that returns the Delaunay mesh of a given sample of points and its dual Voronoi diagram in one single object. |
Authors: | Beatriz Pateiro-Lopez [aut, cre], Alberto Rodriguez-Casal, [aut]. |
Maintainer: | Beatriz Pateiro-Lopez <[email protected]> |
License: | GPL-2 | GPL-3 |
Version: | 2.5 |
Built: | 2024-11-02 06:14:14 UTC |
Source: | https://github.com/beatrizpateiro/alphahull |
Computation of the -shape and
-convex hull of a given sample of points in the plane. The concepts of
-shape and
-convex hull generalize the definition of the convex hull of a finite set of points. The programming is based on the duality between the Voronoi diagram and Delaunay triangulation. The package also includes a function that returns the Delaunay mesh of a given sample of points and its dual Voronoi diagram in one single object.
Package: | alphahull |
Type: | Package |
Version: | 2.5 |
Date: | 2022-06-16 |
License: | R functions: GPL-2 | GPL-3 |
Beatriz Pateiro-Lopez, Alberto Rodriguez-Casal.
Maintainer: Beatriz Pateiro-Lopez <[email protected]>
This function calculates the -convex hull of a given sample of points in the plane for
.
ahull(x, y = NULL, alpha)
ahull(x, y = NULL, alpha)
x , y
|
The |
alpha |
Value of |
An attempt is made to interpret the arguments x and y in a way suitable for computing the -convex hull. Any reasonable way of defining the coordinates is acceptable, see
xy.coords
.
The -convex hull is defined for any finite number of points. However, since the algorithm is based on the Delaunay triangulation, at least three non-collinear points are required.
If y
is NULL and x
is an object of class "delvor"
, then the -convex hull is computed with no need to invoke again the function
delvor
(it reduces the computational cost).
The complement of the -convex hull can be written as the union of
open balls and halfplanes, see
complement
.
The boundary of the -convex hull is formed by arcs of open balls of radius
(besides possible isolated sample points). The arcs are determined by the intersections of some of the balls that define the complement of the
-convex hull. The extremes of an arc are given by
and
where
and
represent the center and radius of the arc, repectively, and
represents the clockwise rotation of angle
of the unitary vector
. Joining the end points of adjacent arcs we can define polygons that help us to determine the area of the estimator , see
areaahull
.
A list with the following components:
arcs |
For each arc in the boundary of the |
xahull |
A 2-column matrix with the coordinates of the original set of points besides possible new end points of the arcs in the boundary of the |
length |
Length of the boundary of the |
complement |
Output matrix from |
alpha |
Value of |
ashape.obj |
Object of class |
Edelsbrunner, H., Kirkpatrick, D.G. and Seidel, R. (1983). On the shape of a set of points in the plane. IEEE Transactions on Information Theory, 29(4), pp.551-559.
Rodriguez-Casal, R. (2007). Set estimation under convexity type assumptions. Annales de l'I.H.P.- Probabilites & Statistiques, 43, pp.763-774.
Pateiro-Lopez, B. (2008). Set estimation under convexity type restrictions. Phd. Thesis. Universidad de Santiago de Compostela. ISBN 978-84-9887-084-8.
## Not run: # Random sample in the unit square x <- matrix(runif(100), nc = 2) # Value of alpha alpha <- 0.2 # Alpha-convex hull ahull.obj <- ahull(x, alpha = alpha) plot(ahull.obj) # Uniform sample of size n=300 in the annulus B(c,0.5)\B(c,0.25), # with c=(0.5,0.5). n <- 300 theta<-runif(n,0,2*pi) r<-sqrt(runif(n,0.25^2,0.5^2)) x<-cbind(0.5+r*cos(theta),0.5+r*sin(theta)) # Value of alpha alpha <- 0.1 # Alpha-convex hull ahull.obj <- ahull(x, alpha = alpha) # The arcs defining the boundary of the alpha-convex hull are ordered plot(x) for (i in 1:dim(ahull.obj$arcs)[1]){ arc(ahull.obj$arcs[i,1:2],ahull.obj$arcs[i,3],ahull.obj$arcs[i,4:5], ahull.obj$arcs[i,6],col=2) Sys.sleep(0.5) } # Random sample from a uniform distribution on a Koch snowflake # with initial side length 1 and 3 iterations x <- rkoch(2000, side = 1, niter = 3) # Value of alpha alpha <- 0.05 # Alpha-convex hull ahull.obj <- ahull(x, alpha = alpha) plot(ahull.obj) ## End(Not run)
## Not run: # Random sample in the unit square x <- matrix(runif(100), nc = 2) # Value of alpha alpha <- 0.2 # Alpha-convex hull ahull.obj <- ahull(x, alpha = alpha) plot(ahull.obj) # Uniform sample of size n=300 in the annulus B(c,0.5)\B(c,0.25), # with c=(0.5,0.5). n <- 300 theta<-runif(n,0,2*pi) r<-sqrt(runif(n,0.25^2,0.5^2)) x<-cbind(0.5+r*cos(theta),0.5+r*sin(theta)) # Value of alpha alpha <- 0.1 # Alpha-convex hull ahull.obj <- ahull(x, alpha = alpha) # The arcs defining the boundary of the alpha-convex hull are ordered plot(x) for (i in 1:dim(ahull.obj$arcs)[1]){ arc(ahull.obj$arcs[i,1:2],ahull.obj$arcs[i,3],ahull.obj$arcs[i,4:5], ahull.obj$arcs[i,6],col=2) Sys.sleep(0.5) } # Random sample from a uniform distribution on a Koch snowflake # with initial side length 1 and 3 iterations x <- rkoch(2000, side = 1, niter = 3) # Value of alpha alpha <- 0.05 # Alpha-convex hull ahull.obj <- ahull(x, alpha = alpha) plot(ahull.obj) ## End(Not run)
This function approximates the -convex hull of tracking data and returns a list of geom_path objects of the boundary.
ahull_track(x, y = NULL, alpha, nps = 10000, sc = 100)
ahull_track(x, y = NULL, alpha, nps = 10000, sc = 100)
x , y
|
The |
alpha |
Value of |
nps |
Number of points to generate in each segment connecting two locations, see Details |
sc |
Scale factor. |
An attempt is made to interpret the arguments x and y in a way suitable for computing the -convex hull. Any reasonable way of defining the coordinates is acceptable, see
xy.coords
.
Increase nps
if the trajectory is not contained in the computed estimator.
A list of geom_path objects defining the boundary of the -convex
Cholaquidis, A., Fraiman, R., Lugosi, G. and Pateiro-Lopez, B. (2014) Set estimation from reflected Brownian motion. arXiv:1411.0433.
Wikelski, M., and Kays, R. (2014). Movebank: archive, analysis and sharing of animal movement data. World Wide Web electronic publication.
## Not run: library(move) library(ggmap) # Data from Movebank # Study Name: Dunn Ranch Bison Tracking Project # Principal Investigator: Stephen Blake, Randy Arndt, Doug Ladd # Max Planck Institute for Ornithology Radolfzell Germany study <- "Dunn Ranch Bison Tracking Project" cainfo <- system.file("CurlSSL", "cacert.pem", package = "RCurl") options(RCurlOptions = list(verbose = FALSE, capath = cainfo, ssl.verifypeer = FALSE)) # Login to movebank (first create the login object) curl <- movebankLogin(username = "xxx", password = "zzz") # Downloads study stored in Movebank track <- getMovebankData(study = study, login = curl) dat <- track@data[track@data[, "deployment_id"] == 13848432,] # Map of animal locations bbox <- ggmap::make_bbox(dat[,"location_long"], dat[,"location_lat"], f = 0.3) map_loc <- get_map(location = bbox, source = "google", maptype = 'satellite') map <- ggmap(map_loc, extent = 'panel', maprange=FALSE) p <- map + geom_path(data = dat, aes(x = location_long, y = location_lat), col=2, size=0.3) p ah_gp <- ahull_track(x = dat[, c("location_long", "location_lat")], alpha = 0.005) p + ah_gp ## End(Not run)
## Not run: library(move) library(ggmap) # Data from Movebank # Study Name: Dunn Ranch Bison Tracking Project # Principal Investigator: Stephen Blake, Randy Arndt, Doug Ladd # Max Planck Institute for Ornithology Radolfzell Germany study <- "Dunn Ranch Bison Tracking Project" cainfo <- system.file("CurlSSL", "cacert.pem", package = "RCurl") options(RCurlOptions = list(verbose = FALSE, capath = cainfo, ssl.verifypeer = FALSE)) # Login to movebank (first create the login object) curl <- movebankLogin(username = "xxx", password = "zzz") # Downloads study stored in Movebank track <- getMovebankData(study = study, login = curl) dat <- track@data[track@data[, "deployment_id"] == 13848432,] # Map of animal locations bbox <- ggmap::make_bbox(dat[,"location_long"], dat[,"location_lat"], f = 0.3) map_loc <- get_map(location = bbox, source = "google", maptype = 'satellite') map <- ggmap(map_loc, extent = 'panel', maprange=FALSE) p <- map + geom_path(data = dat, aes(x = location_long, y = location_lat), col=2, size=0.3) p ah_gp <- ahull_track(x = dat[, c("location_long", "location_lat")], alpha = 0.005) p + ah_gp ## End(Not run)
Given a vector and an angle
,
anglesArc
returns the angles that and
form with the axis OX, where
represents the clockwise rotation of angle
of the vector
.
anglesArc(v, theta)
anglesArc(v, theta)
v |
Vector |
theta |
Angle |
The angle that forms the vector with the axis OX takes its value in
.
angs |
Numeric vector with two components. |
## Not run: # Let v=c(0,1) and theta=pi/4 # Consider the arc such that v is the internal angle bisector that # divides the angle 2*theta into two equal angles # The angles that the arc forms with the OX axis are pi/4 and 3*pi/4 v <- c(0,1) theta <- pi/4 anglesArc(v,theta) ## End(Not run)
## Not run: # Let v=c(0,1) and theta=pi/4 # Consider the arc such that v is the internal angle bisector that # divides the angle 2*theta into two equal angles # The angles that the arc forms with the OX axis are pi/4 and 3*pi/4 v <- c(0,1) theta <- pi/4 anglesArc(v,theta) ## End(Not run)
This function adds the arc of between the angles that
and
form with the axis OX, where
represents the clockwise rotation of angle
of the vector
.
arc(c, r, v, theta, ...)
arc(c, r, v, theta, ...)
c |
Center |
r |
Radius |
v |
Vector |
theta |
Angle |
... |
Arguments to be passed to methods, such as graphical parameters (see |
## Not run: # Plot of the circumference of radius 1 theta <- seq(0, 2*pi, length = 100) r <- 1 plot(r*cos(theta), r*sin(theta), type = "l") # Add in red the arc between pi/4 and 3*pi/4 arc(c(0,0), 1, c(0,1), pi/4, col = 2, lwd = 2) ## End(Not run)
## Not run: # Plot of the circumference of radius 1 theta <- seq(0, 2*pi, length = 100) r <- 1 plot(r*cos(theta), r*sin(theta), type = "l") # Add in red the arc between pi/4 and 3*pi/4 arc(c(0,0), 1, c(0,1), pi/4, col = 2, lwd = 2) ## End(Not run)
This function calculates the area of the -convex hull of a sample of points.
areaahull(x, timeout = 5)
areaahull(x, timeout = 5)
x |
Object of class |
timeout |
A numeric specifying the maximum number of seconds the expression is allowed to run before being interrupted by the timeout. |
area |
Area of the |
## Not run: # Random sample in the unit square x <- matrix(runif(500), nc = 2) # Value of alpha alpha <- 1 # alpha-convex hull ahull.obj <- ahull(x, alpha = alpha) # Area of the alpha-convex hull areaahull(ahull.obj) ## End(Not run)
## Not run: # Random sample in the unit square x <- matrix(runif(500), nc = 2) # Value of alpha alpha <- 1 # alpha-convex hull ahull.obj <- ahull(x, alpha = alpha) # Area of the alpha-convex hull areaahull(ahull.obj) ## End(Not run)
This function calculates the -shape of a given sample for
.
ashape(x, y = NULL, alpha)
ashape(x, y = NULL, alpha)
x , y
|
The |
alpha |
Value of |
An attempt is made to interpret the arguments x and y in a way suitable for computing the -shape, see
xy.coords
.
The -shape is defined for any finite number of points. However, since the algorithm is based on the Delaunay triangulation, at least three non-collinear points are required.
If y
is NULL and x
is an object of class "delvor"
, then the -shape is computed without invoking again the function
delvor
(it reduces the computational cost).
The function ashape
returns (among other values) the matrix edges
. The structure of edges
is that of matrix mesh
returned by the function delvor
. Note that the -shape is a subgraph of the Delaunay triangulation and, therefore,
edges
is a submatrix of mesh
.
A list with the following components:
edges |
A n.seg-row matrix with the coordinates and indexes of the edges of the Delaunay triangulation that form the |
length |
Length of the |
alpha |
Value of |
alpha.extremes |
Vector with the indexes of the sample points that are |
delvor.obj |
Object of class |
x |
A 2-column matrix with the coordinates of the set of points. |
Edelsbrunner, H., Kirkpatrick, D.G. and Seidel, R. (1983). On the shape of a set of points in the plane. IEEE Transactions on Information Theory, 29(4), pp.551-559.
## Not run: # Uniform sample of size n=300 in the annulus B(c,0.5)\B(c,0.25), # with c=(0.5,0.5). n <- 300 theta<-runif(n,0,2*pi) r<-sqrt(runif(n,0.25^2,0.5^2)) x<-cbind(0.5+r*cos(theta),0.5+r*sin(theta)) # Value of alpha alpha <- 0.1 # alpha-shape ashape.obj <- ashape(x, alpha = alpha) # If we change the value of alpha there is no need to compute # again the Delaunay triangulation and Voronoi Diagram alpha <- 0.4 ashape.obj.new <- ashape(ashape.obj$delvor.obj, alpha = alpha) # Random sample from a uniform distribution on a Koch snowflake # with initial side length 1 and 3 iterations x <- rkoch(2000, side = 1, niter = 3) # Value of alpha alpha <- 0.05 # alpha-shape ashape.obj <- ashape(x, alpha = alpha) ## End(Not run)
## Not run: # Uniform sample of size n=300 in the annulus B(c,0.5)\B(c,0.25), # with c=(0.5,0.5). n <- 300 theta<-runif(n,0,2*pi) r<-sqrt(runif(n,0.25^2,0.5^2)) x<-cbind(0.5+r*cos(theta),0.5+r*sin(theta)) # Value of alpha alpha <- 0.1 # alpha-shape ashape.obj <- ashape(x, alpha = alpha) # If we change the value of alpha there is no need to compute # again the Delaunay triangulation and Voronoi Diagram alpha <- 0.4 ashape.obj.new <- ashape(ashape.obj$delvor.obj, alpha = alpha) # Random sample from a uniform distribution on a Koch snowflake # with initial side length 1 and 3 iterations x <- rkoch(2000, side = 1, niter = 3) # Value of alpha alpha <- 0.05 # alpha-shape ashape.obj <- ashape(x, alpha = alpha) ## End(Not run)
This function calculates the complement of the -convex hull of a given sample for
.
complement(x, y = NULL, alpha)
complement(x, y = NULL, alpha)
x , y
|
The |
alpha |
Value of |
An attempt is made to interpret the arguments x and y in a way suitable for computing the -shape. Any reasonable way of defining the coordinates is acceptable, see
xy.coords
.
If y
is NULL and x
is an object of class "delvor"
, then the complement of the -convex hull is computed with no need to invoke again the function
delvor
(it reduces the computational cost).
The complement of the -convex hull is calculated as a union of open balls and halfplanes that do not contain any point of the sample. See Edelsbrunnner et al. (1983) for a basic description of the algorithm. The construction of the complement is based on the Delaunay triangulation and Voronoi diagram of the sample, provided by the function
delvor
. The function complement
returns a matrix compl
. For each row i
, compl[i,]
contains the information relative to an open ball or halfplane of the complement. The first three columns are assigned to the characterization of the ball or halfplane i
. The information relative to the edge of the Delaunay triangulation that generates the ball or halfplane i
is contained in compl[i,4:16]
. Thus, if the row i
refers to an open ball, compl[i,1:3]
contains the center and radius of the ball. Furthermore, compl[i,17:18]
and compl[i,19]
refer to the unitary vector and the angle
that characterize the arc that joins the two sample points that define the ball
i
. If the row i
refers to a halfplane, compl[i,1:3]
determines its equation. For the halfplane ,
compl[i,1:3]=(a,b,-1)
. In the same way, for the halfplane ,
compl[i,1:3]=(a,b,-2)
, for the halfplane ,
compl[i,1:3]=(a,0,-3)
and for the halfplane ,
compl[i,1:3]=(a,0,-4)
.
compl |
Output matrix. For each row |
Edelsbrunner, H., Kirkpatrick, D.G. and Seidel, R. (1983) On the shape of a set of points in the plane. IEEE Transactions on Information Theory, 29(4), pp.551-559.
## Not run: # Random sample in the unit square x <- matrix(runif(100), nc = 2) # Value of alpha alpha <- 0.2 # Complement of the alpha-convex hull compl <- complement(x, alpha = alpha) ## End(Not run)
## Not run: # Random sample in the unit square x <- matrix(runif(100), nc = 2) # Value of alpha alpha <- 0.2 # Complement of the alpha-convex hull compl <- complement(x, alpha = alpha) ## End(Not run)
This function returns a matrix with information about the Delaunay triangulation and Voronoi diagram of a given sample.
delvor(x, y = NULL)
delvor(x, y = NULL)
x , y
|
The |
An attempt is made to interpret the arguments x and y in a way suitable for computing the Delaunay triangulation and Voronoi diagram . Any reasonable way of defining the coordinates is acceptable, see xy.coords
.
The function tri.mesh
from package interp calculates the Delaunay triangulation of at least three non-collinear points. Using the Delaunay triangulation, the function delvor
calculates the correspondig Voronoi diagram. For each edge of the Delaunay triangulation there is a segment in the Voronoi diagram, given by the union of the circumcenters of the two neighbour triangles that share the edge. For those triangles with edges on the convex hull, the corresponding line in the Voronoi diagram is a semi-infinite segment, whose boundless extreme is calculated by the function dummycoor
. The function delvor
returns the sample, the output object of class "triSht"
from the function tri.mesh
and a matrix mesh
with all the necessary information of the Delaunay triangulation and Voronoi diagram. Thus, for each edge of the Delaunay triangulation the output matrix contains the indexes and coordinates of the sample points that form the edge, the indexes and coordinates of the extremes of the corresponding segment in the Voronoi diagram, and an indicator that takes the value 1 for those extremes of the Voronoi diagram that represent a boundless extreme.
A list with the following components:
mesh |
A |
x |
A 2-column matrix with the coordinates of the sample points. |
tri.obj |
Object of class |
Renka, R. J. (1996). Algorithm 751: TRIPACK: a constrained two-dimensional Delaunay triangulation package, ACM Trans. Math. Softw., 22(1), pp.1-8.
## Not run: # Random sample in the unit square x <- matrix(runif(20), nc = 2) # Delaunay triangulation and Voronoi diagram calculation delvor.obj <- delvor(x) ## End(Not run)
## Not run: # Random sample in the unit square x <- matrix(runif(20), nc = 2) # Delaunay triangulation and Voronoi diagram calculation delvor.obj <- delvor(x) ## End(Not run)
This function determines fictitious coordinates for the boundless extreme of a semi-infinite edge of the Voronoi diagram.
dummycoor(tri.obj, l1, l2, m, away)
dummycoor(tri.obj, l1, l2, m, away)
tri.obj |
Object of class |
l1 |
Index of the sample point correponding to one vertex of a triangle of Delaunay that lies on the convex hull, see Details. |
l2 |
Index of the sample point correponding to other vertex of a triangle of Delaunay that lies on the convex hull, see Details. |
m |
Index of the circumcenter of the triangle of Delaunay with one edge on the convex hull. |
away |
Constant that determines how far away the fictitious boundless extreme is located. |
When a triangle of the Delaunay triangulation has one of its edges (given by the segment that joins the sample points with indexes l1
and l2
) on the convex hull, the corresponding segment of the Voronoi diagram is semi-infinite. The finite extreme coincides with the circumcenter of the triangle and the direction of the line is given by the perpendicular bisector of the edge that lies on the convex hull.
dum |
Fictitious coordinates of the boundless extreme. |
This function calculates the Devroye-Wise estimator of a given sample of points in the plane for .
dw(x, y = NULL, eps)
dw(x, y = NULL, eps)
x , y
|
The |
eps |
Value of |
An attempt is made to interpret the arguments x and y in a way suitable for computing the Devroye-Wise estimator. Any reasonable way of defining the coordinates is acceptable, see xy.coords
.
Given a sample of points in the plane, the estimator is defined as union of balls of radius with centers in the sample points. For each arc in the boundary of the Devroye-Wise estimator, the columns of the output matrix store the center
and radius
of the arc, the unitary vector
, the angle
that define the arc and the indices of the end points.
Devroye, L. and Wise, G. (1980) Detection of abnormal behaviour via nonparametric estimation of the support. SIAM J. Appl. Math. 3, pp. 480-488.
## Not run: # Uniform sample of size n = 200 in the annulus B(c, 0.5)\B(c, 0.25), # with c = (0.5, 0.5). n <- 200 theta <- runif(n, 0, 2*pi) r <- sqrt(runif(n, 0.25^2, 0.5^2)) x <- cbind(0.5 + r*cos(theta), 0.5 + r*sin(theta)) eps <- 0.05 dw.obj <- dw(x, eps = eps) plot(x) for(i in 1:dim(dw.obj)[1]){arc(dw.obj[i, 1:2], eps, dw.obj[i, 4:5], dw.obj[i, 6])} ## End(Not run)
## Not run: # Uniform sample of size n = 200 in the annulus B(c, 0.5)\B(c, 0.25), # with c = (0.5, 0.5). n <- 200 theta <- runif(n, 0, 2*pi) r <- sqrt(runif(n, 0.25^2, 0.5^2)) x <- cbind(0.5 + r*cos(theta), 0.5 + r*sin(theta)) eps <- 0.05 dw.obj <- dw(x, eps = eps) plot(x) for(i in 1:dim(dw.obj)[1]){arc(dw.obj[i, 1:2], eps, dw.obj[i, 4:5], dw.obj[i, 6])} ## End(Not run)
This function approximates the RBM-sausage of tracking data and returns a list of geom_path objects of the boundary.
dw_track(x, y = NULL, eps, nps = 20000, sc = 100)
dw_track(x, y = NULL, eps, nps = 20000, sc = 100)
x , y
|
The |
eps |
Value of |
nps |
Number of points to generate in each segment connecting two locations, see Details. |
sc |
Scale factor. |
An attempt is made to interpret the arguments x and y in a way suitable for computing the RBM-sausage. Any reasonable way of defining the coordinates is acceptable, see xy.coords
.
Given a trajectory in the plane, the estimator is defined as the set of points whose distance to the trajectory is less than or equal to (this estimator is analogous to the one of Devroye and Wise (1980) for i.i.d. data). Increase
nps
if the trajectory is not contained in the computed estimator.
A list of geom_path objects defining the boundary of the estimator
Cholaquidis, A., Fraiman, R., Lugosi, G. and Pateiro-Lopez, B. (2014) Set estimation from reflected Brownian motion. arXiv:1411.0433.
Devroye, L. and Wise, G. (1980) Detection of abnormal behaviour via nonparametric estimation of the support. SIAM J. Appl. Math. 3, pp. 480-488.
Wikelski, M., and Kays, R. (2014). Movebank: archive, analysis and sharing of animal movement data. World Wide Web electronic publication.
## Not run: library(move) library(ggmap) # Data from Movebank # Study Name: Dunn Ranch Bison Tracking Project # Principal Investigator: Stephen Blake, Randy Arndt, Doug Ladd # Max Planck Institute for Ornithology Radolfzell Germany study <- "Dunn Ranch Bison Tracking Project" cainfo <- system.file("CurlSSL", "cacert.pem", package = "RCurl") options(RCurlOptions = list(verbose = FALSE, capath = cainfo, ssl.verifypeer = FALSE)) # Login to movebank (first create the login object) curl <- movebankLogin(username = "xxx", password = "zzz") # Downloads study stored in Movebank track <- getMovebankData(study = study, login = curl) dat <- track@data[track@data[, "deployment_id"] == 13848432,] # Map of animal locations bbox <- ggmap::make_bbox(dat[,"location_long"], dat[,"location_lat"], f = 0.3) map_loc <- get_map(location = bbox, source = "google", maptype = 'satellite') map <- ggmap(map_loc, extent = 'panel', maprange=FALSE) p <- map + geom_path(data = dat, aes(x = location_long, y = location_lat), col=2, size=0.3) p ah_dw <- dw_track(x = dat[, c("location_long", "location_lat")], eps = 0.001) p + ah_dw ## End(Not run)
## Not run: library(move) library(ggmap) # Data from Movebank # Study Name: Dunn Ranch Bison Tracking Project # Principal Investigator: Stephen Blake, Randy Arndt, Doug Ladd # Max Planck Institute for Ornithology Radolfzell Germany study <- "Dunn Ranch Bison Tracking Project" cainfo <- system.file("CurlSSL", "cacert.pem", package = "RCurl") options(RCurlOptions = list(verbose = FALSE, capath = cainfo, ssl.verifypeer = FALSE)) # Login to movebank (first create the login object) curl <- movebankLogin(username = "xxx", password = "zzz") # Downloads study stored in Movebank track <- getMovebankData(study = study, login = curl) dat <- track@data[track@data[, "deployment_id"] == 13848432,] # Map of animal locations bbox <- ggmap::make_bbox(dat[,"location_long"], dat[,"location_lat"], f = 0.3) map_loc <- get_map(location = bbox, source = "google", maptype = 'satellite') map <- ggmap(map_loc, extent = 'panel', maprange=FALSE) p <- map + geom_path(data = dat, aes(x = location_long, y = location_lat), col=2, size=0.3) p ah_dw <- dw_track(x = dat[, c("location_long", "location_lat")], eps = 0.001) p + ah_dw ## End(Not run)
This function determines for one or more points whether they belong to the
-convex hull of a sample.
inahull(ahull.obj, p)
inahull(ahull.obj, p)
ahull.obj |
Object of class |
p |
Numeric vector with two components describing a point in the plane or two-column matrix of points. |
The complement of the -convex hull of a sample is calculated by
complement
. The function inahull
checks whether each point in belongs to any of the open balls or halfplanes that define the complement.
in.ahull |
A logical vector specifying whether each point in |
## Not run: # Random sample in the unit square x <- matrix(runif(100), nc = 2) # Value of alpha alpha <- 0.2 # alpha-convex hull ahull.obj <- ahull(x, alpha = alpha) # Check if the point (0.5, 0.5) belongs to the alpha-convex hull inahull(ahull.obj, p = c(0.5, 0.5)) # Check if the points (0.5, 0.5) and (2, 2) belong to the alpha-convex hull inahull(ahull.obj, p = rbind(c(0.5, 0.5), c(2, 2))) ## End(Not run)
## Not run: # Random sample in the unit square x <- matrix(runif(100), nc = 2) # Value of alpha alpha <- 0.2 # alpha-convex hull ahull.obj <- ahull(x, alpha = alpha) # Check if the point (0.5, 0.5) belongs to the alpha-convex hull inahull(ahull.obj, p = c(0.5, 0.5)) # Check if the points (0.5, 0.5) and (2, 2) belong to the alpha-convex hull inahull(ahull.obj, p = rbind(c(0.5, 0.5), c(2, 2))) ## End(Not run)
This function calculates the intersection of two circumferences, given their centers and radius and
, respectively.
inter(c11, c12, r1, c21, c22, r2)
inter(c11, c12, r1, c21, c22, r2)
c11 |
X-coordinate of the center |
c12 |
Y-coordinate of the center |
r1 |
Radius |
c21 |
X-coordinate of the center |
c22 |
Y-coordinate of the center |
r2 |
Radius |
The function inter
is internally called by the function ahull
.
A list with the following components:
n.cut |
Number of intersection points (0,1,2, or Inf). |
v1 |
If there are two intersection points, |
theta1 |
Angle that forms |
v2 |
If there are two intersection points, |
theta2 |
Angle that forms |
This function uses recursion to construct a Kock snowflake curve.
koch(side = 3, niter = 5)
koch(side = 3, niter = 5)
side |
Side length of the initial equilateral triangle. |
niter |
Number of iterations in the development of the snowflake curve. |
The Koch snowflake is a fractal curve described by the Swedish mathematician Helge von Koch in 1904. It is built by starting with an equilateral triangle, removing the inner third of each side, building another equilateral triangle at the location where the side was removed, and then repeating the process.
vertices |
A 2-column matrix with the coordinates of the snowflake vertices. |
von Koch, H. (1904). Sur une courbe continue sans tangente, obtenue par une construction geometrique elementaire. Arkiv for Matematik, 1, pp.681-704.
## Not run: # The first four iterations of a Koch snowflake # with side length of the initial equilateral triangle equal to 3. vertices <- koch(side = 2, niter = 4) plot(vertices[, 1], vertices[, 2], type = "l", asp = TRUE, main = "Koch snowflake", xlab = "", ylab = "", col = 4) polygon(vertices[, 1], vertices[, 2] , col = 4) ## End(Not run)
## Not run: # The first four iterations of a Koch snowflake # with side length of the initial equilateral triangle equal to 3. vertices <- koch(side = 2, niter = 4) plot(vertices[, 1], vertices[, 2], type = "l", asp = TRUE, main = "Koch snowflake", xlab = "", ylab = "", col = 4) polygon(vertices[, 1], vertices[, 2] , col = 4) ## End(Not run)
This function calculates the length of the boundary of the -convex hull of a given sample.
lengthahull(ahull.arcs)
lengthahull(ahull.arcs)
ahull.arcs |
Output matrix of arcs returned by |
The function lengthahull
is internally called by the function ahull
.
length |
Length of the boundary of the |
## Not run: # Random sample in the unit square x <- matrix(runif(100), nc = 2) # Value of alpha alpha <- 0.2 # alpha-convex hull ahull.obj <- ahull(x, alpha = alpha) # Length of the alpha-convex hull ahull.obj$length ## End(Not run)
## Not run: # Random sample in the unit square x <- matrix(runif(100), nc = 2) # Value of alpha alpha <- 0.2 # alpha-convex hull ahull.obj <- ahull(x, alpha = alpha) # Length of the alpha-convex hull ahull.obj$length ## End(Not run)
This function returns a plot of the -convex hull. If desired, it also adds the Delaunay triangulation, Voronoi diagram and
-shape of the sample.
## S3 method for class 'ahull' plot(x, add = FALSE, do.shape = FALSE, wlines = c("none", "both", "del", "vor"), wpoints = TRUE, number = FALSE, col = NULL, xlim = NULL, ylim = NULL, lwd = NULL, ...)
## S3 method for class 'ahull' plot(x, add = FALSE, do.shape = FALSE, wlines = c("none", "both", "del", "vor"), wpoints = TRUE, number = FALSE, col = NULL, xlim = NULL, ylim = NULL, lwd = NULL, ...)
x |
Object of class |
add |
Logical, if TRUE add to a current plot. |
do.shape |
Logical, indicates if the |
wlines |
"Which lines?". I.e. should the Delaunay triangulation be plotted (wlines='del'), should the Voronoi diagram be plotted (wlines='vor'), should both be plotted (wlines='both'), or none (wlines='none', the default)? |
wpoints |
Logical, indicates if sample points should be plotted. |
number |
Logical, defaulting to FALSE; if TRUE then the points plotted will be labelled with their index numbers. |
col |
The colour numbers for plotting the |
xlim |
The limits on the x-axis. |
ylim |
The limits on the y-axis. |
lwd |
The line widths for plotting the tesselations, the |
... |
Arguments to be passed to methods, such as graphical parameters (see |
## Not run: # Random sample in the unit square x <- matrix(runif(100), nc = 2) # Value of alpha alpha <- 0.2 # alpha-convex hull ahull.obj <- ahull(x, alpha = alpha) # Plot including the alpha-convex hull in pink, alpha-shape in blue, # sample points in black, voronoi diagram in green # and Delaunay triangulation in red plot(ahull.obj, do.shape = TRUE, wlines = "both", col = c(6, 4, 1, 2, 3)) # Random sample from a uniform distribution on a Koch snowflake # with initial side length 1 and 3 iterations x <- rkoch(2000, side = 1, niter = 3) # Value of alpha alpha <- 0.05 # Alpha-convex hull ahull.obj <- ahull(x, alpha = alpha) plot(ahull.obj) ## End(Not run)
## Not run: # Random sample in the unit square x <- matrix(runif(100), nc = 2) # Value of alpha alpha <- 0.2 # alpha-convex hull ahull.obj <- ahull(x, alpha = alpha) # Plot including the alpha-convex hull in pink, alpha-shape in blue, # sample points in black, voronoi diagram in green # and Delaunay triangulation in red plot(ahull.obj, do.shape = TRUE, wlines = "both", col = c(6, 4, 1, 2, 3)) # Random sample from a uniform distribution on a Koch snowflake # with initial side length 1 and 3 iterations x <- rkoch(2000, side = 1, niter = 3) # Value of alpha alpha <- 0.05 # Alpha-convex hull ahull.obj <- ahull(x, alpha = alpha) plot(ahull.obj) ## End(Not run)
This function returns a plot of the -shape. If desired, it also adds the Delaunay triangulation and Voronoi diagram of the sample.
## S3 method for class 'ashape' plot(x, add = FALSE, wlines = c("none", "both", "del", "vor"), wpoints = TRUE, number = FALSE, col = NULL, xlim = NULL, ylim = NULL, lwd = NULL, ...)
## S3 method for class 'ashape' plot(x, add = FALSE, wlines = c("none", "both", "del", "vor"), wpoints = TRUE, number = FALSE, col = NULL, xlim = NULL, ylim = NULL, lwd = NULL, ...)
x |
Object of class |
add |
Logical, if TRUE add to a current plot. |
wlines |
"Which lines?". I.e. should the Delaunay triangulation be plotted (wlines='del'), should the Voronoi diagram be plotted (wlines='vor'), should both be plotted (wlines='both'), or none (wlines='none', the default)? |
wpoints |
Logical, indicates if sample points should be plotted. |
number |
Logical, defaulting to FALSE; if TRUE then the points plotted will be labelled with their index numbers. |
col |
The colour numbers for plotting the |
xlim |
The limits on the x-axis. |
ylim |
The limits on the y-axis. |
lwd |
The line widths for plotting the tesselations and the |
... |
Arguments to be passed to methods, such as graphical parameters (see |
## Not run: # Uniform sample of size n=300 in the annulus B(c, 0.5)\B(c, 0.25) # with c=(0.5, 0.5). n <- 300 theta<-runif(n,0,2*pi) r<-sqrt(runif(n,0.25^2,0.5^2)) x<-cbind(0.5+r*cos(theta),0.5+r*sin(theta)) # Value of alpha alpha <- 0.1 # alpha-shape ashape.obj <- ashape(x, alpha = alpha) # Plot alpha-shape in blue, sample points in black, # and Delaunay triangulation in red plot(ashape.obj, wlines= "del", col = c(4, 1, 2)) # Random sample from a uniform distribution on a Koch snowflake # with initial side length 1 and 3 iterations x <- rkoch(2000, side = 1, niter = 3) # Value of alpha alpha <- 0.05 # alpha-shape ashape.obj <- ashape(x, alpha = alpha) # Plot alpha-shape in blue plot(ashape.obj, col = c(4, 1)) ## End(Not run)
## Not run: # Uniform sample of size n=300 in the annulus B(c, 0.5)\B(c, 0.25) # with c=(0.5, 0.5). n <- 300 theta<-runif(n,0,2*pi) r<-sqrt(runif(n,0.25^2,0.5^2)) x<-cbind(0.5+r*cos(theta),0.5+r*sin(theta)) # Value of alpha alpha <- 0.1 # alpha-shape ashape.obj <- ashape(x, alpha = alpha) # Plot alpha-shape in blue, sample points in black, # and Delaunay triangulation in red plot(ashape.obj, wlines= "del", col = c(4, 1, 2)) # Random sample from a uniform distribution on a Koch snowflake # with initial side length 1 and 3 iterations x <- rkoch(2000, side = 1, niter = 3) # Value of alpha alpha <- 0.05 # alpha-shape ashape.obj <- ashape(x, alpha = alpha) # Plot alpha-shape in blue plot(ashape.obj, col = c(4, 1)) ## End(Not run)
This function returns a plot of the Voronoi diagram and Delaunay traingulation.
## S3 method for class 'delvor' plot(x, add = FALSE, wlines = c("both", "del", "vor"), wpoints = TRUE, number = FALSE, col = NULL, xlim = NULL, ylim = NULL, ...)
## S3 method for class 'delvor' plot(x, add = FALSE, wlines = c("both", "del", "vor"), wpoints = TRUE, number = FALSE, col = NULL, xlim = NULL, ylim = NULL, ...)
x |
An object of class |
add |
Logical, if TRUE add to a current plot. |
wlines |
"Which lines?". I.e. should the Delaunay triangulation be plotted (wlines='del'), should the Voronoi diagram be plotted (wlines='vor'), or should both be plotted (wlines='both', the default)? |
wpoints |
Logical, indicates if sample points should be plotted. |
number |
Logical, defaulting to FALSE; if TRUE then the points plotted will be labelled with their index numbers. |
col |
The colour numbers for plotting the data points, Delaunay triangulation, Voronoi diagram, and the point numbers, in that order; defaults to c(1,1,1,1). If fewer than four numbers are given, they are recycled. (If more than four numbers are given, the redundant ones are ignored.) |
xlim |
The limits on the x-axis. |
ylim |
The limits on the y-axis. |
... |
Arguments to be passed to methods, such as graphical parameters (see |
## Not run: # Random sample in the unit square x <- matrix(runif(100), nc = 2) # Delaunay triangulation and Voronoi diagram delvor.obj <- delvor(x) # Plot Voronoi diagram and Delaunay triangulation plot(delvor.obj) ## End(Not run)
## Not run: # Random sample in the unit square x <- matrix(runif(100), nc = 2) # Delaunay triangulation and Voronoi diagram delvor.obj <- delvor(x) # Plot Voronoi diagram and Delaunay triangulation plot(delvor.obj) ## End(Not run)
This function generates ramdom points from a uniform distribution on a Koch snowflake.
rkoch(n, side = 3, niter = 5)
rkoch(n, side = 3, niter = 5)
n |
Number of observations. |
side |
Side length of the initial equilateral triangle of the Koch snowflake curve. |
niter |
Number of iterations in the development of the snowflake curve. |
A 2-column matrix with the coordinates of generated points.
koch
.
## Not run: unifkoch <- rkoch(2000, side = 1, niter = 3) plot(unifkoch, asp = TRUE) ## End(Not run)
## Not run: unifkoch <- rkoch(2000, side = 1, niter = 3) plot(unifkoch, asp = TRUE) ## End(Not run)
This function calculates the clockwise rotation of angle of a given vector
in the plane.
rotation(v, theta)
rotation(v, theta)
v |
Vector |
theta |
Angle |
v.rot |
Vector after rotation. |
## Not run: # Rotation of angle pi/4 of the vector (0,1) rotation(v = c(0, 1), theta = pi/4) ## End(Not run)
## Not run: # Rotation of angle pi/4 of the vector (0,1) rotation(v = c(0, 1), theta = pi/4) ## End(Not run)