pointtree.instance_segmentation

Algorithms for tree instance segmentation.

class pointtree.instance_segmentation.InstanceSegmentationAlgorithm(
trunk_class_id: int,
crown_class_id: int,
branch_class_id: int | None = None,
downsampling_voxel_size: float | None = None,
)[source]

Bases: ABC

Base class for implementing tree instance segmentation algorithms.

Parameters:
  • trunk_class_id – Integer class ID that designates the tree trunk points.

  • crown_class_id – Integer class ID that designates the tree crown points.

  • branch_class_id – Integer class ID that designates the tree branch points. Defaults to None, assuming that branch points are not separately labeled.

downsampling_voxel_size: Voxel size for the voxel-based downsampling of the tree points before performing the

tree instance segmentation. Defaults to None, which means that the tree instance segmentation is performed with the full resolution of the point cloud.

runtime_stats() DataFrame[source]
Returns:

Tracked execution times as pandas.DataFrame with the columns "Description" and "Runtime".

static upsample_instance_ids(
original_coords: ndarray,
downsampled_coords: ndarray,
instance_ids: ndarray,
non_predicted_point_indices: ndarray,
predicted_point_indices: ndarray,
) ndarray[source]

Upsamples instance segmentation labels to a higher point cloud resolution.

Parameters:
  • original_coords – Coordinates of the points from the higher-resolution point cloud

  • downsampled_coords – Coordinates of the points from the lower-resolution point cloud.

  • instance_ids – Instance segmentation labels to be upsampled.

  • non_predicted_point_indices – Indices of the points from the higher-resolution point cloud that were discarded during downsampling.

  • predicted_point_indices – Indices of the points from the higher-resolution point cloud that were kept during downsampling.

Returns:

Upsampled instance segmentation labels.

