cedalion.geometry.segmentation

Functionality to work with segmented MRI scans.

Functions

binaryMaskGenerator(d)

Convert multi-channel probability maps to mutually exclusive binary masks.

cell_coordinates(volume[, flat])

Generate voxel-centre coordinates from a 3D volume.

segmentation_postprocessing(segmentation_dir)

Postprocess SPM12 tissue-probability maps into binary segmentation masks.

sizeOfObject(img[, conn])

Return connected-component sizes for a binary 3-D image, sorted descending.

surface_from_segmentation(segmentation_mask, ...)

Create a surface from a segmentation mask.

voxels_from_segmentation(segmentation_mask, ...)

Generate voxels from a segmentation mask.

cedalion.geometry.segmentation.voxels_from_segmentation(
segmentation_mask: DataArray,
segmentation_types: List[str],
isovalue=0.9,
fill_holes_in_mask=False,
) Voxels[source]

Generate voxels from a segmentation mask.

Parameters:
  • segmentation_mask – xr.DataArray Segmentation mask.

  • segmentation_types – List[str] List of segmentation types.

  • isovalue – float, optional Isovalue for marching cubes, by default 0.9.

  • fill_holes_in_mask – bool, optional Fill holes in the mask, by default False.

Returns:

cdc.Voxels

Voxels in voxel space.

cedalion.geometry.segmentation.surface_from_segmentation(
segmentation_mask: DataArray,
segmentation_types: List[str],
isovalue=0.9,
fill_holes_in_mask=False,
) Surface[source]

Create a surface from a segmentation mask.

Parameters:
  • segmentation_mask (xr.DataArray) – Segmentation mask with dimensions segmentation type, i, j, k.

  • segmentation_types (List[str]) – A list of segmentation types to include in the surface.

  • isovalue (Float) – The isovalue to use for the marching cubes algorithm.

  • fill_holes_in_mask (Bool) – Whether to fill holes in the mask before creating the surface.

Returns:

A cedalion.Surface object.

cedalion.geometry.segmentation.cell_coordinates(volume, flat: bool = False)[source]

Generate voxel-centre coordinates from a 3D volume.

Parameters:
  • volume – 3-D array whose shape defines the grid dimensions.

  • flat – If True, return a flat (N, 3) array with label, i, j, k coordinates; otherwise return a (ni, nj, nk, 3) grid array.

Returns:

xr.DataArray of voxel-centre coordinates in voxel (ijk) space, quantified as dimensionless pint units.

cedalion.geometry.segmentation.segmentation_postprocessing(
segmentation_dir: str,
mask_files: dict[str, str] = {'air': 'c6.nii', 'bone': 'c4.nii', 'csf': 'c3.nii', 'gray': 'c1.nii', 'skin': 'c5.nii', 'white': 'c2.nii'},
isSmooth: bool = True,
fixCSF: bool = True,
removeDisconnected: bool = True,
labelUnassigned: bool = True,
removeAir: bool = True,
subtractTissues: bool = True,
) dict[source]

Postprocess SPM12 tissue-probability maps into binary segmentation masks.

Applies a configurable sequence of steps: Gaussian smoothing, CSF continuity repair, removal of disconnected components, labelling of unassigned voxels, removal of air cavities, and tissue subtraction. Results are saved as NIfTI files alongside the originals.

Parameters:
  • segmentation_dir – Directory containing the SPM12 output files.

  • mask_files – Mapping from tissue name to filename within segmentation_dir. Defaults to the standard SPM12 output names (c1.niic6.nii).

  • isSmooth – Apply Gaussian smoothing to each tissue probability map before thresholding.

  • fixCSF – Repair CSF continuity between brain tissue and bone.

  • removeDisconnected – Remove small disconnected components from each mask.

  • labelUnassigned – Assign unlabelled voxels to the nearest tissue type via Gaussian distance blurring.

  • removeAir – Remove air cavities outside the head.

  • subtractTissues – Enforce non-overlapping masks by subtracting inner tissues from outer ones (inside-out order).

Returns:

Dictionary mapping tissue names to the filenames of the saved postprocessed NIfTI masks.

References

Huang et al. [HDS+13], Harmening and Miklody [HM22]

cedalion.geometry.segmentation.binaryMaskGenerator(d)[source]

Convert multi-channel probability maps to mutually exclusive binary masks.

For each voxel the tissue with the highest probability wins; ties go to the first tissue (index 0 = empty).

Parameters:

d – Array of shape (n_tissues, i, j, k) with tissue probability values.

Returns:

Binary mask array of shape (n_tissues + 1, i, j, k) where index 0 is the “empty/unassigned” class and indices 1…n_tissues correspond to the input tissues.

cedalion.geometry.segmentation.sizeOfObject(img, conn=None)[source]

Return connected-component sizes for a binary 3-D image, sorted descending.

Parameters:
  • img – Binary 3-D (or 2-D) NumPy array.

  • conn – Connectivity (number of neighbours to consider). Defaults to 8 for 2-D and 26 for 3-D images.

Returns:

  • size_descend (list[int]): component sizes sorted in descending order.

  • ind (list[int]): original component indices sorted by descending size.

Return type:

Tuple (size_descend, ind)