segment#

Module: segment.bundles#

RecoBundles(streamlines, *[, greater_than, ...])

Methods

check_range(streamline, gt, lt)

bundle_adjacency(dtracks0, dtracks1, threshold)

Find bundle adjacency between two given tracks/bundles

ba_analysis(recognized_bundle, expert_bundle, *)

Calculates bundle adjacency score between two given bundles

cluster_bundle(bundle, clust_thr, rng, *[, ...])

Clusters bundles

bundle_shape_similarity(bundle1, bundle2, rng, *)

Calculates bundle shape similarity between two given bundles using bundle adjacency (BA) metric

Module: segment.clustering#

Identity()

Provides identity indexing functionality.

Cluster(*[, id, indices, refdata])

Provides functionalities for interacting with a cluster.

ClusterCentroid(centroid, *[, id, indices, ...])

Provides functionalities for interacting with a cluster.

ClusterMap(*[, refdata])

Provides functionalities for interacting with clustering outputs.

ClusterMapCentroid(*[, refdata])

Provides functionalities for interacting with clustering outputs that have centroids.

Clustering()

Methods

QuickBundles(threshold, *[, metric, ...])

Clusters streamlines using QuickBundles.

QuickBundlesX(thresholds, *[, metric])

Clusters streamlines using QuickBundlesX.

TreeCluster(threshold, centroid, *[, indices])

Attributes:

TreeClusterMap(root)

Attributes:

qbx_and_merge(streamlines, thresholds, *[, ...])

Run QuickBundlesX and then run again on the centroids of the last layer.

Module: segment.fss#

FastStreamlineSearch(ref_streamlines, ...[, ...])

Methods

nearest_from_matrix_row(coo_array)

Return the nearest (smallest) for each row given an coup sparse matrix

nearest_from_matrix_col(coo_array)

Return the nearest (smallest) for each column given an coup sparse matrix

Module: segment.mask#

multi_median(data, median_radius, numpass)

Applies median filter multiple times on input data.

applymask(vol, mask)

Mask vol with mask.

bounding_box(vol)

Compute the bounding box of nonzero intensity voxels in the volume.

crop(vol, mins, maxs)

Crops the input volume.

median_otsu(input_volume, *[, vol_idx, ...])

Simple brain extraction tool method for images from DWI data.

segment_from_cfa(tensor_fit, roi, threshold, *)

Segment the cfa inside roi using the values from threshold as bounds.

clean_cc_mask(mask)

Cleans a segmentation of the corpus callosum so no random pixels are included.

Module: segment.metric#

mdf(s1, s2)

Computes the MDF (Minimum average Direct-Flip) distance between two streamlines.

mean_manhattan_distance(a, b)

Compute the average Manhattan-L1 distance (MDF without flip)

mean_euclidean_distance(a, b)

Compute the average Euclidean-L2 distance (MDF without flip)

Module: segment.threshold#

otsu(image, *[, nbins])

Return threshold value based on Otsu's method.

upper_bound_by_rate(data, *[, rate])

Adjusts upper intensity boundary using rates

upper_bound_by_percent(data, *[, percent])

Find the upper bound for visualization of medical images

Module: segment.tissue#

TissueClassifierHMRF(*[, save_history, verbose])

This class contains the methods for tissue classification using the Markov Random Fields modeling approach.

compute_directional_average(data, bvals, *)

Compute the mean signal for each unique b-value shell and fit a linear model.

dam_classifier(data, bvals, wm_threshold, *)

Computes the P-map (fitting slope) on data to extract white and grey matter.

Module: segment.utils#

remove_holes_and_islands(binary_img, *[, ...])

Remove any small mask chunks or holes that could be in the segmentation output.

RecoBundles#

class dipy.segment.bundles.RecoBundles(streamlines, *, greater_than=50, less_than=1000000, cluster_map=None, clust_thr=15, nb_pts=20, rng=None, verbose=False)[source]#

Bases: object

Methods

evaluate_results(model_bundle, ...)

Compare the similarity between two given bundles, model bundle, and extracted bundle.

recognize(model_bundle, model_clust_thr, *)

Recognize the model_bundle in self.streamlines

refine(model_bundle, pruned_streamlines, ...)

Refine and recognize the model_bundle in self.streamlines This method expects once pruned streamlines as input.

evaluate_results(model_bundle, pruned_streamlines, slr_select)[source]#

Compare the similarity between two given bundles, model bundle, and extracted bundle.

Parameters:
model_bundleStreamlines

Model bundle streamlines.

pruned_streamlinesStreamlines

Pruned bundle streamlines.

slr_selecttuple

Select the number of streamlines from model to neighborhood of model to perform the local SLR.

Returns:
ba_valuefloat

bundle adjacency value between model bundle and pruned bundle

bmd_valuefloat

bundle minimum distance value between model bundle and pruned bundle

