| Title: | Binary Image Thinning Algorithms |
|---|---|
| Description: | Thinning (skeletonization) algorithms for binary raster images. Provides seven algorithms behind a single dispatching function: Zhang-Suen (Zhang and Suen 1984) <doi:10.1145/357994.358023>, Guo-Hall (Guo and Hall 1989) <doi:10.1145/62065.62074>, a 2-D adaptation of Lee (Lee, Kashyap, and Chu 1994) <doi:10.1006/cgip.1994.1042>, K3M (Saeed, Tabedzki, Rybnik, and Adamski 2010) <doi:10.2478/v10006-010-0024-4>, the parallel form commonly attributed to Hilditch (1969, in 'Machine Intelligence 4'), OPTA / SPTA (Naccache and Shinghal 1984), and Holt and colleagues (1987) <doi:10.1145/12527.12531>. Also provides the medial axis transform (Blum 1967) and a distance transform implementation following Felzenszwalb and Huttenlocher (2012) <doi:10.4086/toc.2012.v008a019>. The drop-in thinImage() matches the signature of thinImage() in the 'EBImage' package on Bioconductor so existing code can switch parsers without changes. The wider thin() API selects the algorithm by name. |
| Authors: | Bill Denney [aut, cre] (ORCID: <https://orcid.org/0000-0002-5759-428X>, affiliation: Human Predictions, LLC) |
| Maintainer: | Bill Denney <[email protected]> |
| License: | LGPL-3 |
| Version: | 0.2.0 |
| Built: | 2026-05-27 22:44:39 UTC |
| Source: | https://github.com/cran/thinr |
Compute the distance from each foreground pixel to the nearest background pixel, under one of three standard metrics.
distance_transform(image, metric = c("euclidean", "manhattan", "chessboard"))distance_transform(image, metric = c("euclidean", "manhattan", "chessboard"))
image |
A binary image: a matrix where non-zero values are foreground and zero values are background. Logical, integer, and numeric inputs are accepted. |
metric |
Distance metric. One of:
|
A numeric matrix of the same shape as image. Background
pixels are 0; foreground pixels carry their distance to the
nearest background pixel.
Felzenszwalb, P. F., & Huttenlocher, D. P. (2012). Distance transforms of sampled functions. Theory of Computing, 8(19), 415-428. doi:10.4086/toc.2012.v008a019
Rosenfeld, A., & Pfaltz, J. L. (1968). Distance functions on digital pictures. Pattern Recognition, 1(1), 33-61. doi:10.1016/0031-3203(68)90013-7
# A 5x5 image with a single background pixel in the corner. m <- matrix(1L, nrow = 5, ncol = 5) m[1, 1] <- 0L distance_transform(m, metric = "manhattan") distance_transform(m, metric = "chessboard") round(distance_transform(m, metric = "euclidean"), 3)# A 5x5 image with a single background pixel in the corner. m <- matrix(1L, nrow = 5, ncol = 5) m[1, 1] <- 0L distance_transform(m, metric = "manhattan") distance_transform(m, metric = "chessboard") round(distance_transform(m, metric = "euclidean"), 3)
Return the medial axis of a binary image: the locus of foreground pixels that are local maxima of the (squared) Euclidean distance transform in at least one of the four principal directions (horizontal, vertical, NW-SE diagonal, NE-SW diagonal). Each skeleton pixel carries width information via the distance value at that point.
medial_axis(image, return_distance = FALSE)medial_axis(image, return_distance = FALSE)
image |
A binary image: a matrix where non-zero values are foreground and zero values are background. Logical, integer, and numeric inputs are accepted. |
return_distance |
Logical. If |
This is different from thin(): classical thinning algorithms
produce a connected, 1-pixel-wide skeleton without width
information. The medial axis transform (Blum 1967) produces a
skeleton with width information, useful for shape analysis
where local thickness matters.
Either a matrix (when return_distance = FALSE) or a
list(skeleton, distance) (when return_distance = TRUE).
Blum, H. (1967). A transformation for extracting new descriptors of shape. In Models for the Perception of Speech and Visual Form (pp. 362-380). MIT Press.
Felzenszwalb, P. F., & Huttenlocher, D. P. (2012). Distance transforms of sampled functions. Theory of Computing, 8(19), 415-428. doi:10.4086/toc.2012.v008a019
# A 7x9 solid rectangle: the medial axis is the middle row. m <- matrix(0L, nrow = 7, ncol = 9) m[3:5, 3:7] <- 1L medial_axis(m) # Returning width information alongside the skeleton. result <- medial_axis(m, return_distance = TRUE) result$skeleton round(result$distance, 3)# A 7x9 solid rectangle: the medial axis is the middle row. m <- matrix(0L, nrow = 7, ncol = 9) m[3:5, 3:7] <- 1L medial_axis(m) # Returning width information alongside the skeleton. result <- medial_axis(m, return_distance = TRUE) result$skeleton round(result$distance, 3)
Reduce a binary image to its one-pixel-wide skeleton using one of the supported thinning algorithms.
thin( image, method = c("zhang_suen", "guo_hall", "lee", "k3m", "hilditch", "opta", "holt"), max_iter = 1000L )thin( image, method = c("zhang_suen", "guo_hall", "lee", "k3m", "hilditch", "opta", "holt"), max_iter = 1000L )
image |
A binary image: a matrix or array where non-zero values are foreground and zero values are background. Logical, integer, and numeric inputs are all accepted. The image is treated as a 2-D matrix; arrays with more than two dimensions are not yet supported. |
method |
Algorithm to use. One of |
max_iter |
Maximum number of passes. Default 1000. Real binary images of typical sizes converge well under 50 passes; the limit is a safety bound against pathological inputs. |
A matrix of the same shape and storage mode as image, with
foreground pixels marking the thinned skeleton and the rest set to
background.
# A 3x3 solid square thins to a single foreground pixel. m <- matrix(c(0, 0, 0, 0, 0, 0, 1, 1, 1, 0, 0, 1, 1, 1, 0, 0, 1, 1, 1, 0, 0, 0, 0, 0, 0), nrow = 5, byrow = TRUE) thin(m, method = "zhang_suen") thin(m, method = "guo_hall") thin(m, method = "hilditch")# A 3x3 solid square thins to a single foreground pixel. m <- matrix(c(0, 0, 0, 0, 0, 0, 1, 1, 1, 0, 0, 1, 1, 1, 0, 0, 1, 1, 1, 0, 0, 0, 0, 0, 0), nrow = 5, byrow = TRUE) thin(m, method = "zhang_suen") thin(m, method = "guo_hall") thin(m, method = "hilditch")
EBImage::thinImage
Applies Zhang-Suen thinning to a binary image. Provided as a
signature-compatible alternative to EBImage::thinImage() so callers
can switch from EBImage to thinr by changing the namespace prefix
only.
thinImage(x)thinImage(x)
x |
A binary image. Same constraints as |
For access to the other algorithms (Guo-Hall, and eventually Lee /
K3M), use thin().
The thinned skeleton in the same storage mode as x.
m <- matrix(c(0, 1, 1, 1, 0, 0, 1, 1, 1, 0, 0, 1, 1, 1, 0), nrow = 3, byrow = TRUE) thinImage(m)m <- matrix(c(0, 1, 1, 1, 0, 0, 1, 1, 1, 0, 0, 1, 1, 1, 0), nrow = 3, byrow = TRUE) thinImage(m)