| Title: | Global Triangular and Penta-Hexagonal Grids Based on Tessellated Icosahedra |
|---|---|
| Description: | Implementation of icosahedral grids in three dimensions. The spherical-triangular tessellation can be set to create grids with custom resolutions. Both the primary triangular and their inverted penta-hexagonal grids can be calculated. Additional functions are provided that allow plotting of the grids and associated data, the interaction of the grids with other raster and vector objects, and treating the grids as a graphs. |
| Authors: | Adam T. Kocsis [cre, aut] (ORCID: <https://orcid.org/0000-0002-9028-665X>), Deutsche Forschungsgemeinschaft [fnd], FAU GeoZentrum Nordbayern [fnd] |
| Maintainer: | Adam T. Kocsis <[email protected]> |
| License: | GPL-3 |
| Version: | 0.12.0 |
| Built: | 2026-05-26 07:23:21 UTC |
| Source: | https://github.com/icosa-grid/r-icosa |
Shorthand to the subset function.
Function to replace specific elements in a gridlayer object
## S4 method for signature 'gridlayer,ANY,missing' x[i] ## S4 method for signature 'gridlayer,SpatExtent,missing' x[i] ## S4 replacement method for signature 'gridlayer,ANY,ANY' x[i] <- value## S4 method for signature 'gridlayer,ANY,missing' x[i] ## S4 method for signature 'gridlayer,SpatExtent,missing' x[i] ## S4 replacement method for signature 'gridlayer,ANY,ANY' x[i] <- value
x |
( |
i |
( |
value |
The replacement values. |
All these methods are implementing direct replacement in the @values slot of a layer, depending on criteria used for subsetting.
The extraction methods return facelayer-class objects.
This function calculates the shortest arc distance between two points.
arcdist(p1, p2, output = "distance", origin = c(0, 0, 0), radius = authRadius)arcdist(p1, p2, output = "distance", origin = c(0, 0, 0), radius = authRadius)
p1 |
( |
p2 |
( |
output |
( |
origin |
( |
radius |
( |
A single numeric value.
# coordinates of two points point1<- c(0,0) point2<- c(180,0) arcdist(point1,point2,"distance")# coordinates of two points point1<- c(0,0) point2<- c(180,0) arcdist(point1,point2,"distance")
This function calculates the shortest arc distance matrix between two sets of points.
arcdistmat( points1, points2 = NULL, origin = c(0, 0, 0), output = "distance", radius = authRadius )arcdistmat( points1, points2 = NULL, origin = c(0, 0, 0), output = "distance", radius = authRadius )
points1 |
( |
points2 |
( |
origin |
( |
output |
( |
radius |
( |
This function will create all possible shortest arc distances between points in the two sets,
but not between the points within the sets. The function is useful for great circle distance calculations.
For a symmetrical distance matrix leave the points2 argument empty.
A single numeric value.
g <- trigrid(c(4)) res <- arcdistmat(g@vertices) rand<-rpsphere(500) res2 <- arcdistmat(g@vertices, rand)g <- trigrid(c(4)) res <- arcdistmat(g@vertices) rand<-rpsphere(500) res2 <- arcdistmat(g@vertices, rand)
This function calculates points along an arc between two points and a circle center.
arcpoints( p1, p2, breaks = 2, origin = c(0, 0, 0), onlyNew = FALSE, output = "cartesian", radius = authRadius )arcpoints( p1, p2, breaks = 2, origin = c(0, 0, 0), onlyNew = FALSE, output = "cartesian", radius = authRadius )
p1 |
( |
p2 |
( |
breaks |
( |
origin |
( |
onlyNew |
( |
output |
( |
radius |
( |
The function always returns the smaller arc, with angle alpha < pi.
Either an XYZ or a long-lat numeric matrix.
# empty plot plot(NULL, NULL, xlim=c(-180, 180), ylim=c(-90,90)) # then endpoints of the arc point1<-c(-45,-70) point2<-c(130,65) points(arcpoints(point1, point2, breaks=70, output="polar"))# empty plot plot(NULL, NULL, xlim=c(-180, 180), ylim=c(-90,90)) # then endpoints of the arc point1<-c(-45,-70) point2<-c(130,65) points(arcpoints(point1, point2, breaks=70, output="polar"))
Low level plotting of great circle arcs with lines
arcs(x, ...) ## S4 method for signature 'matrix' arcs(x, breaks = 100, breakAtDateline = TRUE, plot = TRUE, ...) ## S4 method for signature 'data.frame' arcs(x, ...)arcs(x, ...) ## S4 method for signature 'matrix' arcs(x, breaks = 100, breakAtDateline = TRUE, plot = TRUE, ...) ## S4 method for signature 'data.frame' arcs(x, ...)
x |
A matrix of longitude and latitude points (WGS 84 longlat) |
... |
Arguments passed to lines (par) |
breaks |
the number of points inserted between every points to draw great circle arcs. |
breakAtDateline |
Logical to indicate whether the lines are to be broken at the dateline. |
plot |
Logical value whether the plotting should be done at all (in case returned values are needed). |
Invisible return of a matrix of coordinates. If breakAtDateline = TRUE, then NA missing values
will be inserted between coordinates where the lines cross the dateline.
# generate random points set.seed(0) example <- rpsphere(10, output="polar") # plotting plot(NULL, NULL, xlim=c(-180, 180), ylim=c(-90,90)) points(example) text(label=1:nrow(example), example, pos=2) arcs(example, col="red", breaks=200)# generate random points set.seed(0) example <- rpsphere(10, output="polar") # plotting plot(NULL, NULL, xlim=c(-180, 180), ylim=c(-90,90)) points(example) text(label=1:nrow(example), example, pos=2) arcs(example, col="red", breaks=200)
The function uses basic trigonometric relationships to transform XYZ coordinates to polar coordinates
CarToPol(x, ...) ## S4 method for signature 'matrix' CarToPol(x, norad = FALSE, origin = c(0, 0, 0)) ## S4 method for signature 'numeric' CarToPol(x, norad = FALSE, origin = c(0, 0, 0)) ## S4 method for signature 'data.frame' CarToPol(x, norad = FALSE, origin = c(0, 0, 0))CarToPol(x, ...) ## S4 method for signature 'matrix' CarToPol(x, norad = FALSE, origin = c(0, 0, 0)) ## S4 method for signature 'numeric' CarToPol(x, norad = FALSE, origin = c(0, 0, 0)) ## S4 method for signature 'data.frame' CarToPol(x, norad = FALSE, origin = c(0, 0, 0))
x |
( |
... |
Arguments passed to class-specific methods. |
norad |
( |
origin |
( |
A 3-column or 2-column numeric, matrix or data.frame with longitude, latitude and, if set accordingly, radius data.
# some random points xyz <- rbind( c(6371, 0,0), c(0, 6371,0), c(1000,1000,1000) ) # conversions CarToPol(xyz)# some random points xyz <- rbind( c(6371, 0,0), c(0, 6371,0), c(1000,1000,1000) ) # conversions CarToPol(xyz)
The function returns which grid faces contain the points clicked in a plot.
cellocator(gridObj, n, output = "faces", ...)cellocator(gridObj, n, output = "faces", ...)
gridObj |
|
n |
( |
output |
( |
... |
Arguments passed to the |
A vector of character values, each corresponding to a face identifier.
Shorthand function to return the @faceCenters slot of an icosahedral grid or a grid linked to a facelayer.
centers(x, ...) ## S4 method for signature 'trigrid' centers(x, output = "polar") ## S4 method for signature 'facelayer' centers(x, output = "polar")centers(x, ...) ## S4 method for signature 'trigrid' centers(x, output = "polar") ## S4 method for signature 'facelayer' centers(x, output = "polar")
x |
( |
... |
Arguments passed to the class specific methods. |
output |
( |
The coordinates of the face centers as a numeric matrix.
a <- trigrid() centers(a)a <- trigrid() centers(a)
This function calculates a possible implementation of the spherical convex hull.
chullsphere( data, center = c(0, 0, 0), radius = authRadius, param = 200, strict = TRUE )chullsphere( data, center = c(0, 0, 0), radius = authRadius, param = 200, strict = TRUE )
data |
( |
center |
( |
radius |
( |
param |
( |
strict |
( |
With the method centroidprojection the function calls the surfacecentroid
function to get the a reference point from the shape. Then all the points are 'projected'
close to this point using the great circles linking them to the reference point.
Each such great circle will be devided to an equal number of points and the closest
will replace the original point coordinates in the convex hull algorithm implemented in chull.
The indices of the data points forming the convex hull as a (numeric) vector.
# generate some random points allData <- rpsphere(1000) # select only a subset points<-allData[allData[,1]>3000,] chullsphere(points)# generate some random points allData <- rpsphere(1000) # select only a subset points<-allData[allData[,1]>3000,] chullsphere(points)
This function will return the length of all edges in the specified grid object.
edgelength(gridObj, ...) ## S4 method for signature 'trigrid' edgelength(gridObj, output = "distance")edgelength(gridObj, ...) ## S4 method for signature 'trigrid' edgelength(gridObj, output = "distance")
gridObj |
( |
... |
Arguments passed to the class specific methods. |
output |
( |
A named numeric vector.
g <- trigrid(3) edges <- edgelength(g, output="deg") edgesg <- trigrid(3) edges <- edgelength(g, output="deg") edges
Shorthand function to get the edges slot of an icosahedral grid or a grid linked to a facelayer.
edges(x) ## S4 method for signature 'obj3d' edges(x) ## S4 method for signature 'facelayer' edges(x)edges(x) ## S4 method for signature 'obj3d' edges(x) ## S4 method for signature 'facelayer' edges(x)
x |
( |
The edges of the grid, as a character matrix.
Make an spdep-style neighbor list for an icosahedral grid
face2nb(x, queen = FALSE)face2nb(x, queen = FALSE)
x |
A trigrid class object |
queen |
Should he queen neighborhood be returned? |
Neighbor-list object such as thoses defined in spdep.
# calculate a grid hex <- hexagrid(deg=5) neighborList <- face2nb(hex) neighborList# calculate a grid hex <- hexagrid(deg=5) neighborList <- face2nb(hex) neighborList
facelayer linked to a trigrid or hexagrid objectThe grids themselves are scaffolds for the assigned data. The data are stored in containers which are linked to the grids.
gridObj |
|
value |
( |
A facelayer class object.
g <- trigrid(c(4,4)) fl <- facelayer(g, 1:length(g)) # faces3d(fl)g <- trigrid(c(4,4)) fl <- facelayer(g, 1:length(g)) # faces3d(fl)
Shorthand function to get the face names of an icosahedral grid or a grid linked to a facelayer.
faces(x) ## S4 method for signature 'trigrid' faces(x) ## S4 method for signature 'gridlayer' faces(x)faces(x) ## S4 method for signature 'trigrid' faces(x) ## S4 method for signature 'gridlayer' faces(x)
x |
( |
The names of the faces of the grid as a character vector.
ball <- hexagrid(2) faces(ball)ball <- hexagrid(2) faces(ball)
This function is used to plot the faces of either a trigrid, hexagrid or facelayer object in 3D space.
faces3d(x, ...) ## S4 method for signature 'trigrid' faces3d(x, ...) ## S4 method for signature 'hexagrid' faces3d(x, ...) ## S4 method for signature 'facelayer' faces3d(x, col = "heat", breaks = NULL, inclusive = TRUE, legend = TRUE, ...)faces3d(x, ...) ## S4 method for signature 'trigrid' faces3d(x, ...) ## S4 method for signature 'hexagrid' faces3d(x, ...) ## S4 method for signature 'facelayer' faces3d(x, col = "heat", breaks = NULL, inclusive = TRUE, legend = TRUE, ...)
x |
|
... |
Further graphical parameters passed to (see |
col |
( |
breaks |
( |
inclusive |
( |
legend |
( |
The function is built on the openGL renderer of the R package rgl.
The function does not return any value.
# create a hexagonal grid g <- hexagrid(c(2,2)) # plot the grid in 3d space # faces3d(g) h <- hexagrid(8) b <- facelayer(h) values(b)<- rnorm(length(b))# create a hexagonal grid g <- hexagrid(c(2,2)) # plot the grid in 3d space # faces3d(g) h <- hexagrid(8) b <- facelayer(h) values(b)<- rnorm(length(b))
Any points set can be binned to an icosahedral grid (i.e. number of incidences can be counted), which will be dependent on the exact positions of grid cells. Rotating the grid in 3d space will result in a different distribution of counts. This distribution can be resampled to a standard orientation structure. The size of the icosahedral grid cells act as a bandwidth parameter.
Discretization with iteratively rotated icosahedral grids
gridensity(x, y, out, trials = 100, FUN = mean) grapply(x, out, ...) ## S4 method for signature 'data.frame,SpatRaster' grapply( x, out, y, coords = c("long", "lat"), iter = 100, FUN = function(x) table(x$cell), APP = mean, miss = NA, APP.args = NULL, counter = TRUE, FUN.args = NULL ) ## S4 method for signature 'data.frame,trigrid' grapply( x, out, y = out, coords = c("long", "lat"), iter = 100, FUN = function(x) table(x$cell), APP = mean, miss = NA, APP.args = NULL, counter = TRUE, FUN.args = NULL ) ## S4 method for signature 'matrix,SpatRaster' grapply(x, out, ...) ## S4 method for signature 'matrix,trigrid' grapply(x, out, ...) ## S4 method for signature 'data.frame,missing' grapply(x, out, y, ...) ## S4 method for signature 'matrix,missing' grapply(x, out, y, ...)gridensity(x, y, out, trials = 100, FUN = mean) grapply(x, out, ...) ## S4 method for signature 'data.frame,SpatRaster' grapply( x, out, y, coords = c("long", "lat"), iter = 100, FUN = function(x) table(x$cell), APP = mean, miss = NA, APP.args = NULL, counter = TRUE, FUN.args = NULL ) ## S4 method for signature 'data.frame,trigrid' grapply( x, out, y = out, coords = c("long", "lat"), iter = 100, FUN = function(x) table(x$cell), APP = mean, miss = NA, APP.args = NULL, counter = TRUE, FUN.args = NULL ) ## S4 method for signature 'matrix,SpatRaster' grapply(x, out, ...) ## S4 method for signature 'matrix,trigrid' grapply(x, out, ...) ## S4 method for signature 'data.frame,missing' grapply(x, out, y, ...) ## S4 method for signature 'matrix,missing' grapply(x, out, y, ...)
x |
Matrix of longitude, latitude data, |
y |
|
out |
|
trials |
|
FUN |
|
... |
Arguments passed to class-specific methods. |
coords |
|
iter |
|
APP |
|
miss |
|
APP.args |
|
counter |
|
FUN.args |
|
The implemented algorithm 1) takes a point cloud (x)) and an icosahedral grid y 2) randomly rotates the icosahedral grid, 3) looks up the points falling on grid cells, 4) resamples the grid to a constant orientation object (either trigrid, hexagrid or SpatRaster). Steps 2-4 are repeated trial times, and then FUN is applied to every vector of values that have same spatial position.
Simple discretization of spatial data is subject to error due to the random assignment to grid cells. grapply offers a framework for the Monte Carlo estimation of the expectation of a function that normally can be applied to
a discretized set of points. This function FUN takes the input the data.frame (or matrix, which will be coerced into one) as argument, with the addition of a new variable cells, which includes
face identifiers for the point set given the random grid rotation (the output of locate). See examples for a step-by-step description.
Either named numeric vector, or a SpatRaster object. If FUN is set to NULL, the output will be either a matrix or SpatRaster.
# example to be run if terra is present if(requireNamespace("terra", quietly=TRUE)){ # randomly generated points x <- rpsphere(100, output="polar") # bandwidth grid gr <- hexagrid(deg=13) # output structure out <- terra::rast(res=5) # Manual example - for understanding what FUN is doing cell<- locate(gr, x) yNew <- cbind(as.data.frame(x), cell=cell) # this is the default function (here named CellCount) CellCount <- function(x) table(x$cell) counts <- CellCount(yNew) # create facelayer fl <- facelayer(gr, counts) # and resample oneOut <- icosa::resample(fl, out) terra::plot(oneOut, main="Density with default grid orientation") points(x, pch=3, col="red") # for density estimation o <- grapply(x=x, out=out,y=gr, iter=7, FUN=CellCount, miss=0) # visualize results terra::plot(o) points(x, pch=3, col="red") }# example to be run if terra is present if(requireNamespace("terra", quietly=TRUE)){ # randomly generated points x <- rpsphere(100, output="polar") # bandwidth grid gr <- hexagrid(deg=13) # output structure out <- terra::rast(res=5) # Manual example - for understanding what FUN is doing cell<- locate(gr, x) yNew <- cbind(as.data.frame(x), cell=cell) # this is the default function (here named CellCount) CellCount <- function(x) table(x$cell) counts <- CellCount(yNew) # create facelayer fl <- facelayer(gr, counts) # and resample oneOut <- icosa::resample(fl, out) terra::plot(oneOut, main="Density with default grid orientation") points(x, pch=3, col="red") # for density estimation o <- grapply(x=x, out=out,y=gr, iter=7, FUN=CellCount, miss=0) # visualize results terra::plot(o) points(x, pch=3, col="red") }
graph class graph from the faces of an icosahedral gridThe function can be applied to both grids and to facelayer-class object of logical values. The resulting graph will have the characteristics of the original grid (directed/undirected etc.).
gridgraph(x, ...) ## S4 method for signature 'trigrid' gridgraph(x, directed = FALSE, distances = FALSE) ## S4 method for signature 'hexagrid' gridgraph(x, directed = FALSE, distances = FALSE) ## S4 method for signature 'facelayer' gridgraph(x)gridgraph(x, ...) ## S4 method for signature 'trigrid' gridgraph(x, directed = FALSE, distances = FALSE) ## S4 method for signature 'hexagrid' gridgraph(x, directed = FALSE, distances = FALSE) ## S4 method for signature 'facelayer' gridgraph(x)
x |
( |
... |
Arguments passed to the class specific methods. |
directed |
|
distances |
|
The function returns an undirected igraph graph.
This function will show where the grid elements are located.
gridlabs(gridObj, type = "f", crs = NULL, ...)gridlabs(gridObj, type = "f", crs = NULL, ...)
gridObj |
|
type |
( |
crs |
( |
... |
Arguments passed to the |
The function has no return value.
gr <- hexagrid(sp=TRUE) plot(gr) gridlabs(gr)gr <- hexagrid(sp=TRUE) plot(gr) gridlabs(gr)
This function will display the names of vertices, faces and edges on 3d plots.
gridlabs3d(gridObj, ...) ## S4 method for signature 'trigrid' gridlabs3d(gridObj, type = "f", ...) ## S4 method for signature 'hexagrid' gridlabs3d(gridObj, type = "f", ...)gridlabs3d(gridObj, ...) ## S4 method for signature 'trigrid' gridlabs3d(gridObj, type = "f", ...) ## S4 method for signature 'hexagrid' gridlabs3d(gridObj, type = "f", ...)
gridObj |
|
... |
Additional arguments passed to |
type |
( |
The function does not return any value.
# create a hexagonal grid g <- hexagrid(c(2,2)) # plot the grid in 3d space # lines3d(g, guides=FALSE) # labels # gridlabs3d(g)# create a hexagonal grid g <- hexagrid(c(2,2)) # plot the grid in 3d space # lines3d(g, guides=FALSE) # labels # gridlabs3d(g)
This function plots 3d guidelines for navigation on the surface of the sphere, includings the rotational axis and a polar coordinate system.
guides3d( axis = 1.5, polgrid = c(30, 30), textPG = FALSE, res = 1, origin = c(0, 0, 0), radius = authRadius, drad = 1.1, ... )guides3d( axis = 1.5, polgrid = c(30, 30), textPG = FALSE, res = 1, origin = c(0, 0, 0), radius = authRadius, drad = 1.1, ... )
axis |
( |
polgrid |
( |
textPG |
( |
res |
( |
origin |
( |
radius |
( |
drad |
( |
... |
Additional arguments passed to |
The function is built on the openGL renderer of the R package rgl.
The function does not return any value.
# create a hexagonal grid g <- hexagrid(c(2,2)) # plot the grid in 3d space # plot3d(g, guides=FALSE) # plot the rotational axis in blue # guides3d(axis=2, polgrid=NULL, col="blue") # plot the polar grid at 10 degree resolution # guides3d(axis=NULL, polgrid=c(10,10), col="red") # plot some coordinates # guides3d(axis=NULL, polgrid=c(30,30), textPG=TRUE, col="orange", cex=1.4)# create a hexagonal grid g <- hexagrid(c(2,2)) # plot the grid in 3d space # plot3d(g, guides=FALSE) # plot the rotational axis in blue # guides3d(axis=2, polgrid=NULL, col="blue") # plot the polar grid at 10 degree resolution # guides3d(axis=NULL, polgrid=c(10,10), col="red") # plot some coordinates # guides3d(axis=NULL, polgrid=c(30,30), textPG=TRUE, col="orange", cex=1.4)
This function will invoke the plot function to draw a heatmap legend.
heatMapLegend( cols, vals, varName, tick.text = NULL, tick.cex = 1.5, barWidth = 3, barHeight = 50, tickLength = 1, xLeft = 88, yBot = 25, add = FALSE, bounds = c(FALSE, FALSE), ... )heatMapLegend( cols, vals, varName, tick.text = NULL, tick.cex = 1.5, barWidth = 3, barHeight = 50, tickLength = 1, xLeft = 88, yBot = 25, add = FALSE, bounds = c(FALSE, FALSE), ... )
cols |
( |
vals |
( |
varName |
( |
tick.text |
( |
tick.cex |
( |
barWidth |
( |
barHeight |
( |
tickLength |
( |
xLeft |
( |
yBot |
( |
add |
( |
bounds |
( |
... |
Arguments passed to the |
The 'percents' refer to the plotting area measured from the lower left corner.
The function does not return any value.
The hexagrid function constrcucts a hexa-pentagonal grid based on the inversion of a
tessellated icosahedron.
tessellation |
( |
deg |
( |
spacing |
( |
sp |
( |
graph |
( |
radius |
( |
center |
( |
verbose |
( |
Inherits from the trigrid class.
The grid structure functions as a frame for data graining, plotting and
calculations. Data can be stored in layers that are linked to the grid object. In the current version only the
facelayer class is implemented which allows the user to render data to the cells
of the grid which are called faces.
The grid 'user interface' is made up of four primary tables: the @vertices table for the coordinates of the vertices,
the faceCenters for the coordinates of the centers of faces,
the faces and the edges tables that contain which vertices form which faces and edges respectively.
In these tables, the faces and vertices are sorted to form spirals that go from the north pole in a counter-clockwise
direction. In case grid subsetting is performed these tables get truncated.
At finer resolutions, the large number of spatial elements render all calculations very resource demanding and slow,
therefore the hierarchical structure created during the tessellation procedure is retained for efficient implementations.
These data are stored in a list in the slot @skeleton and are 0-indexed integer tables for Rccp-based functions. $v
stores vertex, $f the edge, and $e contains the edge data for plotting and calculations. In these tables
the original hierarchy based orderings of the units are retained, during subsetting, additional vectors are used to indicate
deactivation of these units. Any sort of meddling with the @skeleton object will lead to unexpected behavior.
A hexagonal grid object, with class hexagrid.
verticesMatrix of the vertex coordinates.
facesMatrix of the verticies forming the faces
edgesMatrix of the vertices forming the edges.
tessellationContains the tessellation vector.
orientationContains the grid orientation in xyz 3d space, values in radian.
centerThe xyz coordinates of the grid's origin/center.
divContains the number of faces that a single face of the previous tessellation level is decomposed to.
faceCentersContains the xyz coordinates of the centers of the faces on the surface of the sphere.
g <- hexagrid(c(8), sf=TRUE) # based on approximate size (4 degrees edge length) g1 <- hexagrid(deg=4)g <- hexagrid(c(8), sf=TRUE) # based on approximate size (4 degrees edge length) g1 <- hexagrid(deg=4)
hexagrid objectsThe table includes basic properties of hexagrids described with specific tessellation parameters
hexguidehexguide
A data.frame with 120 observations and 20 variables:
totalThe total tessellation of the grid, the number of points inserted between icosahedron vertices along an edge.
level1Level 1 tessellation.
level2Level 2 tessellation - second value of the tessellation vector.
level3Level 3 tessellation - third value of the tessellation vector.
level4Level 4 tessellation - four value of the tessellation vector.
facesThe number of faces in the grid.
verticesThe number of vertices in the grid.
meanEdgeLength_degMean edge length in degrees.
sdEdgeLength_degStandard deviation of edge length in degrees.
meanEdgeLength_kmMean edge length in kilometers.
sdEdgeLength_kmStandard devation of edge length in kilometers.
meanArea_km2Mean face area in square-kilometers.
sdArea_km2Standard deviation of face area in square-kilometers.
timeTime to compute grid with an Intel Xeon E-1650 prcessor.
time_spTime to compute grid with an Intel Xeon E-1650 prcessor, with the 'sp' member.
sizeThe size of the grid in bytes.
size_spThe size of the grid object in bytes, with the 'sp' member.
timeLocate_5000Time to locate 5000 points with an Intel Xeon E-1650 in seconds (wall time).
meanSpacing_degMean spacing in degrees.
sdSpacing_degStandard deviation of spacing in degrees.
This is version 0.2, now including the spacing of the grids.
The function calculates the face names that represent holes in a surface shape
holes(x, ...) ## S4 method for signature 'trigrid' holes(x, y, outside = FALSE) ## S4 method for signature 'facelayer' holes(x)holes(x, ...) ## S4 method for signature 'trigrid' holes(x, y, outside = FALSE) ## S4 method for signature 'facelayer' holes(x)
x |
( |
... |
Arguments passed to class-specific methods. |
y |
( |
outside |
( |
The function uses the horizontal graph of a trigrid-class object, removes the subgraph corresponding to a set of faces (the shape),
and searches for isolated subgraphs. The largest subgraph (highest number of vertices, i.e. faces) is considered to be outside of the shape.
This function relies on the igraph package to run.
A named numeric vector, names correspond to faces, numbers outline the holes. If outside=FALSE and there are no holes in the shape, the function will return NULL.
# create a grid hex <- hexagrid(2, sf=TRUE) # an example shape shape <- paste0("F", c(4, 5, 11, 13, 15, 21, 24, 26, 32, 33, 34, 35, 36)) # visualize basic grid plot(hex) gridlabs(hex) # visualize the shape plot(hex, shape, col="#FF000055", add=TRUE) # calculate holes ho <- holes(x=hex, y=shape) # plot both holes plot(hex, names(ho[ho==1]), add=TRUE, col="#00FF0055") plot(hex, names(ho[ho==2]), add=TRUE, col="#0000FF55")# create a grid hex <- hexagrid(2, sf=TRUE) # an example shape shape <- paste0("F", c(4, 5, 11, 13, 15, 21, 24, 26, 32, 33, 34, 35, 36)) # visualize basic grid plot(hex) gridlabs(hex) # visualize the shape plot(hex, shape, col="#FF000055", add=TRUE) # calculate holes ho <- holes(x=hex, y=shape) # plot both holes plot(hex, names(ho[ho==1]), add=TRUE, col="#00FF0055") plot(hex, names(ho[ho==2]), add=TRUE, col="#0000FF55")
The icosa package provides tools to aggregate and analyze geographic data using grids based on tessellated icosahedra. The procedures can be set to provide a grid with a custom resolution. Both the primary triangular and their inverted penta- hexagonal grids are available for implementation. Additional functions are provided to position points (latitude-longitude data) on the grids, to allow 2D and 3D plotting, use raster and vector spatial data.
Note that similar to R, the package comes with absolutely no warranty. Notes about found bugs and suggestions are more than welcome!
Adam T. Kocsis ([email protected])
Useful links:
# Create a triangular grid tri <- trigrid(c(2,2))# Create a triangular grid tri <- trigrid(c(2,2))
trigrid or hexagrid class object.The length of the object is interpreted as the number of faces it contains.
## S4 method for signature 'trigrid' length(x) ## S4 method for signature 'gridlayer' length(x)## S4 method for signature 'trigrid' length(x) ## S4 method for signature 'gridlayer' length(x)
x |
An integer value.
trigrid and hexagrid classesThis function will invoke the method of the SpatialPolygons class.
This function will invoke the lines method of the sf or the SpatialPolygons class.
## S4 method for signature 'trigrid' lines(x, crs = NULL, col = 1, lwd = 1, lty = 1, ...)## S4 method for signature 'trigrid' lines(x, crs = NULL, col = 1, lwd = 1, lty = 1, ...)
x |
|
crs |
( |
col |
Line colors - as in |
lwd |
Line thickness - as in |
lty |
Line type - as in |
... |
Arguments passed to the |
The function has no return value.
This is a generic function used to plot the edge lines of either a trigrid or a hexagrid object, a facelayer, or Spatial objects in 3d space. The method is also implemented for
the object classes defined by the package 'sp'.
lines3d ## S4 method for signature 'trigrid' lines3d(x, arcs = FALSE, ...) ## S4 method for signature 'Line' lines3d(x, radius = authRadius, ...) ## S4 method for signature 'Lines' lines3d(x, radius = authRadius, ...) ## S4 method for signature 'SpatialLines' lines3d(x, radius = authRadius, ...) ## S4 method for signature 'SpatialLinesDataFrame' lines3d(x, radius = authRadius, ...) ## S4 method for signature 'Polygon' lines3d(x, radius = authRadius, ...) ## S4 method for signature 'Polygons' lines3d(x, radius = authRadius, ...) ## S4 method for signature 'SpatialPolygons' lines3d(x, radius = authRadius, ...) ## S4 method for signature 'SpatialPolygonsDataFrame' lines3d(x, radius = authRadius, ...)lines3d ## S4 method for signature 'trigrid' lines3d(x, arcs = FALSE, ...) ## S4 method for signature 'Line' lines3d(x, radius = authRadius, ...) ## S4 method for signature 'Lines' lines3d(x, radius = authRadius, ...) ## S4 method for signature 'SpatialLines' lines3d(x, radius = authRadius, ...) ## S4 method for signature 'SpatialLinesDataFrame' lines3d(x, radius = authRadius, ...) ## S4 method for signature 'Polygon' lines3d(x, radius = authRadius, ...) ## S4 method for signature 'Polygons' lines3d(x, radius = authRadius, ...) ## S4 method for signature 'SpatialPolygons' lines3d(x, radius = authRadius, ...) ## S4 method for signature 'SpatialPolygonsDataFrame' lines3d(x, radius = authRadius, ...)
x |
|
arcs |
|
... |
Further graphical parameters passed to (see |
radius |
( |
An object of class nonstandardGenericFunction of length 1.
The function is built on the openGL renderer of the R package rgl, which needs to be installed for the function to run. Although the function is works without attaching rgl, note that if you want to attach both icosa and rgl,the rgl package has to be loaded ifrst otherwise the function will not be usable.
The function does not return any value.
# create a hexagonal grid g <- hexagrid(c(2,2)) # plot the grid in 3d space # lines3d(g, col="blue")# create a hexagonal grid g <- hexagrid(c(2,2)) # plot the grid in 3d space # lines3d(g, col="blue")
Basic lookup function of coordinates on an icosahedral grid
locate(x, y, ...) ## S4 method for signature 'trigrid,matrix' locate(x, y, randomborder = FALSE, output = "ui") ## S4 method for signature 'trigrid,numeric' locate(x, y, ...) ## S4 method for signature 'trigrid,data.frame' locate(x, y, ...) ## S4 method for signature 'trigrid,sf' locate(x, y, ...) ## S4 method for signature 'trigrid,SpatialPoints' locate(x, y, ...) ## S4 method for signature 'trigrid,SpatialPointsDataFrame' locate(x, y, ...) ## S4 method for signature 'hexagrid,matrix' locate(x, y, output = "ui", randomborder = FALSE, forceNA = FALSE)locate(x, y, ...) ## S4 method for signature 'trigrid,matrix' locate(x, y, randomborder = FALSE, output = "ui") ## S4 method for signature 'trigrid,numeric' locate(x, y, ...) ## S4 method for signature 'trigrid,data.frame' locate(x, y, ...) ## S4 method for signature 'trigrid,sf' locate(x, y, ...) ## S4 method for signature 'trigrid,SpatialPoints' locate(x, y, ...) ## S4 method for signature 'trigrid,SpatialPointsDataFrame' locate(x, y, ...) ## S4 method for signature 'hexagrid,matrix' locate(x, y, output = "ui", randomborder = FALSE, forceNA = FALSE)
x |
( |
y |
( |
... |
Arguments passed to class specific methods. |
randomborder |
( |
output |
( |
forceNA |
( |
The function returns the cell names (as character) where the input coordinates fall.
# create a grid g <- trigrid(4) # some random points randomPoints<-rpsphere(4, output="polar") # cells locate(g, randomPoints)# create a grid g <- trigrid(4) # some random points randomPoints<-rpsphere(4, output="polar") # cells locate(g, randomPoints)
facelayer class objectFunction to extract the registered face names to which the facelayer renders information.
## S4 method for signature 'gridlayer' names(x)## S4 method for signature 'gridlayer' names(x)
x |
( |
A vector of character values, the names of the faces.
Add an igraph object to a predefined slot in an icosahedral grid
newgraph(gridObj, ...) ## S4 method for signature 'trigrid' newgraph(gridObj, ...)newgraph(gridObj, ...) ## S4 method for signature 'trigrid' newgraph(gridObj, ...)
gridObj |
|
... |
Arguments passed to the |
A new (trigrid or hexagrid) object with the recalculated graph.
#create a grid g<-trigrid(4, graph=FALSE) g<-newgraph(g)#create a grid g<-trigrid(4, graph=FALSE) g<-newgraph(g)
sf object to a predefined slot in a trigrid or hexagrid objectAdd a sf object to a predefined slot in a trigrid or hexagrid object
newsf(x, res = NULL) ## S4 method for signature 'trigrid' newsf(x, res = NULL)newsf(x, res = NULL) ## S4 method for signature 'trigrid' newsf(x, res = NULL)
x |
|
res |
( |
A trigrid or hexagrid object with the new @sf slot.
a<-trigrid(4) a<-newsf(a) plot(a)a<-trigrid(4) a<-newsf(a) plot(a)
SpatialPolygons object to a predefined slot in a trigrid or hexagrid objectAdd a SpatialPolygons object to a predefined slot in a trigrid or hexagrid object
newsp(gridObj, res = NULL) ## S4 method for signature 'trigrid' newsp(gridObj, res = NULL)newsp(gridObj, res = NULL) ## S4 method for signature 'trigrid' newsp(gridObj, res = NULL)
gridObj |
|
res |
( |
A trigrid or hexagrid object with the new @sp slot.
a<-trigrid(4) a<-newsp(a) plot(a)a<-trigrid(4) a<-newsp(a) plot(a)
This function will return an object showing which faces are occupied by the input object.
occupied(gridObj, data, out = "logical", ...)occupied(gridObj, data, out = "logical", ...)
gridObj |
|
data |
( |
out |
( |
... |
Arguments passed to the class specific methods |
This is a wrapper function on the OccupiedFaces methods that are specific to grid class and input data.
The function returns a either a named logical vector or facelayer-class object.
# create a grid g <- trigrid(8, sf=TRUE) # create random points randPoints <- rpsphere(100,output="polar") # the faces occupied by these points occ <- occupied(g, randPoints) # plot using sf slot independently plot(g@sf[occ,"geometry"]) points(randPoints, col="red", pch="+")# create a grid g <- trigrid(8, sf=TRUE) # create random points randPoints <- rpsphere(100,output="polar") # the faces occupied by these points occ <- occupied(g, randPoints) # plot using sf slot independently plot(g@sf[occ,"geometry"]) points(randPoints, col="red", pch="+")
Extracting and setting the grid orientation
orientation(x, ...) ## S4 method for signature 'trigrid' orientation(x, display = "deg", ...) orientation(x) <- value ## S4 replacement method for signature 'trigrid' orientation(x) <- valueorientation(x, ...) ## S4 method for signature 'trigrid' orientation(x, display = "deg", ...) orientation(x) <- value ## S4 replacement method for signature 'trigrid' orientation(x) <- value
x |
|
... |
Values passed on to the |
display |
( |
value |
( |
In case the function returns does, it returns the orientation angles of the grid (as numeric).
The function calculates the face names that represent patches in a surface shape
## S4 method for signature 'trigrid' patches(x, y, ...) ## S4 method for signature 'facelayer' patches(x)## S4 method for signature 'trigrid' patches(x, y, ...) ## S4 method for signature 'facelayer' patches(x)
x |
( |
y |
( |
... |
Arguments passed to class-specific methods. |
The function uses the horizontal graph of a trigrid-class object, and searches for isolated subgraphs.
This function relies on the igraph package to run.
A named numeric vector, names correspond to faces, numbers define the patches.
# create a grid hex <- hexagrid(2, sf=TRUE) # an example shape shape <- paste0("F", c(3,6,7,9, 10, 16, 22, 26)) # visualize basic grid plot(hex) gridlabs(hex) # visualize the shape plot(hex, shape, col="#FF000055", add=TRUE) # calculate holes pa <- patches(x=hex, y=shape) # plot all patches (coloring borders) plot(hex, names(pa[pa==1]), add=TRUE, border="#00FF00", lwd=4) plot(hex, names(pa[pa==2]), add=TRUE, border="#0000FF", lwd=4) plot(hex, names(pa[pa==3]), add=TRUE, border="#00FFFF", lwd=4)# create a grid hex <- hexagrid(2, sf=TRUE) # an example shape shape <- paste0("F", c(3,6,7,9, 10, 16, 22, 26)) # visualize basic grid plot(hex) gridlabs(hex) # visualize the shape plot(hex, shape, col="#FF000055", add=TRUE) # calculate holes pa <- patches(x=hex, y=shape) # plot all patches (coloring borders) plot(hex, names(pa[pa==1]), add=TRUE, border="#00FF00", lwd=4) plot(hex, names(pa[pa==2]), add=TRUE, border="#0000FF", lwd=4) plot(hex, names(pa[pa==3]), add=TRUE, border="#00FFFF", lwd=4)
trigrid, hexagrid or facelayer classesThis function will invoke the plot method of the sf or the SpatialPolygons class.
The function passes arguments to the plot method of the SpatialPolygons class. In case a heatmap is plotted and the plotting device gets resized,
some misalignments can happen. If you want to use a differently sized window, use x11 to set the height and width before running the function.
The function matches data referred to the grid and plots it with sf's plotting methods.
plot ## S4 method for signature 'trigrid,ANY' plot(x, crs = NULL, ...) ## S4 method for signature 'facelayer,ANY' plot( x, crs = NULL, col = "heat", border = NA, alpha = NULL, frame = FALSE, legend = TRUE, breaks = NULL, inclusive = TRUE, discrete = FALSE, ... ) ## S4 method for signature 'trigrid,numeric' plot(x, y, crs = NULL, main = "", ...) ## S4 method for signature 'trigrid,array' plot(x, y, crs = NULL, main = "", ...) ## S4 method for signature 'trigrid,table' plot(x, y, crs = NULL, main = "", ...) ## S4 method for signature 'trigrid,character' plot(x, y, crs = NULL, main = "", ...) ## S4 method for signature 'trigrid,logical' plot(x, y, crs = NULL, main = "", ...)plot ## S4 method for signature 'trigrid,ANY' plot(x, crs = NULL, ...) ## S4 method for signature 'facelayer,ANY' plot( x, crs = NULL, col = "heat", border = NA, alpha = NULL, frame = FALSE, legend = TRUE, breaks = NULL, inclusive = TRUE, discrete = FALSE, ... ) ## S4 method for signature 'trigrid,numeric' plot(x, y, crs = NULL, main = "", ...) ## S4 method for signature 'trigrid,array' plot(x, y, crs = NULL, main = "", ...) ## S4 method for signature 'trigrid,table' plot(x, y, crs = NULL, main = "", ...) ## S4 method for signature 'trigrid,character' plot(x, y, crs = NULL, main = "", ...) ## S4 method for signature 'trigrid,logical' plot(x, y, crs = NULL, main = "", ...)
x |
|
crs |
( |
... |
Arguments passed to the |
col |
( |
border |
( |
alpha |
( |
frame |
( |
legend |
( |
breaks |
( |
inclusive |
( |
discrete |
( |
y |
Data or part of the grid to be plotted. If it is an unnamed character vector, then it is expected to be
a set of faces, which will be treated as a subscript that indicates the faces to be plotted.
If it is a logical vector, then it is expeced to be subscript, indicating a similar operation.
If it is a named logical or character vector, table with names, or single-dimensional named array
then the names are expected to refer to faces of the grid |
main |
The main title of the plot |
An object of class standardGeneric of length 1.
The function has no return value.
# A simple grid, with sf-representation gr <- hexagrid(4, sf=TRUE) dat <- 1:nrow(gr@faces) names(dat) <- paste0("F", dat) plot(x=gr, y=dat)# A simple grid, with sf-representation gr <- hexagrid(4, sf=TRUE) dat <- 1:nrow(gr@faces) names(dat) <- paste0("F", dat) plot(x=gr, y=dat)
The function is built on the openGL renderer of the R package rgl. The default plotting window size is 800x800 pixels. In case you want to override this, please
use the function with defaultPar3d=FALSE after running par3d(windowRect=<>).
plot3d(x,...) ## S3 method for class 'trigrid' plot3d(x, type = c("l"), sphere = NULL, add = FALSE, guides = TRUE, ...) ## S3 method for class 'hexagrid' plot3d( x, type = c("l"), sphere = NULL, color = "gray70", add = FALSE, guides = TRUE, ... ) ## S3 method for class 'facelayer' plot3d(x, type = "f", frame = TRUE, guides = TRUE, defaultPar3d = TRUE, ...)plot3d(x,...) ## S3 method for class 'trigrid' plot3d(x, type = c("l"), sphere = NULL, add = FALSE, guides = TRUE, ...) ## S3 method for class 'hexagrid' plot3d( x, type = c("l"), sphere = NULL, color = "gray70", add = FALSE, guides = TRUE, ... ) ## S3 method for class 'facelayer' plot3d(x, type = "f", frame = TRUE, guides = TRUE, defaultPar3d = TRUE, ...)
x |
|
type |
( |
sphere |
( |
add |
( |
guides |
( |
... |
Further graphical parameters passed to (see |
color |
( |
frame |
( |
defaultPar3d |
( |
An object of class function of length 1.
The function does not return any value.
# create a hexagonal grid g <- hexagrid(c(2,2)) # plot the grid in 3d space # plot3d(g, col="blue") # make a subset to select faces subG <- subset(g, c("F5", "F2")) # plot the subset defined above # plot3d(subG, type="f", col=c("orange"), add=TRUE, lwd=1)# create a hexagonal grid g <- hexagrid(c(2,2)) # plot the grid in 3d space # plot3d(g, col="blue") # make a subset to select faces subG <- subset(g, c("F5", "F2")) # plot the subset defined above # plot3d(subG, type="f", col=c("orange"), add=TRUE, lwd=1)
The function uses basic trigonometric relationships to transform longitude/latitude coordinates on a sphere to xyz Cartesian coordinates.
PolToCar(x, ...) ## S4 method for signature 'matrix' PolToCar(x, radius = authRadius, origin = c(0, 0, 0)) ## S4 method for signature 'numeric' PolToCar(x, radius = authRadius, origin = c(0, 0, 0)) ## S4 method for signature 'data.frame' PolToCar(x, radius = authRadius, origin = c(0, 0, 0), long = NULL, lat = NULL)PolToCar(x, ...) ## S4 method for signature 'matrix' PolToCar(x, radius = authRadius, origin = c(0, 0, 0)) ## S4 method for signature 'numeric' PolToCar(x, radius = authRadius, origin = c(0, 0, 0)) ## S4 method for signature 'data.frame' PolToCar(x, radius = authRadius, origin = c(0, 0, 0), long = NULL, lat = NULL)
x |
( |
... |
Arguments passed to class-specific methods. |
radius |
( |
origin |
( |
long |
( |
lat |
( |
The authalic mean radius of Earth (6371.007 km) is used by this function as a default. The origin is c(0,0,0). The precision of these conversions is not exact (see example c(0,90) below),
but should be considered acceptable when applied at a reasonable scale (e.g. for global analyses using data above 10e-6 meters of resolution).
An xyz 3-column numeric matrix, data.frame or numeric, depending on the class of x.
longLat <- rbind( c(0,0), #note the precision here! c(0, 90), c(-45,12) ) # matrix-method xyz <- PolToCar(longLat) # numeric-method xyz2 <- PolToCar(longLat[1,]) # data.frame method xyz3 <- PolToCar(as.data.frame(longLat))longLat <- rbind( c(0,0), #note the precision here! c(0, 90), c(-45,12) ) # matrix-method xyz <- PolToCar(longLat) # numeric-method xyz2 <- PolToCar(longLat[1,]) # data.frame method xyz3 <- PolToCar(as.data.frame(longLat))
This function will retrieve the position of a vertex or a face on a hexagrid or trigrid object.
pos(gridObj, names, output = "polar")pos(gridObj, names, output = "polar")
gridObj |
|
names |
( |
output |
( |
Vertex and face names can be mixed in a single names argument.
A numeric matrix.
g <- trigrid(c(4,4)) pos(g, c("F2", "P6", "dummyname"))g <- trigrid(c(4,4)) pos(g, c("F2", "P6", "dummyname"))
trigrid or a hexagrid object.The function is used to resolve and resample data stored in SpatRasters and facelayers so they can be fitted to and can be plotted by using trigrid or hexagrid objects.
The function applies different resampling algorithms. Currently there are only two implemented methods, one for upscaling and one for downscaling. The downscaling method "average" will tabluate all face centers from the high resolution grid that fall on a coarse resolution cell and average them. The upscaling method "ebaa" (edge breakpoint area approximation) will estimate the areas covered by the high resolution cells using the number of edge breakpoints.
resample ## S4 method for signature 'SpatRaster,trigrid' resample(x, y, method = "near", na.rm = TRUE, output = "numeric") ## S4 method for signature 'facelayer,trigrid' resample(x, y, method = NULL, res = 5) ## S4 method for signature 'facelayer,SpatRaster' resample(x, y, method = NULL, res = 5) ## S4 method for signature 'SpatRaster,facelayer' resample(x, y, ...)resample ## S4 method for signature 'SpatRaster,trigrid' resample(x, y, method = "near", na.rm = TRUE, output = "numeric") ## S4 method for signature 'facelayer,trigrid' resample(x, y, method = NULL, res = 5) ## S4 method for signature 'facelayer,SpatRaster' resample(x, y, method = NULL, res = 5) ## S4 method for signature 'SpatRaster,facelayer' resample(x, y, ...)
x |
( |
y |
( |
method |
( |
na.rm |
( |
output |
( |
res |
( |
... |
Arguments passed to class-specific methods. |
An object of class standardGeneric of length 1.
This method is necessary to utilize rasterized data in the icosa package. The only method currently implemented upscales the raster data and then resolves the values to the trigrid or hexagrid values, using averages. In the case of resampling SpatRasters, the method argument will be passed to the resample function.
A named numeric vector.
# create a grid g <- trigrid(c(4,4)) # create a data layer fl <- facelayer(g, rnorm(length(faces(g)))) # target structure h <- trigrid(4) # resampling res <- resample(fl, h) fl2<-facelayer(h, res) fl2@values[] <- res# create a grid g <- trigrid(c(4,4)) # create a data layer fl <- facelayer(g, rnorm(length(faces(g)))) # target structure h <- trigrid(4) # resampling res <- resample(fl, h) fl2<-facelayer(h, res) fl2@values[] <- res
Multiple implementations of rotations for coordinate transformation.
## S4 method for signature 'matrix' rotate( x, angles = "random", long = 0, lat = 0, reflong = NULL, pivot = c(0, 0, 0), radius = authRadius, output = "polar" ) ## S4 method for signature 'data.frame' rotate(x, coords = NULL, ...) rotate ## S4 method for signature 'trigrid' rotate(x, angles = "random", pivot = NA, projnote = TRUE)## S4 method for signature 'matrix' rotate( x, angles = "random", long = 0, lat = 0, reflong = NULL, pivot = c(0, 0, 0), radius = authRadius, output = "polar" ) ## S4 method for signature 'data.frame' rotate(x, coords = NULL, ...) rotate ## S4 method for signature 'trigrid' rotate(x, angles = "random", pivot = NA, projnote = TRUE)
x |
|
angles |
( |
long |
( |
lat |
( |
reflong |
( |
pivot |
( |
radius |
The radius of the sphere, relevant only if the |
output |
The output format of the rotations, either |
coords |
( |
... |
Arguments passed to class-specific methods. |
projnote |
( |
An object of class standardGeneric of length 1.
The function implements 3D rotations of various class of objects, that are ultimately reduced to individual points. Internally, point rotation is implemented with 3-axis rotations (Method 1), that are implemented in the X-Y-Z order (note that 3d rotations are not commutative!). For this reason it is not recommended to re-rotate an already rotated object, unless the purpose is to achieve random orientation.
Method 2 parametrizes rotation with three arguments (long, lat, reflong), that can be easier to control. Longitudinal rotations are invariant to the position of the object, but latitudinal
rotation is not. An axis in the equatorial plane will be set perpendicular to the reference longitude (reflong)for latitude-based rotation
(i.e. the latitude difference will equal the latitudinal rotation value only at the reference longitude). If this is not given, then the reference longitude will be that of the centroid of the point cloud.
Note that latitudinal rotation is executed first and only then are the points rotated longitudinally.
Same class object as x.
This function will create a predefined number of points randomly distributed on the surface of a sphere with a given radius.
rpsphere(n = 1, output = "cartesian", radius = authRadius, origin = c(0, 0, 0))rpsphere(n = 1, output = "cartesian", radius = authRadius, origin = c(0, 0, 0))
n |
( |
output |
( |
radius |
( |
origin |
( |
The function uses a three dimensional gaussian distribution to generate points, which are then projected to the surface of the sphere.
A 3-column (XYZ) or a 2-column (long-lat) numeric matrix.
randomPoints <- rpsphere(2000, output="polar") # observe latitudinal pattern plot(randomPoints, xlim=c(-180, 180), ylim=c(-90, 90))randomPoints <- rpsphere(2000, output="polar") # observe latitudinal pattern plot(randomPoints, xlim=c(-180, 180), ylim=c(-90, 90))
The function will take the given trigrid class object and write it's vertex, edge and face information as a .obj file
saveOBJ(x, ...) ## S4 method for signature 'trigrid' saveOBJ(x, file, scale = TRUE) ## S4 method for signature 'hexagrid' saveOBJ(x, file, scale = TRUE)saveOBJ(x, ...) ## S4 method for signature 'trigrid' saveOBJ(x, file, scale = TRUE) ## S4 method for signature 'hexagrid' saveOBJ(x, file, scale = TRUE)
x |
A |
... |
Arguments of class-specific methods. |
file |
A |
scale |
A |
Note that hexagrid class objects are exported in their triangulated form (subfaces). The order of faces for hexagrids is not the natural (UI) order but the internal order of subfaces.
The function has no return value.
gr <- hexagrid(spacing=4) # example written into temporary directroy td <- tempdir() td # actual writing saveOBJ(gr, file=file.path(td, "hexagrid.obj"))gr <- hexagrid(spacing=4) # example written into temporary directroy td <- tempdir() td # actual writing saveOBJ(gr, file=file.path(td, "hexagrid.obj"))
This function will return the distance between neighboring face centers.
spacing(x, ...) ## S4 method for signature 'trigrid' spacing(x, degree = TRUE)spacing(x, ...) ## S4 method for signature 'trigrid' spacing(x, degree = TRUE)
x |
|
... |
Arguments of class-specific methods. |
degree |
( |
The value for every pair is given in either degrees or kilometers depending on degree.
A named numeric vector, one value for every for every neighboring cell pair.
h <- hexagrid(3) spacing(h)h <- hexagrid(3) spacing(h)
SpatialLines class object from an icosahedral gridCreate a SpatialLines class object from an icosahedral grid
SpLines(gridObj, ...) ## S4 method for signature 'trigrid' SpLines(gridObj, dateLine = "break", res = NULL)SpLines(gridObj, ...) ## S4 method for signature 'trigrid' SpLines(gridObj, dateLine = "break", res = NULL)
gridObj |
|
... |
Specific details of the new |
dateLine |
( |
res |
( |
An object of class SpatialLines.
The function will create a SpatialPolygons class 2d representation of the icosahedral grid.
SpPolygons(gridObj, ...) ## S4 method for signature 'trigrid' SpPolygons(gridObj, res = NULL) ## S4 method for signature 'hexagrid' SpPolygons(gridObj, res = NULL)SpPolygons(gridObj, ...) ## S4 method for signature 'trigrid' SpPolygons(gridObj, res = NULL) ## S4 method for signature 'hexagrid' SpPolygons(gridObj, res = NULL)
gridObj |
|
... |
Arguments passed to class-specific methods. |
res |
( |
A SpatialPolygons class object.
a <- trigrid() sp <- SpPolygons(a)a <- trigrid() sp <- SpPolygons(a)
This is a generic function used to access data from either a triangular or hexagonal grid using the names of the faces, integers or logical vectors.
The function extracts subsets of the gridlayer depending on different criteria.
subset ## S4 method for signature 'trigrid' subset(x, i) ## S4 method for signature 'hexagrid' subset(x, i) ## S4 method for signature 'trigrid,ANY,ANY' x[i] ## S4 method for signature 'gridlayer' subset(x, i)subset ## S4 method for signature 'trigrid' subset(x, i) ## S4 method for signature 'hexagrid' subset(x, i) ## S4 method for signature 'trigrid,ANY,ANY' x[i] ## S4 method for signature 'gridlayer' subset(x, i)
x |
( |
i |
( |
An object of class standardGeneric of length 1.
The function returns subsets of the grid pertaining to the specified faces that can be used for additional operations (e.g. plotting).
The subscript vector can be either a logical, character or numeric one. The character vector should contain the names of faces, the logical subscript should have
the same length as the number of faces in the order in which the faces are present in the faces slot.
The numeric vector can either refer to indices to the rownames of faces in the faces slot, or
to surfaces bounded by longitude/latitude data. In the latter case, the the vector should contain an element with a names of at least one of the "lomax", "lamax",
"lomin" or "lamin" strings (lo for longitude, la: latitude, min: minimum, max: maximum). In case a subset around the dateline is needed a larger longitude to a smaller longitude value is needed (e.g. between 150° to -150°).
The following methods are incorporated into the function: If i argument is a vector of integers, they will be interpreted as indices. If the numeric i contains either the lamin, lamax, lomin or lomax names, the subsetting will be done using the latitude-longitude coordinates outlined by these 4 values. Logical subsetting and subsetting by face names are also possible.
Subset of the input grid. The class of the original object is retained, the @skeleton slot contains all previous information.
#create a triangular grid g <- trigrid(c(2,2)) #make a subset pertaining to the faces subG1 <- subset(g, c("F1", "F33")) #additional way of subsetting subG2 <- g[1:15] # selects faces F1 through F15 logicalSub<-sample(c(TRUE,FALSE), nrow(g@faces), replace=TRUE) subG3 <- g[logicalSub] #plot the subset in 3d space # plot3d(subG3) # previously mentioned case around the dateline gDateLine<-g[c(lomax=-150, lomin=150)] # plot3d(gDateLine)#create a triangular grid g <- trigrid(c(2,2)) #make a subset pertaining to the faces subG1 <- subset(g, c("F1", "F33")) #additional way of subsetting subG2 <- g[1:15] # selects faces F1 through F15 logicalSub<-sample(c(TRUE,FALSE), nrow(g@faces), replace=TRUE) subG3 <- g[logicalSub] #plot the subset in 3d space # plot3d(subG3) # previously mentioned case around the dateline gDateLine<-g[c(lomax=-150, lomin=150)] # plot3d(gDateLine)
This function will return the areas of all cells in the specified grid object.
surfacearea(x) ## S4 method for signature 'trigrid' surfacearea(x) ## S4 method for signature 'hexagrid' surfacearea(x)surfacearea(x) ## S4 method for signature 'trigrid' surfacearea(x) ## S4 method for signature 'hexagrid' surfacearea(x)
x |
A named numeric vector, in the metric that was given to the function in the coordinates or the radius of the grid. Default grid configurations yield values in square kilometers.
g <- trigrid(3) surfaces <- surfacearea(g) surfacesg <- trigrid(3) surfaces <- surfacearea(g) surfaces
This function calculated the projected place of the centroid from a pointset on the sphere.
surfacecentroid(x, ...) ## S4 method for signature 'matrix' surfacecentroid( x, output = "polar", center = c(0, 0, 0), radius = authRadius, w = NULL ) ## S4 method for signature 'data.frame' surfacecentroid(x, ...) ## S4 method for signature 'SpatialPoints' surfacecentroid(x, ...)surfacecentroid(x, ...) ## S4 method for signature 'matrix' surfacecentroid( x, output = "polar", center = c(0, 0, 0), radius = authRadius, w = NULL ) ## S4 method for signature 'data.frame' surfacecentroid(x, ...) ## S4 method for signature 'SpatialPoints' surfacecentroid(x, ...)
x |
( |
... |
Arguments passed to the |
output |
( |
center |
( |
radius |
( |
w |
( |
The function calculates the position of the centroid in 3D space (inside the sphere/Earth), which is then projected to the surface.
Either an XYZ or a long-lat numeric vector.
# generate some random points allData <- rpsphere(1000) # select only a subset points<-allData[allData[,3]>1500,] # transform to 2d points2 <- CarToPol(points, norad=TRUE) # the spherical centroid sc <- surfacecentroid(points2, output="polar") sc #3d plot plot(points2, xlim=c(-180, 180), ylim=c(-90, 90)) points(sc[1], sc[2], col="red", cex=5, pch=3)# generate some random points allData <- rpsphere(1000) # select only a subset points<-allData[allData[,3]>1500,] # transform to 2d points2 <- CarToPol(points, norad=TRUE) # the spherical centroid sc <- surfacecentroid(points2, output="polar") sc #3d plot plot(points2, xlim=c(-180, 180), ylim=c(-90, 90)) points(sc[1], sc[2], col="red", cex=5, pch=3)
The function translates the coordinates of a grid object with the specified 3d vector.
translate(gridObj, vec) ## S4 method for signature 'trigrid,numeric' translate(gridObj, vec) ## S4 method for signature 'hexagrid,numeric' translate(gridObj, vec)translate(gridObj, vec) ## S4 method for signature 'trigrid,numeric' translate(gridObj, vec) ## S4 method for signature 'hexagrid,numeric' translate(gridObj, vec)
gridObj |
|
vec |
( |
The same grid structure as the input, but with translated coordinates.
# create a grid and plot it g <- trigrid(3) # lines3d(g) # translate the grid to (15000,15000,15000) g2 <- translate(g, c(15000,15000,15000)) # lines3d(g2)# create a grid and plot it g <- trigrid(3) # lines3d(g) # translate the grid to (15000,15000,15000) g2 <- translate(g, c(15000,15000,15000)) # lines3d(g2)
trigrid() creates a triangular grid based on the
tessellation of an icosahedron.
tessellation |
( |
deg |
( |
spacing |
( |
sp |
( |
graph |
( |
radius |
( |
center |
( |
verbose |
( |
The grid structure functions as a frame for data graining, plotting and spatial
calculations. Data can be stored in layers that are linked to the grid object. In the current version only the
facelayer class is implemented, which allows the user to render data to the cells
of the grid, which are usually referred to as faces.
The grid 'user interface' is made up of four primary tables: the @vertices table for the coordinates of the vertices,
the faceCenters for the coordinates of the centers of faces,
the faces and the edges tables that contain which vertices form which faces and edges respectively.
In these tables, the faces and vertices are sorted to form spirals that go from the north pole in a counter-clockwise
direction. In case grid subsetting is performed these tables get truncated.
At finer resolutions, the large number of spatial elements render all calculations resource demanding and slow,
therefore the hierarchical structure created during the tessellation procedure is retained for efficient implementation.
These data are stored in a list in the slot @skeleton and are 0-indexed integer tables for Rccp-based functions. $v
stores vertex, $f the edge, and $e contains the edge data for plotting and calculations. In these tables
the original hierarchy based orderings of the units are retained, during subsetting, additional vectors are used to indicate
deactivation of these units. Any sort of meddling with the @skeleton object will lead to unexpected behavior.
A triangular grid object, with class trigrid.
verticesMatrix of the vertex XYZ coordinates.
facesMatrix of the verticies forming the faces.
edgesMatrix of the vertices forming the edges.
tessellationContains the tessellation vector.
orientationContains the grid orientation in xyz 3d space, values in radian relative to the (0,1,0) direction.
centeris the xyz coordinates of the grids origin/center.
divvector contains the number of faces that a single face of the previous tessellation level is decomposed to.
faceCenterscontains the xyz coordinates of the centers of the faces on the surface of the sphere.
beltsVector of integers indicating the belt the face belongs to.
edgeLengththe length of an average edge in km and degrees.
graphan 'igraph' class graph object.
lengthinteger vector of length=3. The number of vertices, edges and faces in this order.
crsa CRS class object, by design this is the authalic sphere (ESRI:37008)
rthe radius of the grid
spThe SpatialPolygons representation of the grid. If missing, it can be created with newsp().
sfThe sf representation of the grid. If missing, it can be created with newsf().
skeletondata tables with sequential indexing for the C functions.
# single tessellation value g <- trigrid(c(8)) g # series of tessellations g1 <- trigrid(c(2,3,4)) g1 # based on approximate size (4 degrees edge length) g2 <- trigrid(deg=4)# single tessellation value g <- trigrid(c(8)) g # series of tessellations g1 <- trigrid(c(2,3,4)) g1 # based on approximate size (4 degrees edge length) g2 <- trigrid(deg=4)
trigrid objectsThe table includes basic properties of trigrids described with specific tessellation parameters
triguidetriguide
A data.frame with 120 observations and 20 variables:
totalThe total tessellation of the grid, the number of points inserted between icosahedron vertices along an edge.
level1Level 1 tessellation.
level2Level 2 tessellation - second value of the tessellation vector.
level3Level 3 tessellation - third value of the tessellation vector.
level4Level 4 tessellation - four value of the tessellation vector.
facesThe number of faces in the grid.
verticesThe number of vertices in the grid.
meanEdgeLength_degMean edge length in degrees.
sdEdgeLength_degStandard deviation of edge length in degrees.
meanEdgeLength_kmMean edge length in kilometers.
sdEdgeLength_kmStandard devation of edge length in kilometers.
meanArea_km2Mean face area in square-kilometers.
sdArea_km2Standard deviation of face area in square-kilometers.
timeTime to compute grid with an Intel Xeon E-1650 prcessor.
time_spTime to compute grid with an Intel Xeon E-1650 prcessor, with the 'sp' member.
sizeThe size of the grid in bytes.
size_spThe size of the grid object in bytes, with the 'sp' member.
timeLocate_5000Time to locate 5000 points with an Intel Xeon E-1650 in seconds (wall time).
meanSpacing_degMean spacing in degrees.
sdSpacing_degStandard deviation of spacing in degrees.
This is version 0.2, now including the spacing of the grids.
This function will return a value that is proportional to the irregularity of a triangonal face or subface. The ratio of the lengths of the shortest and the longest edges.
trishape(gridObj) ## S4 method for signature 'trigrid' trishape(gridObj) ## S4 method for signature 'hexagrid' trishape(gridObj)trishape(gridObj) ## S4 method for signature 'trigrid' trishape(gridObj) ## S4 method for signature 'hexagrid' trishape(gridObj)
gridObj |
The value is exactly 1 for an equilateral triangle, and becomes 0 as one of the edges approach 0. The values for hexagrid objects are face-specific means of subface values.
A named numeric vector, one value for every face of the grid.
g <- trigrid(3) shape <- trishape(g)g <- trigrid(3) shape <- trishape(g)
link{facelayer}).The function will get the @values slot of a facelayer object.
values(x,...) ## S4 method for signature 'gridlayer' values(x) values(x) <- value ## S4 replacement method for signature 'gridlayer,ANY' values(x) <- valuevalues(x,...) ## S4 method for signature 'gridlayer' values(x) values(x) <- value ## S4 replacement method for signature 'gridlayer,ANY' values(x) <- value
x |
( |
value |
( |
... |
Arguments passed to class-specific methods. (Not used.) |
An object of class standardGeneric of length 1.
An object of class standardGeneric of length 1.
Great circle distances between face centers and vertices
vertexradius(x, degree = TRUE)vertexradius(x, degree = TRUE)
x |
|
degree |
( |
A numeric matrix that matches in structure with the @faces slot of the provided grid x.
Distances measured on each face are in the same row.
# example grid g <- trigrid(3) # all vertexradius vertrads <- vertexradius(g) # face average averages <- apply(vertrads, 1, mean, na.rm=TRUE)# example grid g <- trigrid(3) # all vertexradius vertrads <- vertexradius(g) # face average averages <- apply(vertrads, 1, mean, na.rm=TRUE)
Shorthand function to return the vertices slot of an icosahedral grid or a grid linked to a facelayer.
vertices(x, ...) ## S4 method for signature 'trigrid' vertices(x, output = "polar") ## S4 method for signature 'facelayer' vertices(x, output = "polar")vertices(x, ...) ## S4 method for signature 'trigrid' vertices(x, output = "polar") ## S4 method for signature 'facelayer' vertices(x, output = "polar")
x |
( |
... |
Additional arguments passed to class-specific methods. |
output |
( |
a <- trigrid(1) vertices(a)a <- trigrid(1) vertices(a)
This function will return neighbouring faces of the input faces.
vicinity(gridObj, faces, ...) ## S4 method for signature 'trigrid,character' vicinity( gridObj, faces, order = 1, output = "vector", self = TRUE, namedorder = FALSE, ... )vicinity(gridObj, faces, ...) ## S4 method for signature 'trigrid,character' vicinity( gridObj, faces, order = 1, output = "vector", self = TRUE, namedorder = FALSE, ... )
gridObj |
|
faces |
( |
... |
Arguments passed to the |
order |
( |
output |
( |
self |
( |
namedorder |
( |
A character vector or a list of character vectors.
g <- trigrid(3) ne <- vicinity(g, c("F4", "F10")) neg <- trigrid(3) ne <- vicinity(g, c("F4", "F10")) ne