recognize(model_bundle, model_clust_thr, *, reduction_thr=10, reduction_distance='mdf', slr=True, num_threads=None, slr_metric=None, slr_x0=None, slr_bounds=None, slr_select=(400, 600), slr_method='L-BFGS-B', pruning_thr=5, pruning_distance='mdf')[source]#

Recognize the model_bundle in self.streamlines

See [1] for further details about the method.

Parameters:
model_bundleStreamlines

model bundle streamlines used as a reference to extract similar streamlines from input tractogram

model_clust_thrfloat

MDF distance threshold for the model bundles

reduction_thrfloat, optional

Reduce search space in the target tractogram by (mm) (default 10)

reduction_distancestring, optional

Reduction distance type can be mdf or mam (default mdf)

slrbool, optional

Use Streamline-based Linear Registration (SLR) locally (default True)

num_threadsint, optional

Number of threads to be used for OpenMP parallelization. If None (default) the value of OMP_NUM_THREADS environment variable is used if it is set, otherwise all available threads are used. If < 0 the maximal number of threads minus \(|num_threads + 1|\) is used (enter -1 to use as many threads as possible). 0 raises an error.

slr_metricBundleMinDistanceMetric
slr_x0array or int or str, optional

Transformation allowed. translation, rigid, similarity or scaling Initial parametrization for the optimization.

If 1D array with:

a) 6 elements then only rigid registration is performed with the 3 first elements for translation and 3 for rotation. b) 7 elements also isotropic scaling is performed (similarity). c) 12 elements then translation, rotation (in degrees), scaling and shearing are performed (affine).

Here is an example of x0 with 12 elements: x0=np.array([0, 10, 0, 40, 0, 0, 2., 1.5, 1, 0.1, -0.5, 0])

This has translation (0, 10, 0), rotation (40, 0, 0) in degrees, scaling (2., 1.5, 1) and shearing (0.1, -0.5, 0).

If int:
  1. 6

    x0 = np.array([0, 0, 0, 0, 0, 0])

  2. 7

    x0 = np.array([0, 0, 0, 0, 0, 0, 1.])

  3. 12

    x0 = np.array([0, 0, 0, 0, 0, 0, 1., 1., 1, 0, 0, 0])

If str:
  1. “rigid”

    x0 = np.array([0, 0, 0, 0, 0, 0])

  2. “similarity”

    x0 = np.array([0, 0, 0, 0, 0, 0, 1.])

  3. “affine”

    x0 = np.array([0, 0, 0, 0, 0, 0, 1., 1., 1, 0, 0, 0])

slr_boundsarray, optional

SLR bounds.

slr_selecttuple, optional

Select the number of streamlines from model to neighborhood of model to perform the local SLR.

slr_methodstring, optional

Optimization method ‘L_BFGS_B’ or ‘Powell’ optimizers can be used. (default ‘L-BFGS-B’)

pruning_thrfloat, optional

Pruning after reducing the search space.

pruning_distancestring, optional

Pruning distance type can be mdf or mam.

Returns:
recognized_transfStreamlines

Recognized bundle in the space of the model tractogram

recognized_labelsarray

Indices of recognized bundle in the original tractogram

References

refine(model_bundle, pruned_streamlines, model_clust_thr, *, reduction_thr=14, reduction_distance='mdf', slr=True, slr_metric=None, slr_x0=None, slr_bounds=None, slr_select=(400, 600), slr_method='L-BFGS-B', pruning_thr=6, pruning_distance='mdf')[source]#

Refine and recognize the model_bundle in self.streamlines This method expects once pruned streamlines as input. It refines the first output of RecoBundles by applying second local slr (optional), and second pruning. This method is useful when we are dealing with noisy data or when we want to extract small tracks from tractograms. This time, search space is created using pruned bundle and not model bundle.

See [1], [2] for further details about the method.

Parameters:
model_bundleStreamlines

model bundle streamlines used as a reference to extract similar streamlines from input tractogram

pruned_streamlinesStreamlines

Recognized bundle from target tractogram by RecoBundles.

model_clust_thrfloat

MDF distance threshold for the model bundles

reduction_thrfloat

Reduce search space by (mm) (default 14)

reduction_distancestring

Reduction distance type can be mdf or mam (default mdf)

slrbool

Use Streamline-based Linear Registration (SLR) locally.

slr_metricBundleMinDistanceMetric

Bundle distance metric.

slr_x0array or int or str

Transformation allowed. translation, rigid, similarity or scaling Initial parametrization for the optimization.

If 1D array with:
  1. 6 elements then only rigid registration is performed with the 3 first elements for translation and 3 for rotation.

  2. 7 elements also isotropic scaling is performed (similarity).

  3. 12 elements then translation, rotation (in degrees), scaling and shearing are performed (affine).

Here is an example of x0 with 12 elements: x0=np.array([0, 10, 0, 40, 0, 0, 2., 1.5, 1, 0.1, -0.5, 0])

