Title: | Convex Hull |
---|---|
Description: | Computes the convex hull in arbitrary dimension, based on the Qhull library (<http://www.qhull.org>). The package provides a complete description of the convex hull: edges, ridges, facets, adjacencies. Triangulation is optional. |
Authors: | C. B. Barber [cph] (author of the Qhull library), The Geometry Center [cph], Stéphane Laurent [cph, aut, cre] |
Maintainer: | Stéphane Laurent <[email protected]> |
License: | GPL-3 |
Version: | 0.7.4 |
Built: | 2024-11-16 05:35:44 UTC |
Source: | https://github.com/stla/cxhull |
Computes the convex hull of a set of points.
cxhull(points, triangulate = FALSE)
cxhull(points, triangulate = FALSE)
points |
numeric matrix, one point per row |
triangulate |
logical, whether to triangulate the convex hull |
A list providing a lot of information about the convex hull. See the README file for details.
library(cxhull) points <- rbind( c(0.5,0.5,0.5), c(0,0,0), c(0,0,1), c(0,1,0), c(0,1,1), c(1,0,0), c(1,0,1), c(1,1,0), c(1,1,1) ) cxhull(points)
library(cxhull) points <- rbind( c(0.5,0.5,0.5), c(0,0,0), c(0,0,1), c(0,1,0), c(0,1,1), c(1,0,0), c(1,0,1), c(1,1,0), c(1,1,1) ) cxhull(points)
Computes the vertices and the edges of the convex hull of a set of points.
cxhullEdges(points, adjacencies = FALSE, orderEdges = FALSE)
cxhullEdges(points, adjacencies = FALSE, orderEdges = FALSE)
points |
numeric matrix, one point per row; it must contain at least three columns (the two-dimensional case is not implemented yet) |
adjacencies |
Boolean, whether to return the vertex adjacencies |
orderEdges |
Boolean, whether to order the edges in the output |
A list with two fields: vertices
and edges
. The
vertices
field is a list which provides an id for each vertex and
its coordinates. If adjacencies=TRUE
, it provides in addition the
ids of the adjacent vertices for each vertex. The edges
fields is
an integer matrix with two columns. Each row provides the two ids of the
vertices of the corresponding edge.
library(cxhull) # let's try with the hexacosichoron (see `?hexacosichoron`) # it is convex so its convex hull is itself VE <- cxhullEdges(hexacosichoron) edges <- VE[["edges"]] random_edge <- edges[sample.int(720L, 1L), ] A <- hexacosichoron[random_edge[1L], ] B <- hexacosichoron[random_edge[2L], ] sqrt(c(crossprod(A - B))) # this is 2/phi # Now let's project the polytope to the H4 Coxeter plane phi <- (1 + sqrt(5)) / 2 u1 <- c( 0, 2*phi*sin(pi/30), 0, 1 ) u2 <- c( 2*phi*sin(pi/15), 0, 2*sin(2*pi/15), 0 ) u1 <- u1 / sqrt(c(crossprod(u1))) u2 <- u2 / sqrt(c(crossprod(u2))) # projections to the Coxeter plane proj <- function(v){ c(c(crossprod(v, u1)), c(crossprod(v, u2))) } points <- t(apply(hexacosichoron, 1L, proj)) # we will assign a color to each edge # according to the norms of its two vertices norms2 <- round(apply(points, 1L, crossprod), 1L) ( tbl <- table(norms2) ) #> 0.4 1.6 2.4 3.6 #> 30 30 30 30 values <- as.numeric(names(tbl)) grd <- as.matrix(expand.grid(values, values)) grd <- grd[grd[, 1L] <= grd[, 2L], ] pairs <- apply(grd, 1L, paste0, collapse = "-") colors <- hcl.colors(nrow(grd), palette = "Hawaii", rev = TRUE) if(require("colorspace")) { colors <- colorspace::darken(colors, amount = 0.3) } names(colors) <- pairs # plot #### opar <- par(mar = c(0, 0, 0, 0)) plot( points[!duplicated(points), ], pch = 19, cex = 0.3, asp = 1, axes = FALSE, xlab = NA, ylab = NA ) for(i in 1L:nrow(edges)){ twopoints <- points[edges[i, ], ] nrms2 <- round(sort(apply(twopoints, 1L, crossprod)), 1L) pair <- paste0(nrms2, collapse = "-") lines(twopoints, lwd = 0.5, col = colors[pair]) } par(opar)
library(cxhull) # let's try with the hexacosichoron (see `?hexacosichoron`) # it is convex so its convex hull is itself VE <- cxhullEdges(hexacosichoron) edges <- VE[["edges"]] random_edge <- edges[sample.int(720L, 1L), ] A <- hexacosichoron[random_edge[1L], ] B <- hexacosichoron[random_edge[2L], ] sqrt(c(crossprod(A - B))) # this is 2/phi # Now let's project the polytope to the H4 Coxeter plane phi <- (1 + sqrt(5)) / 2 u1 <- c( 0, 2*phi*sin(pi/30), 0, 1 ) u2 <- c( 2*phi*sin(pi/15), 0, 2*sin(2*pi/15), 0 ) u1 <- u1 / sqrt(c(crossprod(u1))) u2 <- u2 / sqrt(c(crossprod(u2))) # projections to the Coxeter plane proj <- function(v){ c(c(crossprod(v, u1)), c(crossprod(v, u2))) } points <- t(apply(hexacosichoron, 1L, proj)) # we will assign a color to each edge # according to the norms of its two vertices norms2 <- round(apply(points, 1L, crossprod), 1L) ( tbl <- table(norms2) ) #> 0.4 1.6 2.4 3.6 #> 30 30 30 30 values <- as.numeric(names(tbl)) grd <- as.matrix(expand.grid(values, values)) grd <- grd[grd[, 1L] <= grd[, 2L], ] pairs <- apply(grd, 1L, paste0, collapse = "-") colors <- hcl.colors(nrow(grd), palette = "Hawaii", rev = TRUE) if(require("colorspace")) { colors <- colorspace::darken(colors, amount = 0.3) } names(colors) <- pairs # plot #### opar <- par(mar = c(0, 0, 0, 0)) plot( points[!duplicated(points), ], pch = 19, cex = 0.3, asp = 1, axes = FALSE, xlab = NA, ylab = NA ) for(i in 1L:nrow(edges)){ twopoints <- points[edges[i, ], ] nrms2 <- round(sort(apply(twopoints, 1L, crossprod)), 1L) pair <- paste0(nrms2, collapse = "-") lines(twopoints, lwd = 0.5, col = colors[pair]) } par(opar)
A matrix giving the 62 vertices of da Vinci's 72-sided sphere, a convex polyhedra with 72 faces.
daVinciSphere
daVinciSphere
A matrix with 62 rows and 3 columns.
http://www.matematicasvisuales.com/english/html/geometry/space/sphereCampanus.html
Dihedral angles of a convex hull.
dihedralAngles(hull)
dihedralAngles(hull)
hull |
an output of |
A dataframe with three columns. The two first columns represent
the edges, given as a pair of vertex indices. The third column provides
the dihedral angle in degrees corresponding to the edge, that is the
angle between the two faces incident to this edge. This is useful to find
edges between two coplanar faces: if the faces are exactly coplanar then
the dihedral angle is 180, but because of numerical approximation one can
consider that there is coplanarity when the dihedral angle is greater
than 179, for example. This function is used in
plotConvexHull3d
to get rid of such edges (if the user
sets a value to the argument angleThreshold
).
# a cube #### library(cxhull) points <- rbind( c(0.5,0.5,0.5), c(0,0,0), c(0,0,1), c(0,1,0), c(0,1,1), c(1,0,0), c(1,0,1), c(1,1,0), c(1,1,1) ) hull <- cxhull(points) dihedralAngles(hull)
# a cube #### library(cxhull) points <- rbind( c(0.5,0.5,0.5), c(0,0,0), c(0,0,1), c(0,1,0), c(0,1,1), c(1,0,0), c(1,0,1), c(1,1,0), c(1,1,1) ) hull <- cxhull(points) dihedralAngles(hull)
Edges of a triangulated 3D convex hull given by the ids of the vertices in a matrix, plus a column indicating the border edges.
EdgesAB(hull)
EdgesAB(hull)
hull |
an output of |
A character matrix with three columns. Each row provides the ids of the two vertices of an edge, and a yes/no indicator of whether the edge is a border edge.
library(cxhull) library(rgl) dodecahedron <- t(dodecahedron3d()$vb[-4L, ]) hull <- cxhull(dodecahedron, triangulate = TRUE) triangles <- TrianglesXYZ(hull) triangles3d(triangles, color = "yellow") edges <- EdgesAB(hull) trueEdges <- edges[edges[, 3L] == "yes", c(1L, 2L)] otherEdges <- edges[edges[, 3L] == "no", c(1L, 2L)] vertices <- VerticesXYZ(hull) for(i in 1:nrow(trueEdges)){ lines3d(vertices[trueEdges[i, ], ], color = "blue", lwd = 3) } for(i in 1:nrow(otherEdges)){ lines3d(vertices[otherEdges[i, ], ], color = "red", lwd = 3) }
library(cxhull) library(rgl) dodecahedron <- t(dodecahedron3d()$vb[-4L, ]) hull <- cxhull(dodecahedron, triangulate = TRUE) triangles <- TrianglesXYZ(hull) triangles3d(triangles, color = "yellow") edges <- EdgesAB(hull) trueEdges <- edges[edges[, 3L] == "yes", c(1L, 2L)] otherEdges <- edges[edges[, 3L] == "no", c(1L, 2L)] vertices <- VerticesXYZ(hull) for(i in 1:nrow(trueEdges)){ lines3d(vertices[trueEdges[i, ], ], color = "blue", lwd = 3) } for(i in 1:nrow(otherEdges)){ lines3d(vertices[otherEdges[i, ], ], color = "red", lwd = 3) }
The coordinates of the extremities of the edges in a matrix, plus a column indicating which edges are border edges.
EdgesXYZ(hull)
EdgesXYZ(hull)
hull |
an output of |
A numeric matrix with four columns. The first three values of a row are the coordinates of a vertex at the extremity of an edge, and the fourth column indicates whether the edge is a border edge.
A matrix giving the 120 vertices of the hexacosichoron, a
regular convex 4D polytope also known as the "600-cell", with edge length
2/phi
, where phi
is the golden number. It has 720 edges.
hexacosichoron
hexacosichoron
A matrix with 120 rows and 4 columns.
https://www.qfbox.info/4d/600-cell
Extract the vertices and the faces from a 3d convex hull.
hullMesh(hull, simplify = TRUE, rgl = FALSE)
hullMesh(hull, simplify = TRUE, rgl = FALSE)
hull |
a 3d convex hull, output of |
simplify |
Boolean, whether to return the faces as a matrix instead
of a list if possible, i.e. if all faces have the same number of edges;
this argument is possibly ignored if |
rgl |
Boolean, whether to return a rgl mesh (class
|
A list giving the vertices and the faces, or a rgl mesh.
Unless all faces are triangular, the output does not define a mesh with coherently oriented faces.
library(cxhull) hull <- cxhull(daVinciSphere) septuaginta <- hullMesh(hull, rgl = TRUE) library(rgl) open3d(windowRect = c(50, 50, 562, 562)) shade3d(septuaginta, color = "darkred") # some quad faces are misoriented: open3d(windowRect = c(50, 50, 562, 562)) shade3d(septuaginta, color = "tomato", back = "culled")
library(cxhull) hull <- cxhull(daVinciSphere) septuaginta <- hullMesh(hull, rgl = TRUE) library(rgl) open3d(windowRect = c(50, 50, 562, 562)) shade3d(septuaginta, color = "darkred") # some quad faces are misoriented: open3d(windowRect = c(50, 50, 562, 562)) shade3d(septuaginta, color = "tomato", back = "culled")
Summary of a triangulated 3D convex hull
hullSummary(hull)
hullSummary(hull)
hull |
an output of |
A list with the vertices and the facets.
library(cxhull) # pyramid pts <- rbind( c(0, 0, 0), c(1, 0, 0), c(1, 1, 0), c(0.5, 0.5, 1), c(0.5, 0.5, 0.9), c(0, 1, 0) ) hull <- cxhull(pts, triangulate = TRUE) hullSummary(hull)
library(cxhull) # pyramid pts <- rbind( c(0, 0, 0), c(1, 0, 0), c(1, 1, 0), c(0.5, 0.5, 1), c(0.5, 0.5, 0.9), c(0, 1, 0) ) hull <- cxhull(pts, triangulate = TRUE) hullSummary(hull)
Plot a triangulated 3d convex hull with rgl.
plotConvexHull3d( hull, angleThreshold = NULL, edgesAsTubes = TRUE, verticesAsSpheres = TRUE, palette = NULL, bias = 1, interpolate = "linear", g = identity, facesColor = "navy", edgesColor = "gold", tubesRadius = 0.03, spheresRadius = 0.05, spheresColor = edgesColor, alpha = 1 )
plotConvexHull3d( hull, angleThreshold = NULL, edgesAsTubes = TRUE, verticesAsSpheres = TRUE, palette = NULL, bias = 1, interpolate = "linear", g = identity, facesColor = "navy", edgesColor = "gold", tubesRadius = 0.03, spheresRadius = 0.05, spheresColor = edgesColor, alpha = 1 )
hull |
an output of |
angleThreshold |
a threshold angle in degrees, typically |
edgesAsTubes |
Boolean, whether to draw the edges as tubes |
verticesAsSpheres |
Boolean, whether to draw the vertices as spheres |
palette |
a vector of colors to make a color gradient for the faces;
if |
bias , interpolate
|
if |
g |
a function defined on [0, 1] and taking its values in [0, 1]; it is
composed with the function created by |
facesColor |
the color(s) for the faces; this argument is ignored if
the argument |
edgesColor |
the color for the edges |
tubesRadius |
the radius of the tubes when |
spheresRadius |
the radius of the spheres when
|
spheresColor |
the color of the spheres when
|
alpha |
number between 0 and 1 controlling the opacity of the faces |
No value.
library(cxhull) library(rgl) cuboctahedron <- t(cuboctahedron3d()$vb[-4L, ]) hull <- cxhull(cuboctahedron, triangulate = TRUE) # single color #### open3d(windowRect = c(50, 50, 562, 562)) plotConvexHull3d(hull) # gradient #### open3d(windowRect = c(50, 50, 562, 562)) if(getRversion() < "4.1.0"){ palette <- "Viridis" }else{ palette <- "Rocket" } plotConvexHull3d(hull, palette = hcl.colors(256, palette), bias = 0.5) library(cxhull) library(rgl) # Leonardo da Vinci's 72-sided sphere #### hull <- cxhull(daVinciSphere, triangulate = TRUE) # there are some undesirable edges: plotConvexHull3d( hull, tubesRadius = 0.07, spheresRadius = 0.1 ) # => use `angleThreshold` to get rid of these edges: plotConvexHull3d( hull, angleThreshold = 179, tubesRadius = 0.07, spheresRadius = 0.1 )
library(cxhull) library(rgl) cuboctahedron <- t(cuboctahedron3d()$vb[-4L, ]) hull <- cxhull(cuboctahedron, triangulate = TRUE) # single color #### open3d(windowRect = c(50, 50, 562, 562)) plotConvexHull3d(hull) # gradient #### open3d(windowRect = c(50, 50, 562, 562)) if(getRversion() < "4.1.0"){ palette <- "Viridis" }else{ palette <- "Rocket" } plotConvexHull3d(hull, palette = hcl.colors(256, palette), bias = 0.5) library(cxhull) library(rgl) # Leonardo da Vinci's 72-sided sphere #### hull <- cxhull(daVinciSphere, triangulate = TRUE) # there are some undesirable edges: plotConvexHull3d( hull, tubesRadius = 0.07, spheresRadius = 0.1 ) # => use `angleThreshold` to get rid of these edges: plotConvexHull3d( hull, angleThreshold = 179, tubesRadius = 0.07, spheresRadius = 0.1 )
Coordinates of the vertices of the triangles of a triangulated 3D convex hull.
TrianglesXYZ(hull)
TrianglesXYZ(hull)
hull |
an output of |
A matrix with three columns. Each row represents the coordinates of a vertex of a triangle.
library(cxhull) library(rgl) dodecahedron <- t(dodecahedron3d()$vb[-4L, ]) hull <- cxhull(dodecahedron, triangulate = TRUE) triangles <- TrianglesXYZ(hull) triangles3d(triangles, color = "firebrick")
library(cxhull) library(rgl) dodecahedron <- t(dodecahedron3d()$vb[-4L, ]) hull <- cxhull(dodecahedron, triangulate = TRUE) triangles <- TrianglesXYZ(hull) triangles3d(triangles, color = "firebrick")
The coordinates of the vertices of a 3D convex hull.
VerticesXYZ(hull)
VerticesXYZ(hull)
hull |
an output of |
A matrix with three columns. Each row represents the coordinates of a vertex and the row names are the ids of the vertices.