gala.morpho: Morphological operations

gala.morpho.damify(a, in_place=False)[source]

Add dams to a borderless segmentation.

gala.morpho.get_neighbor_idxs(ar, idxs, connectivity=1)[source]

Return indices of neighboring voxels given array, indices, connectivity.

Parameters:

ar : ndarray

The array in which neighbors are to be found.

idxs : int or container of int

The indices for which to find neighbors.

connectivity : int in {1, 2, …, ar.ndim}

The number of orthogonal steps allowed to be considered a neighbor.

Returns:

neighbor_idxs : 2D array, shape (nidxs, nneighbors)

The neighbor indices for each index passed.

Examples

>>> ar = np.arange(16).reshape((4, 4))
>>> ar
array([[ 0,  1,  2,  3],
       [ 4,  5,  6,  7],
       [ 8,  9, 10, 11],
       [12, 13, 14, 15]])
>>> get_neighbor_idxs(ar, [5, 10], connectivity=1)
array([[ 9,  6,  1,  4],
       [14, 11,  6,  9]])
>>> get_neighbor_idxs(ar, 9, connectivity=2)
array([[13, 10,  5,  8, 14, 12,  6,  4]])
gala.morpho.hminima(a, thresh)[source]

Suppress all minima that are shallower than thresh.

Parameters:

a : array

The input array on which to perform hminima.

thresh : float

Any local minima shallower than this will be flattened.

Returns:

out : array

A copy of the input array with shallow minima suppressed.

gala.morpho.hollowed(ar, skinsize=1)[source]

Return a copy of ar with the center zeroed out.

‘skinsize’ determines how thick of a crust to leave untouched.

gala.morpho.imhmin(a, thresh)

Suppress all minima that are shallower than thresh.

Parameters:

a : array

The input array on which to perform hminima.

thresh : float

Any local minima shallower than this will be flattened.

Returns:

out : array

A copy of the input array with shallow minima suppressed.

gala.morpho.impose_minima(a, minima, connectivity=1)[source]

Transform ‘a’ so that its only regional minima are those in ‘minima’.

Parameters:
‘a’: an ndarray ‘minima’: a boolean array of same shape as ‘a’ ‘connectivity’: the connectivity of the structuring element used in morphological reconstruction.
Value:
an ndarray of same shape as a with unmarked local minima paved over.
gala.morpho.manual_split(probs, seg, body, seeds, connectivity=1, boundary_seeds=None)[source]

Manually split a body from a segmentation using seeded watershed.

Input:
  • probs: the probability of boundary in the volume given.
  • seg: the current segmentation.
  • body: the label to be split.
  • seeds: the seeds for the splitting (should be just two labels).

[-connectivity: the connectivity to use for watershed.] [-boundary_seeds: if not None, these locations become inf in probs.]

Value:
  • the segmentation with the selected body split.
gala.morpho.minimum_seeds(current_seeds, min_seed_coordinates, connectivity=1)[source]

Ensure that each point in given coordinates has its own seed.

gala.morpho.morphological_reconstruction(marker, mask, connectivity=1)[source]

Perform morphological reconstruction of the marker into the mask.

See the Matlab image processing toolbox documentation for details: http://www.mathworks.com/help/toolbox/images/f18-16264.html

gala.morpho.multiscale_regular_seeds(off_limits, num_seeds)[source]

Return evenly-spaced seeds, but thinned in areas with no boundaries.

Parameters:

off_limits : array of bool, shape (M, N)

A binary array where True indicates the position of a boundary, and thus where we don’t want to place seeds.

num_seeds : int

The desired number of seeds.

Returns:

seeds : array of int, shape (M, N)

An array of seed points. Each seed gets its own integer ID, starting from 1.

gala.morpho.non_traversing_segments(a)[source]

Find segments that enter the volume but do not leave it elsewhere.

Parameters:

a : array of int

A segmented volume.

Returns:

nt : 1D array of int

The IDs of any segments not traversing the volume.

Examples

>>> segs = np.array([[1, 2, 3, 3, 4],
...                  [1, 2, 2, 3, 4],
...                  [1, 5, 5, 3, 4],
...                  [1, 1, 5, 3, 4]], int)
>>> non_traversing_segments(segs)
array([1, 2, 4, 5])
gala.morpho.orphans(a)[source]

Find all the segments that do not touch the volume boundary.

This function differs from agglo.Rag.orphans() in that it does not use the graph, but rather computes orphans directly from a volume.

Parameters:

a : array of int

A segmented volume.

Returns:

orph : 1D array of int

The IDs of any segments not touching the volume boundary.

Examples

>>> segs = np.array([[1, 1, 1, 2],
...                  [1, 3, 4, 2],
...                  [1, 2, 2, 2]], int)
>>> orphans(segs)
array([3, 4])
>>> orphans(segs[:2])
array([], dtype=int64)
gala.morpho.raveled_steps_to_neighbors(shape, connectivity=1)[source]

Compute the stepsize along all axes for given connectivity and shape.

Parameters:

shape : tuple of int

The shape of the array along which we are stepping.

connectivity : int in {1, 2, …, len(shape)}

The number of orthogonal steps we can take to reach a “neighbor”.