This has translation (0, 10, 0), rotation (40, 0, 0) in degrees, scaling (2., 1.5, 1) and shearing (0.1, -0.5, 0).

If int:
  1. 6

    x0 = np.array([0, 0, 0, 0, 0, 0])

  2. 7

    x0 = np.array([0, 0, 0, 0, 0, 0, 1.])

  3. 12

    x0 = np.array([0, 0, 0, 0, 0, 0, 1., 1., 1, 0, 0, 0])

If str:
  1. “rigid”

    x0 = np.array([0, 0, 0, 0, 0, 0])

  2. “similarity”

    x0 = np.array([0, 0, 0, 0, 0, 0, 1.])

  3. “affine”

    x0 = np.array([0, 0, 0, 0, 0, 0, 1., 1., 1, 0, 0, 0])

slr_boundsarray

SLR bounds.

slr_selecttuple

Select the number of streamlines from model to neighborhood of model to perform the local SLR.

slr_methodstring

Optimization method ‘L_BFGS_B’ or ‘Powell’ optimizers can be used.

pruning_thrfloat

Pruning after reducing the search space.

pruning_distancestring

Pruning distance type can be mdf or mam.

Returns:
recognized_transfStreamlines

Recognized bundle in the space of the model tractogram

recognized_labelsarray

Indices of recognized bundle in the original tractogram

References

check_range#

dipy.segment.bundles.check_range(streamline, gt, lt)[source]#

bundle_adjacency#

dipy.segment.bundles.bundle_adjacency(dtracks0, dtracks1, threshold)[source]#

Find bundle adjacency between two given tracks/bundles

See [3] for further details about the method.

Parameters:
dtracks0Streamlines

White matter tract from one subject

dtracks1Streamlines

White matter tract from another subject

thresholdfloat

Threshold controls how much strictness user wants while calculating bundle adjacency between two bundles. Smaller threshold means bundles should be strictly adjacent to get higher BA score.

Returns:
resFloat

Bundle adjacency score between two tracts

References

ba_analysis#

dipy.segment.bundles.ba_analysis(recognized_bundle, expert_bundle, *, nb_pts=20, threshold=6.0)[source]#

Calculates bundle adjacency score between two given bundles

See [3] for further details about the method.

Parameters:
recognized_bundleStreamlines

Extracted bundle from the whole brain tractogram (eg: AF_L)

expert_bundleStreamlines

Model bundle used as reference while extracting similar type bundle from input tractogram

nb_ptsinteger, optional

Discretizing streamlines to have nb_pts number of points

thresholdfloat, optional

Threshold used for in computing bundle adjacency. Threshold controls how much strictness user wants while calculating bundle adjacency between two bundles. Smaller threshold means bundles should be strictly adjacent to get higher BA score.

Returns:
Bundle adjacency score between two tracts

References

cluster_bundle#

dipy.segment.bundles.cluster_bundle(bundle, clust_thr, rng, *, nb_pts=20, select_randomly=500000)[source]#

Clusters bundles

See [3] for further details about the method.

Parameters:
bundleStreamlines

White matter tract

clust_thrfloat

clustering threshold used in quickbundlesX

rngnp.random.Generator

numpy’s random generator for generating random values.

nb_pts: integer, optional

Discretizing streamlines to have nb_points number of points

select_randomly: integer, optional

Randomly select streamlines from the input bundle

Returns:
centroidsStreamlines

clustered centroids of the input bundle

References

bundle_shape_similarity#

dipy.segment.bundles.bundle_shape_similarity(bundle1, bundle2, rng, *, clust_thr=(5, 3, 1.5), threshold=6)[source]#

Calculates bundle shape similarity between two given bundles using bundle adjacency (BA) metric

See [3], [2] for further details about the method.

Parameters:
bundle1Streamlines

White matter tract from one subject (eg: AF_L)

bundle2Streamlines

White matter tract from another subject (eg: AF_L)

rngnp.random.Generator

Random number generator.

clust_thrarray-like, optional

list of clustering thresholds used in quickbundlesX

thresholdfloat, optional

Threshold used for in computing bundle adjacency. Threshold controls how much strictness user wants while calculating shape similarity between two bundles. Smaller threshold means bundles should be strictly similar to get higher shape similarity score.

Returns:
ba_valuefloat

Bundle similarity score between two tracts

References

Identity#

class dipy.segment.clustering.Identity[source]#

Bases: object

Provides identity indexing functionality.

This can replace any class supporting indexing used for referencing (e.g. list, tuple). Indexing an instance of this class will return the index provided instead of the element. It does not support slicing.

Cluster#

class dipy.segment.clustering.Cluster(*, id=0, indices=None, refdata=None)[source]#

Bases: object

Provides functionalities for interacting with a cluster.

Useful container to retrieve index of elements grouped together. If a reference to the data is provided to cluster_map, elements will be returned instead of their index when possible.

Parameters:
cluster_mapClusterMap object

Reference to the set of clusters this cluster is being part of.

