Title: | Implementation of the 3D Alpha-Shape for the Reconstruction of 3D Sets from a Point Cloud |
---|---|
Description: | Implementation in R of the alpha-shape of a finite set of points in the three-dimensional space. The alpha-shape generalizes the convex hull and allows to recover the shape of non-convex and even non-connected sets in 3D, given a random sample of points taken into it. Besides the computation of the alpha-shape, this package provides users with functions to compute the volume of the alpha-shape, identify the connected components and facilitate the three-dimensional graphical visualization of the estimated set. |
Authors: | Thomas Lafarge, Beatriz Pateiro-Lopez |
Maintainer: | Beatriz Pateiro-Lopez <[email protected]> |
License: | GPL-2 |
Version: | 1.3.2 |
Built: | 2024-11-11 07:06:21 UTC |
Source: | CRAN |
-shape computationThis function calculates the 3D -shape of a given sample of
points in the three-dimensional space for
.
ashape3d(x, alpha, pert = FALSE, eps = 1e-09)
ashape3d(x, alpha, pert = FALSE, eps = 1e-09)
x |
A 3-column matrix with the coordinates of the input points.
Alternatively, an object of class |
alpha |
A single value or vector of values for |
pert |
Logical. If the input points are not in general position and
|
eps |
Scaling factor used for data perturbation when the input points are not in general position, see Details. |
If x
is an object of class "ashape3d"
, then ashape3d
does not recompute the 3D Delaunay triangulation (it reduces the
computational cost).
If the input points x
are not in general position and pert
is
set to TRUE, the function adds random noise to the data. The noise is
generated from a normal distribution with mean zero and standard deviation
eps*sd(x)
.
An object of class "ashape3d"
with the following components
(see Edelsbrunner and Mucke (1994) for notation):
tetra |
For each
tetrahedron of the 3D Delaunay triangulation, the matrix |
triang |
For each triangle of the 3D Delaunay
triangulation, the matrix |
edge |
For each edge of the 3D Delaunay
triangulation, the matrix |
vertex |
For each sample point, the matrix |
x |
A 3-column matrix with the coordinates of the original sample points. |
alpha |
A single value or vector of values of |
xpert |
A 3-column matrix with the coordinates of the perturbated
sample points (only when the input points are not in general position and
|
Edelsbrunner, H., Mucke, E. P. (1994). Three-Dimensional Alpha Shapes. ACM Transactions on Graphics, 13(1), pp.43-72.
T1 <- rtorus(1000, 0.5, 2) T2 <- rtorus(1000, 0.5, 2, ct = c(2, 0, 0), rotx = pi/2) x <- rbind(T1, T2) # Value of alpha alpha <- 0.25 # 3D alpha-shape ashape3d.obj <- ashape3d(x, alpha = alpha) plot(ashape3d.obj) # For new values of alpha, we can use ashape3d.obj as input (faster) alpha <- c(0.15, 1) ashape3d.obj <- ashape3d(ashape3d.obj, alpha = alpha) plot(ashape3d.obj, indexAlpha = 2:3)
T1 <- rtorus(1000, 0.5, 2) T2 <- rtorus(1000, 0.5, 2, ct = c(2, 0, 0), rotx = pi/2) x <- rbind(T1, T2) # Value of alpha alpha <- 0.25 # 3D alpha-shape ashape3d.obj <- ashape3d(x, alpha = alpha) plot(ashape3d.obj) # For new values of alpha, we can use ashape3d.obj as input (faster) alpha <- c(0.15, 1) ashape3d.obj <- ashape3d(ashape3d.obj, alpha = alpha) plot(ashape3d.obj, indexAlpha = 2:3)
This function calculates and clusters the different connected components of
the -shape of a given sample of points in the three-dimensional
space.
components_ashape3d(as3d, indexAlpha = 1)
components_ashape3d(as3d, indexAlpha = 1)
as3d |
An object of class |
indexAlpha |
A single value or vector with the indexes of
|
The function components_ashape3d
computes the connected components of
the -shape for each value of
in
as3d$alpha[indexAlpha]
when indexAlpha
is numeric.
If indexAlpha="all"
or indexAlpha="ALL"
then the function
computes the connected components of the -shape for all values
of
in
as3d$alpha
.
If indexAlpha
is a single value then the function returns a
vector v
of length equal to the sample size. For each sample point
i
, v[i]
represents the label of the connected component to
which the point belongs (for isolated points, v[i]=-1
). The labels of
the connected components are ordered by size where the largest one (in
number of vertices) gets the smallest label which is one.
Otherwise components_ashape3d
returns a list of vectors describing
the connected components of the -shape for each selected value
of
.
T1 <- rtorus(1000, 0.5, 2) T2 <- rtorus(1000, 0.5, 2, ct = c(2, 0, 0), rotx = pi/2) x <- rbind(T1, T2) alpha <- c(0.25, 2) ashape3d.obj <- ashape3d(x, alpha = alpha) plot(ashape3d.obj, indexAlpha = "all") # Connected components of the alpha-shape for both values of alpha comp <- components_ashape3d(ashape3d.obj, indexAlpha = "all") class(comp) # Number of components and points in each component for alpha=0.25 table(comp[[1]]) # Number of components and points in each component for alpha=2 table(comp[[2]]) # Plot the connected components for alpha=0.25 plot(ashape3d.obj, byComponents = TRUE, indexAlpha = 1)
T1 <- rtorus(1000, 0.5, 2) T2 <- rtorus(1000, 0.5, 2, ct = c(2, 0, 0), rotx = pi/2) x <- rbind(T1, T2) alpha <- c(0.25, 2) ashape3d.obj <- ashape3d(x, alpha = alpha) plot(ashape3d.obj, indexAlpha = "all") # Connected components of the alpha-shape for both values of alpha comp <- components_ashape3d(ashape3d.obj, indexAlpha = "all") class(comp) # Number of components and points in each component for alpha=0.25 table(comp[[1]]) # Number of components and points in each component for alpha=2 table(comp[[2]]) # Plot the connected components for alpha=0.25 plot(ashape3d.obj, byComponents = TRUE, indexAlpha = 1)
-shapeThis function checks whether points are inside an -shape.
inashape3d(as3d, indexAlpha = 1, points)
inashape3d(as3d, indexAlpha = 1, points)
as3d |
An object of class |
indexAlpha |
A single value or vector with the indexes of
|
points |
A 3-column matrix with the coordinates of the input points. |
The function inashape3d
checks whether each point in points
is
inside the -shape for each value of
in
as3d$alpha[indexAlpha]
.
If indexAlpha="all"
or indexAlpha="ALL"
then the function
checks whether each point in points
is inside the -shape
for all values of
in
as3d$alpha
.
If indexAlpha
is a single value then the function returns a
vector of boolean of length the number of input points. The element at
position i
is TRUE
if the point in points[i,]
is inside
the -shape.
Otherwise inashape3d
returns a list of vectors of boolean values
(each object in the list as described above).
T1 <- rtorus(2000, 0.5, 2) T2 <- rtorus(2000, 0.5, 2, ct = c(2, 0, 0), rotx = pi/2) x <- rbind(T1, T2) ashape3d.obj <- ashape3d(x, alpha = 0.4) # Random sample of points in a plane points <- matrix(c(5*runif(10000) - 2.5, rep(0.01, 5000)), nc = 3) in3d <- inashape3d(ashape3d.obj, points = points) plot(ashape3d.obj, transparency = 0.2) colors <- ifelse(in3d, "blue", "green") points3d(points, col = colors)
T1 <- rtorus(2000, 0.5, 2) T2 <- rtorus(2000, 0.5, 2, ct = c(2, 0, 0), rotx = pi/2) x <- rbind(T1, T2) ashape3d.obj <- ashape3d(x, alpha = 0.4) # Random sample of points in a plane points <- matrix(c(5*runif(10000) - 2.5, rep(0.01, 5000)), nc = 3) in3d <- inashape3d(ashape3d.obj, points = points) plot(ashape3d.obj, transparency = 0.2) colors <- ifelse(in3d, "blue", "green") points3d(points, col = colors)
-shape in 3DThis function plots the -shape in 3D using the package
rgl
.
## S3 method for class 'ashape3d' plot(x, clear = TRUE, col = c(2, 2, 2), byComponents = FALSE, indexAlpha = 1, transparency = 1, walpha = FALSE, triangles = TRUE, edges = TRUE, vertices = TRUE, ...)
## S3 method for class 'ashape3d' plot(x, clear = TRUE, col = c(2, 2, 2), byComponents = FALSE, indexAlpha = 1, transparency = 1, walpha = FALSE, triangles = TRUE, edges = TRUE, vertices = TRUE, ...)
x |
An object of class |
clear |
Logical, specifying whether the current rgl device should be cleared. |
col |
A vector of length three specifying the colors of the triangles,
edges and vertices composing the |
byComponents |
Logical, if TRUE the connected components of the
|
indexAlpha |
A single value or vector with the indexes of
|
transparency |
The coefficient of transparency, from 0 (transparent) to
1 (opaque), used to plot the |
walpha |
Logical, if TRUE the value of |
triangles |
Logical, if TRUE triangles are plotted. |
edges |
Logical, if TRUE edges are plotted. |
vertices |
Logical, if TRUE vertices are plotted. |
... |
Material properties. See |
The function plot.ashape3d
opens a rgl device for each value of
in
x$alpha[indexAlpha]
. Device information is displayed
in the console.
If indexAlpha="all"
or indexAlpha="ALL"
then the function
represents the -shape for all values of
in
as3d$alpha
.
T1 <- rtorus(1000, 0.5, 2) T2 <- rtorus(1000, 0.5, 2, ct = c(2, 0, 0), rotx = pi/2) x <- rbind(T1, T2) alpha <- c(0.15, 0.25, 1) ashape3d.obj <- ashape3d(x, alpha = alpha) # Plot the alpha-shape for all values of alpha plot(ashape3d.obj, indexAlpha = "all") # Plot the connected components of the alpha-shape for alpha=0.25 plot(ashape3d.obj, byComponents = TRUE, indexAlpha = 2)
T1 <- rtorus(1000, 0.5, 2) T2 <- rtorus(1000, 0.5, 2, ct = c(2, 0, 0), rotx = pi/2) x <- rbind(T1, T2) alpha <- c(0.15, 0.25, 1) ashape3d.obj <- ashape3d(x, alpha = alpha) # Plot the alpha-shape for all values of alpha plot(ashape3d.obj, indexAlpha = "all") # Plot the connected components of the alpha-shape for alpha=0.25 plot(ashape3d.obj, byComponents = TRUE, indexAlpha = 2)
This function generates random points within the torus whose minor
radius is
, major radius is
and center is
.
rtorus(n, r, R, ct = c(0, 0, 0), rotx = NULL)
rtorus(n, r, R, ct = c(0, 0, 0), rotx = NULL)
n |
Number of observations. |
r |
Minor radius (radius of the tube). |
R |
Major radius (distance from the center of the tube to the center of the torus). |
ct |
A vector with the coordinates of the center of the torus. |
rotx |
If not NULL, a rotation through an angle |
T1 <- rtorus(2000, 0.5, 2.5) bbox3d(color = c("white", "black")) points3d(T1, col = 4) T2 <- rtorus(2000, 0.5, 2.5, ct = c(2, 0, 0.5), rotx = pi/2) points3d(T2, col = 2)
T1 <- rtorus(2000, 0.5, 2.5) bbox3d(color = c("white", "black")) points3d(T1, col = 4) T2 <- rtorus(2000, 0.5, 2.5, ct = c(2, 0, 0.5), rotx = pi/2) points3d(T2, col = 2)
This function calculates the normal vectors of all the triangles which
belong to the boundary of the -shape.
surfaceNormals(x, indexAlpha = 1, display = FALSE, col = 3, scale = 1, ...)
surfaceNormals(x, indexAlpha = 1, display = FALSE, col = 3, scale = 1, ...)
x |
An object of class |
indexAlpha |
A single value or vector with the indexes of
|
display |
Logical, if TRUE, |
col |
Color of the normal vectors. |
scale |
Scale parameter to control the length of the surface normals, only affect display. |
... |
Material properties. See |
The function surfaceNormals
computes the normal vectors of all the
triangles which belong to the boundary of the -shape for each
value of
in
x$alpha[indexAlpha]
. The magnitude of each
vector is equal to the area of its associated triangle.
If indexAlpha="all"
or indexAlpha="ALL"
then the function
computes the normal vectors for all values of in
as3d$alpha
.
If indexAlpha
is a single value then the function returns an
object of class "normals"
with the following components:
normals |
Three-column matrix with the euclidean coordinates of the normal vectors. |
centers |
Three-column matrix with the euclidean
coordinates of the centers of the triangles that form the
|
Otherwise surfaceNormals
returns a list of class
"normals-List"
(each object in the list as described above).
x <- rtorus(1000, 0.5, 1) alpha <- 0.25 ashape3d.obj <- ashape3d(x, alpha = alpha) surfaceNormals(ashape3d.obj, display = TRUE)
x <- rtorus(1000, 0.5, 1) alpha <- 0.25 ashape3d.obj <- ashape3d(x, alpha = alpha) surfaceNormals(ashape3d.obj, display = TRUE)
This function calculates the volume of the -shape of a given
sample of points in the three-dimensional space.
volume_ashape3d(as3d, byComponents = FALSE, indexAlpha = 1)
volume_ashape3d(as3d, byComponents = FALSE, indexAlpha = 1)
as3d |
An object of class |
byComponents |
Logical, if FALSE (default) |
indexAlpha |
A single value or vector with the indexes of
|
The function volume_ashape3d
computes the volume of the
-shape for each value of
in
as3d$alpha[indexAlpha]
when indexAlpha
is numeric.
If indexAlpha="all"
or indexAlpha="ALL"
then the function
computes the volume of the -shape for all values of
in
as3d$alpha
.
If indexAlpha
is a single value then the function returns the
volume of the -shape (if the argument
byComponents
is set
to FALSE) or a vector with the volumes of each connected component of the
-shape (if the argument
byComponents
is set to TRUE).
Otherwise volume_ashape3d
returns a list (each object in the list as
described above).
C1 <- matrix(runif(6000), nc = 3) C2 <- matrix(runif(6000), nc = 3) + 2 x <- rbind(C1, C2) ashape3d.obj <- ashape3d(x, alpha = 0.75) plot(ashape3d.obj, byComponents = TRUE) # Compute the volume of the alpha-shape volume_ashape3d(ashape3d.obj) # Compute the volumes of the connected components of the alpha-shape volume_ashape3d(ashape3d.obj, byComponents = TRUE)
C1 <- matrix(runif(6000), nc = 3) C2 <- matrix(runif(6000), nc = 3) + 2 x <- rbind(C1, C2) ashape3d.obj <- ashape3d(x, alpha = 0.75) plot(ashape3d.obj, byComponents = TRUE) # Compute the volume of the alpha-shape volume_ashape3d(ashape3d.obj) # Compute the volumes of the connected components of the alpha-shape volume_ashape3d(ashape3d.obj, byComponents = TRUE)