Returns:

steps : array of int64

The steps needed to get to neighbors from a particular raveled index.

Examples

>>> shape = (5, 4, 9)
>>> steps = raveled_steps_to_neighbors(shape)
>>> sorted(steps)
[-36, -9, -1, 1, 9, 36]
>>> steps2 = raveled_steps_to_neighbors(shape, 2)
>>> sorted(steps2)
[-45, -37, -36, -35, -27, -10, -9, -8, -1, 1, 8, 9, 10, 27, 35, 36, 37, 45]
gala.morpho.regional_minima(a, connectivity=1)[source]

Find the regional minima in an ndarray.

gala.morpho.relabel_connected(im, connectivity=1)[source]

Ensure all labels in im are connected.

Parameters:

im : array of int

The input label image.

connectivity : int in {1, …, im.ndim}, optional

The connectivity used to determine if two voxels are neighbors.

Returns:

im_out : array of int

The relabeled image.

Examples

>>> image = np.array([[1, 1, 2],
...                   [2, 1, 1]])
>>> im_out = relabel_connected(image)
>>> im_out
array([[1, 1, 2],
       [3, 1, 1]])
gala.morpho.remove_merged_boundaries(labels, connectivity=1)[source]

Remove boundaries in a label field when they separate the same region.

By convention, the boundary label is 0, and labels are positive.

Parameters:

labels : array of int

The label field to be processed.

connectivity : int in {1, …, labels.ndim}, optional

The morphological connectivity for considering neighboring voxels.

Returns:

labels_out : array of int

The same label field, with unnecessary boundaries removed.

Examples

>>> labels = np.array([[1, 0, 1], [0, 1, 0], [2, 0, 3]], np.int)
>>> remove_merged_boundaries(labels)
array([[1, 1, 1],
       [0, 1, 0],
       [2, 0, 3]])
gala.morpho.seg_to_bdry(seg, connectivity=1)[source]

Given a borderless segmentation, return the boundary map.

gala.morpho.split_exclusions(image, labels, exclusions, dilation=0, connectivity=1, standard_seeds=False)[source]

Ensure that no segment in ‘labels’ overlaps more than one exclusion.

gala.morpho.undam(seg)[source]

Assign zero-dams to nearest non-zero region.

gala.morpho.watershed(a, seeds=None, connectivity=1, mask=None, smooth_thresh=0.0, smooth_seeds=False, minimum_seed_size=0, dams=False, override_skimage=False, show_progress=False)[source]

Perform the watershed algorithm of Vincent & Soille (1991).

Parameters:

a : np.ndarray, arbitrary shape and type

The input image on which to perform the watershed transform.

seeds : np.ndarray, int or bool type, same shape as a (optional)

The seeds for the watershed. If provided, these are the only basins allowed, and the algorithm proceeds by flooding from the seeds. Otherwise, every local minimum is used as a seed.

connectivity : int, {1, …, a.ndim} (optional, default 1)

The neighborhood of each pixel, defined as in scipy.ndimage.

mask : np.ndarray, type bool, same shape as a. (optional)

If provided, perform watershed only in the parts of a that are set to True in mask.

smooth_thresh : float (optional, default 0.0)

Local minima that are less deep than this threshold are suppressed, using hminima.

smooth_seeds : bool (optional, default False)

Perform binary opening on the seeds, using the same connectivity as the watershed.

minimum_seed_size : int (optional, default 0)

Remove seed regions smaller than this size.

dams : bool (optional, default False)

Place a dam where two basins meet. Set this to True if you require 0-labeled boundaries between different regions.

override_skimage : bool (optional, default False)

skimage.morphology.watershed is used to implement the main part of the algorithm when dams=False. Use this flag to use the separate pure Python implementation instead.

show_progress : bool (optional, default False)

Show a cute little ASCII progress bar (using the progressbar package)

Returns:

ws : np.ndarray, same shape as a, int type.

The watershed transform of the input image.

gala.morpho.watershed_sequence(a, seeds=None, mask=None, axis=0, n_jobs=1, **kwargs)[source]

Perform a watershed on a plane-by-plane basis.

See documentation for watershed for available kwargs.

The watershed algorithm views image intensity as “height” and finds flood basins within it. These basins are then viewed as the different labeled regions of an image.

This function performs watershed on an ndarray on each plane separately, then concatenate the results.

Parameters:

a : numpy ndarray, arbitrary type or shape.

The input image on which to perform the watershed transform.

seeds : bool/int numpy.ndarray, same shape as a (optional, default None)

The seeds for the watershed.

mask : bool numpy.ndarray, same shape as a (optional, default None)

If provided, perform watershed only over voxels that are True in the mask.

axis : int, {1, …, a.ndim} (optional, default: 0)

Which axis defines the plane sequence. For example, if the input image is 3D and axis=1, then the output will be the watershed on a[:, 0, :], a[:, 1, :], a[:, 2, :], … and so on.

n_jobs : int, optional

Use joblib to distribute each plane over given number of processing cores. If -1, multiprocessing.cpu_count is used.

Returns:

ws : numpy ndarray, int type

The labeled watershed basins.

Other Parameters:
 

**kwargs : keyword arguments passed through to the watershed function.