idint, optional

Id of this cluster in its associated cluster_map object.

refdatalist, optional

Actual elements that clustered indices refer to.

Methods

assign(*indices)

Assigns indices to this cluster.

Notes

A cluster does not contain actual data but instead knows how to retrieve them using its ClusterMap object.

assign(*indices)[source]#

Assigns indices to this cluster.

Parameters:
*indiceslist of indices

Indices to add to this cluster.

ClusterCentroid#

class dipy.segment.clustering.ClusterCentroid(centroid, *, id=0, indices=None, refdata=None)[source]#

Bases: Cluster

Provides functionalities for interacting with a cluster.

Useful container to retrieve the indices of elements grouped together and the cluster’s centroid. If a reference to the data is provided to cluster_map, elements will be returned instead of their index when possible.

Parameters:
cluster_mapClusterMapCentroid object

Reference to the set of clusters this cluster is being part of.

idint, optional

Id of this cluster in its associated cluster_map object.

refdatalist, optional

Actual elements that clustered indices refer to.

Methods

assign(id_datum, features)

Assigns a data point to this cluster.

update()

Update centroid of this cluster.

Notes

A cluster does not contain actual data but instead knows how to retrieve them using its ClusterMapCentroid object.

assign(id_datum, features)[source]#

Assigns a data point to this cluster.

Parameters:
id_datumint

Index of the data point to add to this cluster.

features2D array

Data point’s features to modify this cluster’s centroid.

update()[source]#

Update centroid of this cluster.

Returns:
convergedbool

Tells if the centroid has moved.

ClusterMap#

class dipy.segment.clustering.ClusterMap(*, refdata=None)[source]#

Bases: object

Provides functionalities for interacting with clustering outputs.

Useful container to create, remove, retrieve and filter clusters. If refdata is given, elements will be returned instead of their index when using Cluster objects.

Parameters:
refdatalist

Actual elements that clustered indices refer to.

Attributes:
clusters
refdata

Methods

add_cluster(*clusters)

Adds one or multiple clusters to this cluster map.

clear()

Remove all clusters from this cluster map.

clusters_sizes()

Gets the size of every cluster contained in this cluster map.

get_large_clusters(min_size)

Gets clusters which contains at least min_size elements.

get_small_clusters(max_size)

Gets clusters which contains at most max_size elements.

remove_cluster(*clusters)

Remove one or multiple clusters from this cluster map.

size()

Gets number of clusters contained in this cluster map.

add_cluster(*clusters)[source]#

Adds one or multiple clusters to this cluster map.

Parameters:
*clustersCluster object, …

Cluster(s) to be added in this cluster map.

clear()[source]#

Remove all clusters from this cluster map.

property clusters#
clusters_sizes()[source]#

Gets the size of every cluster contained in this cluster map.

Returns:
list of int

Sizes of every cluster in this cluster map.

get_large_clusters(min_size)[source]#

Gets clusters which contains at least min_size elements.

Parameters:
min_sizeint

Minimum number of elements a cluster needs to have to be selected.

Returns:
list of Cluster objects

Clusters having at least min_size elements.

get_small_clusters(max_size)[source]#

Gets clusters which contains at most max_size elements.

Parameters:
max_sizeint

Maximum number of elements a cluster can have to be selected.

Returns:
list of Cluster objects

Clusters having at most max_size elements.

property refdata#
remove_cluster(*clusters)[source]#

Remove one or multiple clusters from this cluster map.

Parameters:
*clustersCluster object, …

Cluster(s) to be removed from this cluster map.

size()[source]#

Gets number of clusters contained in this cluster map.

ClusterMapCentroid#

class dipy.segment.clustering.ClusterMapCentroid(*, refdata=None)[source]#

Bases: ClusterMap

Provides functionalities for interacting with clustering outputs that have centroids.

Allows to retrieve easily the centroid of every cluster. Also, it is a useful container to create, remove, retrieve and filter clusters. If refdata is given, elements will be returned instead of their index when using ClusterCentroid objects.

Parameters:
refdatalist

Actual elements that clustered indices refer to.

Attributes:
centroids
clusters
refdata

Methods

add_cluster(*clusters)

Adds one or multiple clusters to this cluster map.

clear()

Remove all clusters from this cluster map.

clusters_sizes()

Gets the size of every cluster contained in this cluster map.

get_large_clusters(min_size)

Gets clusters which contains at least min_size elements.

get_small_clusters(max_size)

Gets clusters which contains at most max_size elements.

remove_cluster(*clusters)

Remove one or multiple clusters from this cluster map.

size()

Gets number of clusters contained in this cluster map.

property centroids#

Clustering#

class dipy.segment.clustering.Clustering[source]#

Bases: object

Methods

cluster(data, *[, ordering])

Clusters data.

abstract cluster(data, *, ordering=None)[source]#

Clusters data.

Subclasses will perform their clustering algorithm here.

