pycroglia.core.slimskel3d package

Subpackages

Submodules

pycroglia.core.slimskel3d.graph2skel.graph2skel(nodes, links, shape)[source]

Reconstruct a 3D skeleton (binary mask) from node/link graph. It’s not an exact recreation of the original skeleton

Parameters:
  • node (list[dict]) – List of node dicts with keys: - “idx”: list/array of voxel indices (flattened) - “links”: list of link indices connected to this node

  • link (list[dict]) – List of link dicts with key: - “point”: list/array of voxel indices (flattened) for edges

  • w (int) – Dimensions of the output volume (z, y, x)

  • l (int) – Dimensions of the output volume (z, y, x)

  • h (int) – Dimensions of the output volume (z, y, x)

  • nodes (list[Node])

  • links (list[Link])

  • shape (tuple[int, int, int])

Returns:

3D boolean skeleton array.

Return type:

np.ndarray

class pycroglia.core.slimskel3d.skel2graph.CenterOfMass(x=0, y=0, z=0)[source]

Bases: object

Represents the geometric center of mass (mean voxel position) of a node.

Parameters:
x

Mean x-coordinate of all voxels belonging to the node.

Type:

float

y

Mean y-coordinate of all voxels belonging to the node.

Type:

float

z

Mean z-coordinate of all voxels belonging to the node.

Type:

float

x: float = 0
y: float = 0
z: float = 0

Bases: object

Graph edge corresponding to a sequence of canal voxels between two nodes.

A link represents a continuous path through the skeleton, connecting two nodes (branch points or endpoints). It stores both voxel indices and references to its endpoints.

Parameters:
start

Starting node of the link.

Type:

Node

end

Ending node of the link.

Type:

Node

points

1D array of voxel indices (flattened) along the link path.

Type:

NDArray

label

Connected-component label of the skeleton the link belongs to.

Type:

int

end: Node
label: int
points: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
start: Node
class pycroglia.core.slimskel3d.skel2graph.Node(id, indices, center_of_mass, label, links=<factory>, connections=<factory>, is_endpoint=False)[source]

Bases: object

Graph node corresponding to a cluster of voxels in a skeletonized volume.

A node is created from one or more skeleton voxels that represent either a branching point or an endpoint in the skeleton.

Parameters:
id

Unique identifier of the node in the graph.

Type:

int

indices

1D array of voxel indices (flattened) belonging to the node.

Type:

NDArray

center_of_mass

Average (x, y, z) location of node voxels.

Type:

CenterOfMass

label

Connected-component label of the skeleton the node belongs to.

Type:

int

List of Link objects connecting this node to others.

Type:

list[Link]

connections

List of Node objects directly connected to this node.

Type:

list[Node]

is_endpoint

True if this node is an endpoint (degree 1), False otherwise.

Type:

bool

center_of_mass: CenterOfMass
connections: list[Node]
id: int
indices: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
is_endpoint: bool = False
label: int
links: list[Link]
pycroglia.core.slimskel3d.skel2graph.skel2graph(skeleton, threshold=0.0)[source]

“Convert a 3D skeletonized volume into a graph representation.

This function implements the logic of the MATLAB function Skel2Graph3D, producing a network graph from a binary skeleton volume. The skeleton is represented as a set of Node objects (branching points and endpoints) and Link objects (paths connecting nodes). An adjacency matrix is also returned.

Parameters:
  • skeleton (NDArray) – Binary 3D skeleton volume (z, y, x). Typically the output of skeleton3D. Foreground voxels (True/1) are treated as skeleton voxels.

  • threshold (float, optional) – Minimum length of a branch to be included. If 0, short branches that terminate in endpoints are also kept. Default is 100.0.

Returns:

Sparse adjacency matrix (n_nodes × n_nodes), where entry (i, j)

contains the number of voxels along the link between nodes i and j.

nodes (list[Node]):
List of Node objects describing graph nodes. Each node contains:
  • indices: Linear voxel indices belonging to the node.

  • center_of_mass: Average (x, y, z) position of voxels.

  • label: Connected-component label from the original skeleton.

  • links: List of Link objects connected to this node.

  • connections: Neighboring Node objects.

  • is_endpoint: True if this node is a terminal endpoint.

links (list[Link]):
List of Link objects describing graph edges. Each link contains:
  • start: Starting Node.

  • end: Ending Node.

  • points: Linear voxel indices along the path.

  • label: Connected-component label from the original skeleton.