Shape:
  • original_coords: \((N, 3)\)

  • downsampled_coords: \((N', 3)\)

  • instance_ids: \((N')\)

  • instance_ids: \((N - N')\)

  • instance_ids: \((N')\)

  • Output: \((N)\).

where

\(N = \text{ number of points before in the higher-resolution point cloud}\)
\(N' = \text{ number of points before in the lower-resolution point cloud}\)
class pointtree.instance_segmentation.MultiStageAlgorithm(
trunk_class_id: int,
crown_class_id: int,
*,
branch_class_id: int | None = None,
algorithm: Literal['full', 'watershed_crown_top_positions', 'watershed_matched_tree_positions'] = 'full',
downsampling_voxel_size: float | None = None,
visualization_folder: str | None = None,
eps_trunk_clustering: float = 2.5,
min_samples_trunk_clustering: int = 1,
min_trunk_points: int = 100,
grid_size_canopy_height_model: float = 0.5,
min_distance_crown_tops: float = 7,
min_points_crown_detection: float = 100,
min_tree_height: float = 2.5,
smooth_canopy_height_model: bool = True,
smooth_sigma: float = 1,
distance_match_trunk_and_crown_top: float = 5,
correct_watershed: bool = True,
max_point_spacing_region_growing: float = 0.08,
max_radius_region_growing: float = 1,
multiplier_outside_coarse_border: float = 2,
num_neighbors_region_growing: int = 27,
z_scale: float = 0.5,
)[source]

Bases: InstanceSegmentationAlgorithm

Multi-stage algorithm for tree instance segmentation.

Parameters:
  • trunk_class_id – Integer class ID that designates the tree trunk points.

  • crown_class_id – Integer class ID that designates the tree crown points.

  • branch_class_id – Integer class ID that designates the tree branch points. Defaults to None, assuming that branch points are not separately labeled. If branches are separately labeled, the branch points are treated as trunk points by this algorithm.

  • algorithm – Variant of the algorithm to be used: "full": The full algorithm is used. "watershed_crown_top_positions": A marker-controlled Watershed segmentation is performed, using the crown top positions as markers. "watershed_matched_tree_positions": A marker-controlled Watershed segmentation is performed, using the tree positions as markers, resulting from the matching of crown top positions and trunk positions.

  • downsampling_voxel_size – Voxel size for the voxel-based downsampling of the tree points before performing the tree instance segmentation. Defaults to None, which means that the tree instance segmentation is performed with the full resolution of the point cloud.

  • visualization_folder – Path of a directory in which to store visualizations of intermediate results of the algorithm. Defaults to None, which means that no visualizations are stored.

Parameters for the DBSCAN clustering of trunk points:
  • eps_trunk_clustering – Parameter \(\epsilon\) for clustering the trunk points using the DBSCAN algorithm. Further details on the meaning of the parameter can be found in the documentation of sklearn.cluster.DBSCAN. Defaults to 2.5 m.

  • min_samples_trunk_clustering – Parameter \(MinSamples\) for clustering the trunk points using the DBSCAN algorithm. Further details on the meaning of the parameter can be found in the documentation of sklearn.cluster.DBSCAN. Defaults to 1.

  • min_trunk_points – Minimum number of points a cluster of trunk points must contain to be considered as a trunk. Defaults to 100.

Parameters for the construction and maximum filtering of the canopy height model:
  • grid_size_canopy_height_model – Width of the 2D grid cells used to create the canopy height model. Defaults to 0.5 m.

  • min_distance_crown_tops – Minimum horizontal distance that local maxima of the canopy height model must have in order to be considered as separate crown tops. Defaults to 7 m.

  • min_points_crown_detection – Minimum number of points that must be contained in a cell of the canopy height model for the cell to be considered a possible crown top. Defaults to 100.

  • min_tree_height – Minimum height that a local maximum of the canopy height model must have to be considered as a crown top. Defaults to 2.5 m.

  • smooth_canopy_height_model – Whether to smooth the canopy height model using a Gaussian blur filter. Defaults to True.

  • smooth_sigma – Parameter \(\sigma\) of the Gaussian blur filter used to smooth the canopy height model. Defaults to 1.

Parameters for the matching of trunk positions and crown top positions:

distance_match_trunk_and_crown_top – Maximum horizontal distance between trunk positions and crown top positions up to which both are considered to belong to the same tree. Defaults to 5 m.

Parameters for the Watershed segmentation:

correct_watershed – Whether erroneous parts of the watershed segmentation should be replaced by a Voronoi segmentation. Defaults to True.

Parameters for the region growing segmentation:
  • max_point_spacing_region_growing – The results of the Watershed segmentation are only refined in sufficiently dense point cloud regions. For this purpose, the average distance of the points to their nearest neighbor is determined for each tree. If this average distance is less than max_point_spacing_region_growing the tree is considered for refining its segmentation using the region growing approach. Defaults to 0.08 m.

  • max_radius_region_growing – Maximum radius in which to search for neighboring points during region growing. Defaults to 1m .

  • multiplier_outside_coarse_border – In our region growing approach, the points are processed in a sorted order, with points with the smallest distance to an already assigned tree point being processed first. For points that lie outside the crown boundary, the distance is multiplied by a constant factor to restrict growth in areas outside the crown boundary. This parameter defines this constant factor. The larger the factor, the more growth is restricted in areas outside the crown boundaries. Defaults to 2.

  • num_neighbors_region_growing – In our region growing approach, the k-nearest neighbors of each seed point a are searched and may be added to the same tree instance. This parameter specifies the number of neighbors to search. Defaults to 27.

  • z_scale – Factor by which the z-coordinates are multiplied in region growing. Using a value between zero and one favors upward growth. Defaults to 0.5.

cluster_trunks(
tree_coords: ndarray,
classification: ndarray,
) Tuple[ndarray, ndarray][source]

Clusters tree trunk points using the DBSCAN algorithm and assigns the same unique instance ID to the points belonging to the same cluster.

Parameters:
  • tree_coords – Coordinates of all tree points.

  • classification – Class ID of all tree points.

Returns:

Tuple of two arrays. The first contains the instance ID of each tree point. Points that do not belong to any

instance are assigned the ID \(-1\). The second contains the unique instance IDs.

Shape:
  • tree_coords: \((N, 3)\)

  • classification: \((N)\)

  • Output: Tuple of two arrays. The first has shape \((N)\) and the second \((T)\)

where

\(N = \text{ number of tree points}\)
\(T = \text{ number of trunk instances}\)
coarse_segmentation(
tree_coords: ndarray,
instance_ids: ndarray,
tree_positions_grid: ndarray,
canopy_height_model: ndarray,
grid_origin: ndarray,
*,
point_cloud_id: str | None = None,
trunk_positions_grid: ndarray | None = None,
) Tuple[ndarray, ndarray, ndarray, ndarray][source]

Coarse tree instance segmentation using the marker-controlled Watershed algorithm.

Parameters:
  • tree_coords – Coordinates of all tree points.

  • instance_ids – Instance IDs of each tree point.

  • tree_positions_grid – The 2D positions of all tree instances as in integer coordinates in the grid coordinate system of the canopy height model.

  • canopy_height_model – The canopy height model to segment.

  • grid_origin – 2D coordinates of the origin of the canopy height model.

  • point_cloud_id – ID of the point cloud to be used in the file names of the visualization results. Defaults to None, which means that no visualizations are created.

  • trunk_positions_grid – The 2D positions of tree trunks that were not matched with a tree crown as in integer coordinates in the grid coordinate system of the canopy height model. Only used for visualization purposes.

Returns:

Tuple of four arrays. The first contains the instance ID of each tree point. Points that do not belong to any instance are assigned the ID \(-1\). The second contains the unique instance IDs. The third contains the 2D segmentation mask with the pixels at the boundary lines between neighboring trees are assigned the background value of 0. The fourth contains the same segmentation mask without boundary lines.

Shape:
  • tree_coords: \((N, 3)\)

  • instance_ids: \((N)\)

  • tree_positions_grid: \((N, 2)\)

  • canopy_height_model: \((W, H)\)

  • grid_origin: \((2)\)

  • Output: Tuple of five arrays. The first has shape \((N)\) and the second \((I)\). The third and fourth have shape \((W, H)\). The fifth has shape \((N, 2)\).

where

\(N = \text{ number of tree points}\)
\(I = \text{ number of tree instances}\)
\(W = \text{ extent of the canopy height model in x-direction}\)
\(H = \text{ extent of the canopy height model in y-direction}\)
compute_crown_distance_fields(
watershed_labels_without_border: ndarray,
tree_positions_grid: ndarray,
) ndarray[source]

Calculates signed distance fields from the 2D segmentation mask of each tree. The distance values specify the 2D distance to the boundary of the segmentation mask. The distance value is negative for pixels that belong to the tree and positive for pixels that do not belong to the tree.

Parameters:
  • watershed_labels_without_border – Watershed labels without boundary lines.

  • tree_positions_grid – Indices of the grid cells of the canopy height model in which the tree positions are located.

Returns:

Signed distance masks for all trees.

Shape:
  • watershed_labels_without_border: \((W, H)\)

  • tree_positions_grid: \((N, 2)\)

  • Output: \((I, W, H)\).

where

\(N = \text{ number of tree points}\)
\(I = \text{ number of tree instances}\)
\(W = \text{ extent of the distance fields in x-direction}\)
\(H = \text{ extent of the distance fields in y-direction}\)
compute_crown_top_positions(
tree_coords: ndarray,
classification: ndarray,
) Tuple[ndarray, ndarray, ndarray, ndarray][source]

Constructs a 2D canopy height model, identifies the local maxima corresponding to the crown tops and calculates their 2D position.

Parameters:
  • tree_coords – Coordinates of all tree points.

  • classification – Class IDs of each tree point.

Returns:

Tuple of four arrays. The first contains the 2D position of each crown top as floating point coordinates. The second contains the tree positions as integer coordinates in the grid coordinate system of the canopy height model. The third contains the canopy height model. The fourth contains the position of the grid origin.

Shape:
  • tree_coords: \((N, 3)\)

  • classification: \((N)\)

  • Output: Tuple of three tensors. The first has shape \((C, 2)\). The second has shape \((W, H)\). The third has shape \((2)\).

where

\(N = \text{ number of tree points}\)
\(C = \text{ number of crown instances}\)
\(W = \text{ extent of the canopy height model in x-direction}\)
\(H = \text{ extent of the canopy height model in y-direction}\)
compute_trunk_positions(
tree_coords: ndarray,
instance_ids: ndarray,
unique_instance_ids: ndarray,
) ndarray[source]

Computes the 2D position of each trunk.

Parameters:
  • tree_coords – Coordinates of all tree points.

  • instance_ids – Instance IDs of each tree point.

  • unique_instance_ids – Unique instance IDs.

  • min_points – Minimum number of points an instance must have to not be discarded.

Returns:

2D position of each trunk.

Shape:
  • tree_coords: \((N, 3)\)

  • instance_ids: \((N)\)

  • unique_instance_ids: \((T)\)

  • Output: \((T, 2)\)

where

\(N = \text{ number of tree points}\)
\(T = \text{ number of trunk instances}\)
static create_height_map(
points: ndarray,
grid_size: float,
bounding_box: ndarray | None = None,
) Tuple[ndarray, ndarray, ndarray][source]

Creates a 2D height map from a given point cloud. For this purpose, the 3D point cloud is projected onto a 2D grid and the maximum z-coordinate within each grid cell is recorded. The value of grid cells that do not contain any point is set to zero. Before creating the height map, the point cloud is normalized by subtracting the minimum z-coordinate.

Parameters:
  • points – Coordinates of each point.

  • grid_size – Grid size of the height map.

  • bounding_box – Bounding box coordinates specifying the area for which to compute the height map. The first element is expected to be the minimum xy-coordinate of the bounding box and the second the maximum xy-coordinate.

Returns:

Tuple of three arrays. The first is the height map. The second is a map of the same size containing the number of points within each grid cell. The third contains the position of the grid origin.

Shape:
  • points: \((N, 3)\)

  • bounding_box: \((2, 2)\)

  • Output: Tuple of three arrays. The first and second have shape \((W, H)\). The third has shape \((2)\).

where

\(N = \text{ number of points}\)
\(W = \text{ extent of the height map in x-direction}\)
\(H = \text{ extent of the height map in y-direction}\)
determine_overlapping_crowns(
tree_coords: ndarray,
classification: ndarray,
instance_ids: ndarray,
unique_instance_ids: ndarray,
canopy_height_model: ndarray,
watershed_labels_with_border: ndarray,
watershed_labels_without_border: ndarray,
) Tuple[ndarray, ndarray][source]

Identifies trees whose crowns overlap with neighboring trees. If such trees have sufficient point density and contain trunk points, their segmentation is refined. For this purpose, the trunk points are selected as seed points for region growing and the instance IDs of the crown points are reset to -1.

Parameters:
  • tree_coords – Coordinates of all tree points.

  • classification – Class ID of all tree points.

  • instance_ids – Instance IDs of each tree point.

  • unique_instance_ids – Unique instance IDs.

  • canopy_height_model – The canopy height model.

  • watershed_labels_with_border – Watershed labels with boundary lines between neighboring trees are assigned the background value of 0.

Returns:

Tuple of two arrays. The first contains a boolean mask indicating which points were selected as seed points

for region growing. The second contains the updated instance IDs for each tree points. For trees whose segmentation is to be refined the instance IDs of the crown points is set to -1.

Shape:
  • tree_coords: \((N, 3)\)

  • classification: \((N)\)

  • instance_ids: \((N)\)

  • unique_instance_ids: \((I)\)

  • canopy_height_model: \((W, H)\)

  • watershed_labels_with_border: \((W, H)\)

  • Output: Tuple of two arrays. Both have shape \((N)\).

where

\(N = \text{ number of tree points}\)
\(I = \text{ number of tree instances}\)
\(W = \text{ extent of the canopy height model in x-direction}\)
\(H = \text{ extent of the canopy height model in y-direction}\)
filter_trunk_clusters(
instance_ids: ndarray,
unique_instance_ids: ndarray,
min_points: int,
) Tuple[ndarray, ndarray][source]

Removes trunk instances with less than min_points from the set of trunk instances.

Parameters:
  • instance_ids – Instance IDs of each tree point.

  • unique_instance_ids – Unique instance IDs.

  • min_points – Minimum number of points an instance must have to not be discarded.

Returns:

Tuple of two arrays. The first contains the updated instance ID of each tree point. Points that do not belong to any instance are assigned the ID \(-1\). The second contains the unique instance IDs.

Shape:
  • instance_ids: \((N)\)

  • unique_instance_ids: \((T)\)

  • Output: Tuple of two arrays. The first has shape \((N)\) and the second \((T')\)

where

\(N = \text{ number of tree points}\)
\(T = \text{ number of trunk instances}\)
\(T' = \text{ number of trunk instances after filtering}\)
grow_trees(
tree_coords: ndarray,
instance_ids: ndarray,
unique_instances_ids: ndarray,
grid_origin: ndarray,
crown_distance_fields,
seed_mask: ndarray,
) ndarray[source]

Region growing segmentation.

Parameters:
  • tree_coords – Coordinates of all tree points.

  • instance_ids – Instance IDs of each tree point.

  • grid_origin – 2D coordinates of the origin of the crown distance fields.

  • crown_distance_fields – 2D signed distance fields idicating the distance to the crown boundary of each tree to segment.

  • seed_mask – Boolean mask indicating which points were selected as seed points for region growing.

Returns:

Updated instance IDs.

Shape:
  • tree_coords: \((N, 3)\)

  • instance_ids: \((N)\)

  • grid_origin: \((2)\)

  • crown_distance_fields: \((I', W, H)\)

  • seed_mask: \((N)\)

  • Output: \((N)\).

where

\(N = \text{ number of tree points}\)
\(I' = \text{ number of tree instances to segment}\)
\(W = \text{ extent of the distance fields in x-direction}\)
\(H = \text{ extent of the distance fields in y-direction}\)
match_trunk_and_crown_tops(
trunk_positions: ndarray,
crown_top_positions: ndarray,
) ndarray[source]

Merges trunk and crown topm positions corresponding to the same tree.

Parameters:
  • trunk_positions – 2D position of each trunk.

  • crown_top_positions – 2D position of each crown top.

Returns:

Merged tree positions.

Shape:
  • trunk_positions: \((T, 2)\)

  • crown_top_positions: \((C, 2)\)

  • Output: \((I, 2)\)

where

\(T = \text{ number of trunk instances}\)
\(C = \text{ number of crown instances}\)
\(I = \text{ number of merged tree instances}\)
watershed_correction(
watershed_labels_with_border: ndarray,
watershed_labels_without_border: ndarray,
tree_positions_grid: ndarray,
point_cloud_id: str | None = None,
) Tuple[ndarray, ndarray][source]

Detects erroneous parts within a Watershed segmentation mask and replaces them by a Voronoi segmentation.

Parameters:
  • watershed_labels_with_border – Uncorrected Watershed labels with boundary lines between neighboring trees are assigned the background value of 0.

  • watershed_labels_without_border – Uncorrected Watershed labels without boundary lines.

  • tree_positions_grid – Indices of the grid cells of the canopy height model in which the tree positions are located.

  • point_cloud_id – ID of the point cloud to be used in the file names of the visualization results. Defaults to None, which means that no visualizations are created.

Returns:

Tuple of two arrays. The first contains the corrected Watershed segmentation mask with the pixels at the boundary lines between neighboring trees are assigned the background value of 0. The second contains the same Watershed segmentation mask without boundary lines.

Shape:
  • watershed_labels_with_border: \((W, H)\)

  • watershed_labels_without_border: \((W, H)\)

  • tree_positions_grid: \((I, 2)\)

  • Output: Tuple of two arrays. Both have shape \((W, H)\).

where

\(I = \text{ number of tree instances}\)
\(W = \text{ extent of the canopy height model in x-direction}\)
\(H = \text{ extent of the canopy height model in y-direction}\)
watershed_segmentation(
canopy_height_model: ndarray,
tree_positions_grid: ndarray,
) Tuple[ndarray, ndarray][source]

Marker-controlled Watershed segmentation of the canopy height model.

Parameters:
  • canopy_height_model – The canopy height model to segment.

  • tree_positions_grid – Indices of the grid cells of the canopy height model in which the tree positions are located.

Returns:

Tuple of two arrays. The first contains the Watershed segmentation mask with the pixels at the boundary lines between neighboring trees are assigned the background value of 0. The second contains the same Watershed segmentation mask without boundary lines.

Shape:
  • canopy_height_model: \((W, H)\)

  • tree_positions_grid: \((I, 2)\)

  • Output: Tuple of two arrays. Both have shape \((W, H)\).

where

\(I = \text{ number of tree instances}\)
\(W = \text{ extent of the canopy height model in x-direction}\)
\(H = \text{ extent of the canopy height model in y-direction}\)
class pointtree.instance_segmentation.PriorityQueue[source]

Bases: object

Priority queue that allows to update the priority of items in the queue.

REMOVED = '_REMOVED_TASK'
add(key: Hashable, entry: Any, priority: float) None[source]

Adds a new item to the queue or updates the priority of an existing item.

Parameters:
  • key – The key of the item to add or update.

  • entry – The data entry of the item.

  • priority – The priority of the item.

get(key: Hashable) Tuple[float, Any] | None[source]

Retrieves an item without removing it from the queue.

Parameters:

key – The key of the item to retrieve.

Returns:

Priority of the retrieved item and its data entry if the key exists, None otherwise.

pop() Tuple[float, Hashable, Any][source]

Removes and returns the lowest priority item.

Returns:

The priority, key, and data entry of the lowest priority item.

Raises:

KeyError – If the queue is empty.

remove(key: Hashable) None[source]

Marks an existing item as REMOVED.

Parameters:

key – The key of the item to be removed.

Raises:

KeyError – If the key does not exist.

update(key: Hashable, new_entry: Any) None[source]

Updates the data entry of an item in the queue without changing its priority.

Parameters:
  • key – The key of the item to update.

  • new_entry – New data entry for the item.

Raises:

KeyError – If the key does not exist.

pointtree.instance_segmentation.remap_instance_ids(
instance_ids: ndarray,
) Tuple[ndarray, ndarray][source]

Remaps instance IDs that may not be continuous to a continuous range.

Parameters:

instance_ids – Instance IDs to remap.

Returns:

Tuple of two arrays. The first contains the remapped instance IDs. The second contains the unique instance IDs after remapping.

Shape:
  • instance_ids: \((N)\)

  • Output: Tuple of two arrays. The first has shape \((N)\) and the second \((I)\)

where

\(N = \text{ number of points}\)
\(I = \text{ number of instances}\)