Parameters:
datalist of N-dimensional arrays

Each array represents a data point.

orderingiterable of indices, optional

Specifies the order in which data points will be clustered.

Returns:
ClusterMap object

Result of the clustering.

QuickBundles#

class dipy.segment.clustering.QuickBundles(threshold, *, metric='MDF_12points', max_nb_clusters=None)[source]#

Bases: Clustering

Clusters streamlines using QuickBundles.

Given a list of streamlines, the QuickBundles algorithm [3] sequentially assigns each streamline to its closest bundle in \(\mathcal{O}(Nk)\) where \(N\) is the number of streamlines and \(k\) is the final number of bundles. If for a given streamline its closest bundle is farther than threshold, a new bundle is created and the streamline is assigned to it except if the number of bundles has already exceeded max_nb_clusters.

Parameters:
thresholdfloat

The maximum distance from a bundle for a streamline to be still considered as part of it.

metricstr or Metric object, optional

The distance metric to use when comparing two streamlines. By default, the Minimum average Direct-Flip (MDF) distance [3] is used and streamlines are automatically resampled so they have 12 points.

max_nb_clustersint, optional

Limits the creation of bundles.

Methods

cluster(streamlines, *[, ordering])

Clusters streamlines into bundles.

References

Examples

>>> from dipy.segment.clustering import QuickBundles
>>> from dipy.data import get_fnames
>>> from dipy.io.streamline import load_tractogram
>>> from dipy.tracking.streamline import Streamlines
>>> fname = get_fnames(name='fornix')
>>> fornix = load_tractogram(fname, 'same',
...                          bbox_valid_check=False).streamlines
>>> streamlines = Streamlines(fornix)
>>> # Segment fornix with a threshold of 10mm and streamlines resampled
>>> # to 12 points.
>>> qb = QuickBundles(threshold=10.)
>>> clusters = qb.cluster(streamlines)
>>> len(clusters)
4
>>> list(map(len, clusters))
[61, 191, 47, 1]
>>> # Resampling streamlines differently is done explicitly as follows.
>>> # Note this has an impact on the speed and the accuracy (tradeoff).
>>> from dipy.segment.featurespeed import ResampleFeature
>>> from dipy.segment.metricspeed import AveragePointwiseEuclideanMetric
>>> feature = ResampleFeature(nb_points=2)
>>> metric = AveragePointwiseEuclideanMetric(feature)
>>> qb = QuickBundles(threshold=10., metric=metric)
>>> clusters = qb.cluster(streamlines)
>>> len(clusters)
4
>>> list(map(len, clusters))
[58, 142, 72, 28]
cluster(streamlines, *, ordering=None)[source]#

Clusters streamlines into bundles.

Performs quickbundles algorithm using predefined metric and threshold.

Parameters:
streamlineslist of 2D arrays

Each 2D array represents a sequence of 3D points (points, 3).

orderingiterable of indices, optional

Specifies the order in which data points will be clustered.

Returns:
ClusterMapCentroid object

Result of the clustering.

QuickBundlesX#

class dipy.segment.clustering.QuickBundlesX(thresholds, *, metric='MDF_12points')[source]#

Bases: Clustering

Clusters streamlines using QuickBundlesX.

See [4] for further details about the method.

Parameters:
thresholdslist of float

Thresholds to use for each clustering layer. A threshold represents the maximum distance from a cluster for a streamline to be still considered as part of it.

metricstr or Metric object, optional

The distance metric to use when comparing two streamlines. By default, the Minimum average Direct-Flip (MDF) distance [3] is used and streamlines are automatically resampled so they have 12 points.

Methods

cluster(streamlines, *[, ordering])

Clusters streamlines into bundles.

References

cluster(streamlines, *, ordering=None)[source]#

Clusters streamlines into bundles.

Performs QuickbundleX using a predefined metric and thresholds.

Parameters:
streamlineslist of 2D arrays

Each 2D array represents a sequence of 3D points (points, 3).

orderingiterable of indices

Specifies the order in which data points will be clustered.

Returns:
TreeClusterMap object

Result of the clustering.

TreeCluster#

class dipy.segment.clustering.TreeCluster(threshold, centroid, *, indices=None)[source]#

Bases: ClusterCentroid

Attributes:
is_leaf

Methods

assign(id_datum, features)

Assigns a data point to this cluster.

update()

Update centroid of this cluster.

add

return_indices

add(child)[source]#
property is_leaf#
return_indices()[source]#

TreeClusterMap#

class dipy.segment.clustering.TreeClusterMap(root)[source]#

Bases: ClusterMap

Attributes:
clusters
refdata

Methods

add_cluster(*clusters)

Adds one or multiple clusters to this cluster map.

clear()

Remove all clusters from this cluster map.

clusters_sizes()

Gets the size of every cluster contained in this cluster map.

get_large_clusters(min_size)

Gets clusters which contains at least min_size elements.

get_small_clusters(max_size)