Return type:

adjacency (csr_matrix)

Processing Steps:
  1. Pad the skeleton by 1 voxel to avoid boundary issues.

  2. Compute 26-neighborhoods for all skeleton voxels.

  3. Classify voxels as:
    • Nodes: >2 neighbors

    • Endpoints: exactly 1 neighbor

    • Canals: exactly 2 neighbors (forming links)

  4. Cluster adjacent node voxels into connected components (`Node`s).

  5. Follow canal voxels to connect nodes into `Link`s.

  6. Remove interior voxels of links to avoid re-traversal.

  7. Optionally discard links shorter than threshold unless they connect to endpoints and threshold == 0.

  8. Mark 1-degree nodes as endpoints.

  9. Build the sparse adjacency matrix.

  10. Convert voxel indices and centers of mass back to non-padded coordinates.

Notes

  • The voxel indexing is 0-based and follows NumPy’s (z, y, x) order.

  • Padding ensures correct neighborhood extraction near boundaries.

  • The graph is undirected: adjacency[i, j] = adjacency[j, i].

  • The function is designed to mirror the behavior of the original MATLAB implementation for reproducibility.

pycroglia.core.slimskel3d.skeleton3D.skeleton3D(volume)[source]

3D skeletonization of a binary volume using the Lee–Kashyap–Chu medial axis thinning algorithm.

This implementation follows the algorithm from:

Lee, T.C., Kashyap, R.L., & Chu, C.N. (1994). “Building skeleton models via 3-D medial surface/axis thinning algorithms.” Computer Vision, Graphics, and Image Processing, 56(6): 462–478.

The algorithm iteratively deletes border voxels from the foreground volume, ensuring that deletions do not alter the topology of the object. Candidate voxels are removed only if they satisfy three conditions:

  1. Not endpoints: The voxel has more than one neighbor.

  2. Euler invariant: Removal preserves the Euler characteristic, verified with a precomputed look-up table (LUT).

  3. Simple point: Removal does not break connectivity, verified using recursive octant labeling.

The thinning continues until no voxels can be removed from any of the six sweep directions.

Parameters:

volume (NDArray) – Input 3D binary volume of shape (z, y, x). Must be boolean or convertible to boolean.

Returns:

Skeletonized 3D binary volume with the same shape and dtype as the input.

Return type:

NDArray

Notes

  • A 1-voxel border of zeros is temporarily padded around the input to avoid edge effects; this is removed before returning the result.

  • Neighborhoods are evaluated in 3×3×3 cubes flattened in (z, y, x) order.

  • The Euler LUT is precomputed once via _fill_euler_lut().

  • Simple point detection is performed with _is_simple_point(), which recursively labels connected components of the neighborhood.

References

  • Lee, T.C., Kashyap, R.L., & Chu, C.N. (1994). Building skeleton models via 3-D medial surface/axis thinning algorithms. CVGIP: Image Understanding, 56(6): 462–478.

pycroglia.core.slimskel3d.slimskel3d.slimskel3d(vol, threshold)[source]

Skeletonize and iteratively slim a 3D binary image.

This function mirrors the behavior of the MATLAB SlimSkel3D: it skeletonizes a binary volume, converts it into a graph representation, prunes spurious branches, and iterates until the network length stabilizes.

Parameters:
  • cell_to_skel (np.ndarray) – 3D binary array representing the object (True/1 = foreground).

  • thr (int) – Minimum branch length (in voxels). Branches shorter than thr are removed during skeleton-to-graph conversion.

  • vol (ndarray[tuple[int, ...], dtype[_ScalarType_co]])

  • threshold (int)

Returns:

A 3D binary array of the slimmed skeleton, same shape as input.

Return type:

np.ndarray

Notes

  • Uses skeletonize_3d (equivalent to MATLAB Skeleton3D).

  • Uses skel2graph to build graph representation (nodes + links).

  • Uses graph2skel to reconstruct skeleton from graph.

  • Iterates until the total skeleton link length no longer changes.

Example

>>> import numpy as np
>>> from pycroglia.core.slimSkel3D.slim_skel import slim_skel_3d
>>> vol = np.zeros((30, 30, 30), dtype=ool)
>>> vol[5:25, 15, 15] = 1   # simple line
>>> slim = slim_skel_3d(vol, thr=5)
>>> slim.sum()  # number of skeleton voxels
20

Module contents