View source: R/mesh_components.R
fm_mesh_components | R Documentation |
Compute subsets of vertices and triangles/tetrahedrons in an fm_mesh_2d or fm_mesh_3d object that are connected by edges/triangles.
fm_mesh_components(mesh)
## S3 method for class 'fm_mesh_2d'
fm_mesh_components(mesh)
## S3 method for class 'fm_mesh_3d'
fm_mesh_components(mesh)
mesh |
An fm_mesh_2d or fm_mesh_3d object |
A list with elements vertex
and triangle
/tetra
, vectors of
integer labels for which connected component they belong, and info
, a
data.frame
with columns
component |
Connected component integer label. |
nV |
The number of vertices in the component. |
nT |
The number of triangles/tetrahedrons in the component. |
area/volume |
The surface area or volume associated with the component. Component labels are not comparable across different meshes, but some ordering stability is guaranteed by initiating each component from the lowest numbered triangle whenever a new component is initiated. |
Finn Lindgren Finn.Lindgren@gmail.com
fm_mesh_2d()
, fm_rcdt_2d()
, fm_mesh_3d()
# Construct two simple meshes:
loc <- matrix(c(0, 1, 0, 1), 2, 2)
mesh1 <- fm_mesh_2d(loc = loc, max.edge = 0.1)
bnd <- fm_nonconvex_hull(loc, 0.3)
mesh2 <- fm_mesh_2d(boundary = bnd, max.edge = 0.1)
# Compute connectivity information:
conn1 <- fm_mesh_components(mesh1)
conn2 <- fm_mesh_components(mesh2)
# One component, simply connected mesh
conn1$info
# Two disconnected components
conn2$info
# Extract the subset mesh for each component:
# (Note: some information is lost, such as fixed segments,
# and boundary edge labels.)
mesh3_1 <- fm_rcdt_2d_inla(
loc = mesh2$loc,
tv = mesh2$graph$tv[conn2$triangle == 1, , drop = FALSE],
delaunay = FALSE
)
mesh3_2 <- fm_rcdt_2d_inla(
loc = mesh2$loc,
tv = mesh2$graph$tv[conn2$triangle == 2, , drop = FALSE],
delaunay = FALSE
)
if (require("ggplot2")) {
ggplot() +
geom_fm(data = mesh3_1, fill = "red", alpha = 0.5) +
geom_fm(data = mesh3_2, fill = "blue", alpha = 0.5)
}
(m <- fm_mesh_3d(
matrix(c(1, 0, 0, 0, 1, 0, 0, 0, 1, 0, 0, 0), 4, 3, byrow = TRUE),
matrix(c(1, 2, 3, 4), 1, 4, byrow = TRUE)
))
# Compute connectivity information:
(conn <- fm_mesh_components(m))
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.