Gets clusters which contains at most max_size elements.

remove_cluster(*clusters)

Remove one or multiple clusters from this cluster map.

size()

Gets number of clusters contained in this cluster map.

get_clusters

iter_preorder

traverse_postorder

get_clusters(wanted_level)[source]#
iter_preorder(node)[source]#
property refdata#
traverse_postorder(node, visit)[source]#

qbx_and_merge#

dipy.segment.clustering.qbx_and_merge(streamlines, thresholds, *, nb_pts=20, select_randomly=None, rng=None, verbose=False)[source]#

Run QuickBundlesX and then run again on the centroids of the last layer.

Running again QuickBundles at a layer has the effect of merging some of the clusters that may be originally divided because of branching. This function help obtain a result at a QuickBundles quality but with QuickBundlesX speed. The merging phase has low cost because it is applied only on the centroids rather than the entire dataset.

See [3] and [4] for further details about the method.

Parameters:
streamlinesStreamlines

Streamlines.

thresholdssequence

List of distance thresholds for QuickBundlesX.

nb_ptsint

Number of points for discretizing each streamline

select_randomlyint

Randomly select a specific number of streamlines. If None all the streamlines are used.

rngnumpy.random.Generator

If None then generator is initialized internally.

verbosebool, optional.

If True, log information. Default False.

Returns:
clustersobj

Contains the clusters of the last layer of QuickBundlesX after merging.

References

FastStreamlineSearch#

class dipy.segment.fss.FastStreamlineSearch(ref_streamlines, max_radius, *, nb_mpts=4, bin_size=20.0, resampling=24, bidirectional=True)[source]#

Bases: object

Methods

radius_search(streamlines, radius, *[, ...])

Radius Search using Fast Streamline Search

Radius Search using Fast Streamline Search

For each given streamlines, return all reference streamlines within the given radius. See [5] for further details.

Parameters:
streamlinesStreamlines

Streamlines to generate the tree structure.

radiusfloat

Search radius (with MDF / average L2 distance) must be smaller than max_radius when FFS was initialized.

use_negativebool, optional

When used with bidirectional, negative values are returned for reversed order neighbors.

Returns:
resscipy COOrdinates sparse matrix (nb_slines x nb_slines_ref)

Adjacency matrix containing all neighbors within the given radius

Notes

Given streamlines should be already aligned with ref streamlines. Preferably in millimeter space (voxmm or rasmm).

References

nearest_from_matrix_row#

dipy.segment.fss.nearest_from_matrix_row(coo_array)[source]#

Return the nearest (smallest) for each row given an coup sparse matrix

Parameters:
coo_arrayscipy COOrdinates sparse array (nb_slines x nb_slines_ref)

Adjacency matrix containing all neighbors within the given radius

Returns:
non_zero_idsnumpy array (nb_non_empty_row x 1)

Indices of each non-empty slines (row)

nearest_idnumpy array (nb_non_empty_row x 1)

Indices of the nearest reference match (column)

nearest_distnumpy array (nb_non_empty_row x 1)

Distance for each nearest match

nearest_from_matrix_col#

dipy.segment.fss.nearest_from_matrix_col(coo_array)[source]#

Return the nearest (smallest) for each column given an coup sparse matrix

Parameters:
coo_arrayscipy COOrdinates sparse matrix (nb_slines x nb_slines_ref)

Adjacency matrix containing all neighbors within the given radius

Returns:
non_zero_idsnumpy array (nb_non_empty_col x 1)

Indices of each non-empty reference (column)

nearest_idnumpy array (nb_non_empty_col x 1)

Indices of the nearest slines match (row)

nearest_distnumpy array (nb_non_empty_col x 1)

Distance for each nearest match

multi_median#

dipy.segment.mask.multi_median(data, median_radius, numpass)[source]#

Applies median filter multiple times on input data.

Parameters:
datandarray

The input volume to apply filter on.

median_radiusint

Radius (in voxels) of the applied median filter

numpass: int

Number of pass of the median filter

Returns:
datandarray

Filtered input volume.

applymask#

dipy.segment.mask.applymask(vol, mask)[source]#

Mask vol with mask.

Parameters:
volndarray

Array with \(V\) dimensions

maskndarray

Binary mask. Has \(M\) dimensions where \(M <= V\). When \(M < V\), we append \(V - M\) dimensions with axis length 1 to mask so that mask will broadcast against vol. In the typical case vol can be 4D, mask can be 3D, and we append a 1 to the mask shape which (via numpy broadcasting) has the effect of applying the 3D mask to each 3D slice in vol (vol[..., 0] to vol[..., -1).

Returns:
masked_volndarray

vol multiplied by mask where mask may have been extended to match extra dimensions in vol

bounding_box#

dipy.segment.mask.bounding_box(vol)[source]#

Compute the bounding box of nonzero intensity voxels in the volume.

Parameters:
volndarray

Volume to compute bounding box on.

Returns:
npminslist

Array containing minimum index of each dimension

npmaxslist

Array containing maximum index of each dimension

crop#

dipy.segment.mask.crop(vol, mins, maxs)[source]#

Crops the input volume.

Parameters:
volndarray

Volume to crop.

minsarray

Array containing minimum index of each dimension.

maxsarray

Array containing maximum index of each dimension.

Returns:
volndarray

The cropped volume.

median_otsu#

dipy.segment.mask.median_otsu(input_volume, *, vol_idx=None, median_radius=4, numpass=4, autocrop=False, dilate=None, finalize_mask=False)[source]#

Simple brain extraction tool method for images from DWI data.

It uses a median filter smoothing of the input_volumes vol_idx and an automatic histogram Otsu thresholding technique, hence the name median_otsu.

This function is inspired from Mrtrix’s bet which has default values median_radius=3, numpass=2. However, from tests on multiple 1.5T and 3T data from GE, Philips, Siemens, the most robust choice is median_radius=4, numpass=4.

Parameters:
input_volumendarray

3D or 4D array of the brain volume.

vol_idxNone or array, optional

1D array representing indices of axis=3 of a 4D input_volume. None is only an acceptable input if input_volume is 3D.

median_radiusint, optional

Radius (in voxels) of the applied median filter.

numpass: int, optional

Number of pass of the median filter.

autocrop: bool, optional

if True, the masked input_volume will also be cropped using the bounding box defined by the masked data. Should be on if DWI is upsampled to 1x1x1 resolution.

dilateNone or int, optional

number of iterations for binary dilation

finalize_maskbool, optional

Whether to remove potential holes or islands. Useful for solving minor errors.

Returns:
maskedvolumendarray

Masked input_volume

mask3D ndarray

The binary brain mask

Notes

Copyright (C) 2011, the scikit-image team All rights reserved.

Redistribution and use in source and binary forms, with or without modification, are permitted provided that the following conditions are met:

  1. Redistributions of source code must retain the above copyright notice, this list of conditions and the following disclaimer.

  2. Redistributions in binary form must reproduce the above copyright notice, this list of conditions and the following disclaimer in the documentation and/or other materials provided with the distribution.

  3. Neither the name of skimage nor the names of its contributors may be used to endorse or promote products derived from this software without specific prior written permission.

THIS SOFTWARE IS PROVIDED BY THE AUTHOR “AS IS” AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE AUTHOR BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.

segment_from_cfa#

dipy.segment.mask.segment_from_cfa(tensor_fit, roi, threshold, *, return_cfa=False)[source]#

Segment the cfa inside roi using the values from threshold as bounds.

Parameters:
tensor_fitTensorFit object

TensorFit object

roindarray

A binary mask, which contains the bounding box for the segmentation.

thresholdarray-like

An iterable that defines the min and max values to use for the thresholding. The values are specified as (R_min, R_max, G_min, G_max, B_min, B_max)

return_cfabool, optional

If True, the cfa is also returned.

Returns:
maskndarray

Binary mask of the segmentation.

cfandarray, optional

Array with shape = (…, 3), where … is the shape of tensor_fit. The color fractional anisotropy, ordered as a nd array with the last dimension of size 3 for the R, G and B channels.

clean_cc_mask#

dipy.segment.mask.clean_cc_mask(mask)[source]#

Cleans a segmentation of the corpus callosum so no random pixels are included.

Parameters:
maskndarray

Binary mask of the coarse segmentation.

Returns:
new_cc_maskndarray

Binary mask of the cleaned segmentation.

mdf#

dipy.segment.metric.mdf(s1, s2)[source]#

Computes the MDF (Minimum average Direct-Flip) distance between two streamlines.

Streamlines must have the same number of points.

See [3] for a definition of the distance.

Parameters:
s12D array

A streamline (sequence of N-dimensional points).

s22D array

A streamline (sequence of N-dimensional points).

Returns:
double

Distance between two streamlines.

References

mean_manhattan_distance#

dipy.segment.metric.mean_manhattan_distance(a, b)[source]#

Compute the average Manhattan-L1 distance (MDF without flip)

Arrays are representing a single streamline or a list of streamlines that have the same number of N-dimensional points (two last axis).

Parameters:
a2D or 3D array

A streamline or concatenated streamlines (array of S streamlines by P points in N dimension).

b2D or 3D array

A streamline or concatenated streamlines (array of S streamlines by P points in N dimension).

Returns:
1D array

Distance between each S streamlines

mean_euclidean_distance#

dipy.segment.metric.mean_euclidean_distance(a, b)[source]#

Compute the average Euclidean-L2 distance (MDF without flip)

Arrays are representing a single streamline or a list of streamlines that have the same number of N-dimensional points (two last axis).

Parameters:
a2D or 3D array

A streamline or concatenated streamlines (array of S streamlines by P points in N dimension).

b2D or 3D array

A streamline or concatenated streamlines (array of S streamlines by P points in N dimension).

Returns:
1D array

Distance between each S streamlines

otsu#

dipy.segment.threshold.otsu(image, *, nbins=256)[source]#

Return threshold value based on Otsu’s method. Copied from scikit-image to remove dependency.

Parameters:
imagearray

Input image.

nbinsint

Number of bins used to calculate histogram. This value is ignored for integer arrays.

Returns:
thresholdfloat

Threshold value.

upper_bound_by_rate#

dipy.segment.threshold.upper_bound_by_rate(data, *, rate=0.05)[source]#

Adjusts upper intensity boundary using rates

It calculates the image intensity histogram, and based on the rate value it decide what is the upperbound value for intensity normalization, usually lower bound is 0. The rate is the ratio between the amount of pixels in every bins and the bins with highest pixel amount

Parameters:
datafloat

Input intensity value data

ratefloat

representing the threshold whether a specific histogram bin that should be count in the normalization range

Returns:
highfloat

the upper_bound value for normalization

upper_bound_by_percent#

dipy.segment.threshold.upper_bound_by_percent(data, *, percent=1)[source]#

Find the upper bound for visualization of medical images

Calculate the histogram of the image and go right to left until you find the bound that contains more than a percentage of the image.

Parameters:
datandarray
percentfloat
Returns:
upper_boundfloat

TissueClassifierHMRF#

class dipy.segment.tissue.TissueClassifierHMRF(*, save_history=False, verbose=True)[source]#

Bases: object

This class contains the methods for tissue classification using the Markov Random Fields modeling approach.

Methods

classify(image, nclasses, beta, *[, ...])

This method uses the Maximum a posteriori - Markov Random Field approach for segmentation by using the Iterative Conditional Modes and Expectation Maximization to estimate the parameters.

classify(image, nclasses, beta, *, tolerance=1e-05, max_iter=100)[source]#

This method uses the Maximum a posteriori - Markov Random Field approach for segmentation by using the Iterative Conditional Modes and Expectation Maximization to estimate the parameters.

Parameters:
imagendarray,

3D structural image.

nclassesint,

Number of desired classes.

betafloat,

Smoothing parameter, the higher this number the smoother the output will be.

tolerance: float, optional

Value that defines the percentage of change tolerated to prevent the ICM loop to stop. Default is 1e-05. If you want tolerance check to be disabled put ‘tolerance = 0’.

max_iterint, optional

Fixed number of desired iterations. Default is 100. This parameter defines the maximum number of iterations the algorithm will perform. The loop may terminate early if the change in energy sum between iterations falls below the threshold defined by tolerance. However, if tolerance is explicitly set to 0, this early stopping mechanism is disabled, and the algorithm will run for the specified number of iterations unless another stopping criterion is met.

Returns:
initial_segmentationndarray,

3D segmented image with all tissue types specified in nclasses.

final_segmentationndarray,

3D final refined segmentation containing all tissue types.

PVEndarray,

3D probability map of each tissue type.

compute_directional_average#

dipy.segment.tissue.compute_directional_average(data, bvals, *, s0_map=None, masks=None, b0_mask=None, b0_threshold=50, low_signal_threshold=50)[source]#

Compute the mean signal for each unique b-value shell and fit a linear model.

Parameters:
datandarray

The diffusion MRI data.

bvalsndarray

The b-values corresponding to the diffusion data.

s0_mapndarray, optional

Precomputed mean signal map for b=0 images.

masksndarray, optional

Precomputed masks for each unique b-value shell.

b0_maskndarray, optional

Precomputed mask for b=0 images.

b0_thresholdfloat, optional

The intensity threshold for a b=0 image.

low_signal_thresholdfloat, optional

The threshold below which a voxel is considered to have low signal.

Returns:
Pfloat

The slope of the linear model.

Vfloat

The intercept of the linear model.

dam_classifier#

dipy.segment.tissue.dam_classifier(data, bvals, wm_threshold, *, b0_threshold=50, low_signal_threshold=50)[source]#

Computes the P-map (fitting slope) on data to extract white and grey matter.

See [6] for further details about the method.

Parameters:
datandarray

The diffusion MRI data.

bvalsndarray

The b-values corresponding to the diffusion data.

wm_thresholdfloat

The threshold below which a voxel is considered white matter.

b0_thresholdfloat, optional

The intensity threshold for a b=0 image.

low_signal_thresholdfloat, optional

The threshold below which a voxel is considered to have low signal.

Returns:
wm_maskndarray

A binary mask for white matter.

gm_maskndarray

A binary mask for grey matter.

References

remove_holes_and_islands#

dipy.segment.utils.remove_holes_and_islands(binary_img, *, slice_wise=False)[source]#

Remove any small mask chunks or holes that could be in the segmentation output.

Parameters:
binary_imgnp.ndarray

Binary image

slice_wisebool, optional

Whether to run slice wise background correction as well

Returns:
largest_imgnp.ndarray