scilpy.tractanalysis package

scilpy.tractanalysis.afd_along_streamlines module

scilpy.tractanalysis.afd_along_streamlines.afd_map_along_streamlines(sft, fodf, fodf_basis, length_weighting, is_legacy=True)[source]

Compute the mean Apparent Fiber Density (AFD) [1] and mean Radial fODF (radfODF) maps along a bundle.

[1] Raffelt et. al (2012). Apparent fibre density: a novel measure for the analysis of diffusion-weighted magnetic resonance images. Neuroimage, 59(4), 3976-3994.

Parameters:
  • sft (StatefulTractogram) – StatefulTractogram containing the streamlines needed.

  • fodf (nibabel.image) – fODF with shape (X, Y, Z, #coeffs) coeffs depending on the sh_order

  • fodf_basis (string) – Has to be descoteaux07 or tournier07

  • length_weighting (bool) – If set, will weigh the AFD values according to segment lengths

  • is_legacy (bool) – If true, uses legacy basis.

Returns:

  • afd_sum (np.array) – AFD map (weighted if length_weighting)

  • rd_sum (np.array) – rdAFD map (weighted if length_weighting)

scilpy.tractanalysis.bingham_metric_along_streamlines module

scilpy.tractanalysis.bingham_metric_along_streamlines.bingham_metric_map_along_streamlines(sft, bingham_coeffs, metric, max_theta, length_weighting)[source]

Compute mean map for a given Bingham metric along streamlines.

Parameters:
  • sft (StatefulTractogram) – StatefulTractogram containing the streamlines needed.

  • bingham_coeffs (ndarray) – Array of shape (X, Y, Z, N_LOBES, NB_PARAMS) containing the Bingham distributions parameters. Note, NB_PARAMS is usually 7.

  • metric (ndarray) – Array of shape (X, Y, Z) containing the Bingham metric of interest.

  • max_theta (float) – Maximum angle in degrees between the fiber direction and the Bingham peak direction.

  • length_weighting (bool) – If True, will weigh the metric values according to segment lengths.

scilpy.tractanalysis.bingham_metric_along_streamlines.bingham_metric_sum_along_streamlines(sft, bingham_coeffs, metric, max_theta, length_weighting)[source]

Compute a sum map along a bundle for a given Bingham metric.

Parameters:
  • sft (StatefulTractogram) – StatefulTractogram containing the streamlines needed.

  • bingham_coeffs (ndarray (X, Y, Z, N_LOBES, NB_PARAMS)) – Bingham distributions parameters volume.

  • metric (ndarray (X, Y, Z)) – The Bingham metric of interest.

  • max_theta (float) – Maximum angle in degrees between the fiber direction and the Bingham peak direction.

  • length_weighting (bool) – If True, will weight the metric values according to segment lengths.

Returns:

  • metric_sum_map (np.array) – Bingham metric sum map.

  • weight_map (np.array) – Segment lengths.

scilpy.tractanalysis.bundle_operations module

scilpy.tractanalysis.bundle_operations.compute_bundle_diameter(sft, data_labels, fitting_func)[source]

Compute the bundle diameter, in each section.

Parameters:
  • sft (StatefulTractogram) – The bundle.

  • data_labels (np.ndarray) – The labels maps, with one label per section of the bundle.

  • fitting_func (str) – One of [‘lin_up’, ‘lin_down’, ‘exp’, ‘inv’, ‘log’], or None.

Returns:

  • bundle_dict (dict) – The bundle diameter per section, where keys are the sections’ labels, and bundle_dict[label1] = {‘mean’: the_mean, ‘std’: the_std} (STD is in fact an error measure, NOT the standard deviation).

  • centroid_smooth (np.ndarray) – The centroids per section. size [nb_labels, 3]

  • radius (np.ndarray) – The radius per section. size [nb_labels, ]

  • error (np.ndarray) – The error (mean squared error) per section. size [nb_labels, ]

  • pts_labels (list[list]) – The labels associated with each streamline.

scilpy.tractanalysis.bundle_operations.detect_ushape(sft, minU, maxU)[source]

Extract streamlines depending on their “u-shapeness”.

Parameters:
  • sft (Statefull tractogram) – Tractogram used to extract streamlines depending on their ushapeness.

  • minU (Float) – Minimum ufactor of a streamline.

  • maxU (Float) – Maximum ufactor of a streamline.

Returns:

list – Only the ids are returned so proper filtering can be done afterwards.

Return type:

the ids of u-shaped streamlines

scilpy.tractanalysis.bundle_operations.get_streamlines_centroid(streamlines, nb_points)[source]

Compute centroid from streamlines using QuickBundles.

Parameters:
  • streamlines (list of ndarray) – The list of streamlines from which we compute the centroid.

  • nb_points (int) – Number of points defining the centroid streamline.

Return type:

List of length one, containing a np.ndarray of shape (nb_points, 3)

scilpy.tractanalysis.bundle_operations.outliers_removal_using_hierarchical_quickbundles(streamlines, nb_points=12, min_threshold=0.5, nb_samplings_max=30, sampling_seed=1234, fast_approx=False)[source]

Classify inliers and outliers from a list of streamlines.

Parameters:
  • streamlines (list of ndarray) – The list of streamlines from which inliers and outliers are separated.

  • min_threshold (float) – Quickbundles distance threshold for the last threshold.

  • nb_samplings_max (int) – Number of run executed to explore the search space. A different sampling is used each time.

  • sampling_seed (int) – Random number generation initialization seed.

Returns:

ndarray

Return type:

Float value representing the 0-1 score for each streamline

scilpy.tractanalysis.bundle_operations.prune(streamlines, threshold, features)[source]

Discriminate streamlines based on a metrics, usually summary from function outliers_removal_using_hierarchical_quickbundles. :param streamlines: The list of streamlines from which inliers and outliers are separated. :type streamlines: list of ndarray :param threshold: Threshold use to discriminate streamlines using the feature. :type threshold: float :param features: Values that represent a relevant metric to disciminate streamlines. :type features: ndarray

Returns:

Indices for outliers (below threshold), indices for inliers (above threshold).

Return type:

tuple

scilpy.tractanalysis.bundle_operations.remove_outliers_qb(streamlines, threshold, nb_points=12, nb_samplings=30, fast_approx=False)[source]

Wrapper to classify inliers and outliers from a list of streamlines. Uses Quickbundles to separate streamlines that are different.

Parameters:
  • streamlines (list of ndarray) – The list of streamlines from which inliers and outliers are separated.

  • threshold (float) – Quickbundles distance threshold for the last threshold.

  • nb_points (int)

  • nb_samplings (int)

  • fast_approx (bool)

Returns:

list: streamlines considered inliers list: streamlines considered outliers

Return type:

A tuple containing

scilpy.tractanalysis.bundle_operations.uniformize_bundle_sft(sft, axis=None, ref_bundle=None, swap=False)[source]

Uniformize the streamlines in the given tractogram.

Parameters:
  • sft (StatefulTractogram) – The tractogram that contains the list of streamlines to be uniformized

  • axis (int, optional) – Orient endpoints in the given axis

  • ref_bundle (streamlines) – Orient endpoints the same way as this bundle (or centroid)

  • swap (boolean, optional) – Swap the orientation of streamlines

scilpy.tractanalysis.bundle_operations.uniformize_bundle_sft_using_mask(sft, mask, swap=False)[source]

Uniformize the streamlines in the given tractogram so head is closer to to a region of interest.

Parameters:
  • sft (StatefulTractogram) – The tractogram that contains the list of streamlines to be uniformized

  • mask (np.ndarray) – Mask to use as a reference for the ROI.

  • swap (boolean, optional) – Swap the orientation of streamlines

scilpy.tractanalysis.connectivity_segmentation module

scilpy.tractanalysis.connectivity_segmentation.compute_connectivity(indices, atlas_data, real_labels, segmenting_func)[source]

Segments a tractogram into “bundles”, or “connections” between all pairs of labels.

Parameters:
  • indices (ArraySequence) – The list of 3D indices [i, j, k] of all voxels traversed by all streamlines. This is the output of the streamlines_to_voxel_coordinates function.

  • atlas_data (np.ndarray) – The loaded image containing the labels.

  • real_labels (np.ndarray) – The list of labels of interest in the image.

  • segmenting_func (Callable) – The function used for segmentation. Ex: extract_longest_segments_from_profile

Returns:

connectivity – A dict containing one key per real_labels (ex, 1, 2) (starting point).

  • The value of connectivity[1] is again a dict with again the

    real_labels as keys.

  • The value of connectivity[1][2] is a list of length n, where n is

    the number of streamlines starting in 1 and finishing in 2. Each value is a dict of the following shape:

    >>> 'strl_idx': int --> The idex of the streamline in the raw data.
    >>> 'in_idx:    int -->
    >>> 'out_idx': int  -->
    

Return type:

dict

scilpy.tractanalysis.connectivity_segmentation.construct_hdf5_from_connectivity(sft, indices, points_to_idx, real_labels, con_info, hdf5_file, saving_options, out_paths, prune_from_length, min_length, max_length, remove_loops, loop_max_angle, remove_outliers, outlier_threshold, remove_curv_dev, curv_qb_distance, nbr_cpu)[source]
Parameters:
  • sft (StatefulTractogram) – The tractogram.

  • indices (ArraySequence) – Results from streamlines_to_voxel_coordinates.

  • points_to_idx (ArraySequence) – Results from streamlines_to_voxel_coordinates.

  • real_labels (np.ndarray) – The labels.

  • con_info (dict) – The result from compute_connectivity.

  • hdf5_file (hdf5 file) – The opened hdf5_file to which to add the bundles (as groups).

  • saving_options (dict) – Steps for which intermediate files should be saved on disk (not in the hdf5). Keys are: ‘raw’, ‘intermediate’, ‘discarded’, ‘final’. Values are True or False.

  • out_paths (dict) – Name of the intermediate files. Keys are: ‘raw’, ‘invalid_length’, ‘valid_length’, ‘loops’, ‘outliers’, ‘qb_curv’, ‘no_loops’, ‘inliers’, ‘final’. They will be saved if not stated otherwise in saving_options.

  • prune_from_length (bool) – If true, limit length between [min_length, max_length]. Else, skip pruning (step 1).

  • min_length (float)

  • max_length (float)

  • remove_loops (bool) – If true, remove looping streamlines. Else skip step 2.

  • loop_max_angle (float)

  • remove_outliers (bool) – If true, remove outliers using Quickbundles. Else skip step 3.

  • outlier_threshold (float)

  • remove_curv_dev (bool) – If true, remove sharp turns base on Quickbundles. Else skip step 4.

  • curv_qb_distance (float)

  • nbr_cpu (int) – Number of cpu for steps allowing multiprocessing.

scilpy.tractanalysis.connectivity_segmentation.extract_longest_segments_from_profile(strl_indices, atlas_data)[source]

For one given streamline, find the labels at both ends.

Parameters:
  • strl_indices (np.ndarray) – The indices of all voxels traversed by this streamline.

  • atlas_data (np.ndarray) – The loaded image containing the labels.

Returns:

segments_info – A list of length 1 with the information dict if , else, an empty list.

Return type:

list[dict]

scilpy.tractanalysis.distance_to_centroid module

scilpy.tractanalysis.distance_to_centroid.associate_labels(target_sft, min_label=1, max_label=20)[source]

Associate labels to the streamlines in a target SFT using their lengths. Even if unequal distance between points, the labels are interpolated linearly so all the points are labeled according to their position.

min and max labels are used in case there is a cut in the bundle.

Parameters:
  • target_sft (StatefulTractogram) – The target SFT to label, streamlines can be in any space.

  • min_label (int) – Minimum label to use.

  • max_label (int) – Maximum label to use.

Returns:

Array – Labels for each point along the streamlines.

Return type:

np.uint16

scilpy.tractanalysis.distance_to_centroid.closest_match_to_centroid(bundle_pts, centroid_pts, nb_pts)[source]

Assign a label to each point in the bundle_pts based on the closest centroid_pts. The labels are between 1 and nb_pts, where nb_pts is the number of points in the centroid_pts. The labels are assigned based on the order of the centroid_pts.

The 3D points are expected to be in the same space.

Typically the bundle_pts will be voxel indices (from argwhere) and the centroid_pts will be the 3D positions of a single streamline.

Parameters:
  • bundle_pts (np.array) – Coordinates of all streamlines (N x nb_pts x 3)

  • centroid_pts (np.array) – Coordinates of all streamlines (nb_pts x 3)

  • nb_pts (int) – Number of point for the association to centroids

Returns:

Labels (between 1 and nb_pts) for all bundle points

Return type:

Array

scilpy.tractanalysis.distance_to_centroid.compute_distance_map(labels_map, binary_mask, nb_pts, use_manhattan=False)[source]

Computes the distance map for each label in the labels_map.

Parameters:
  • (numpy.ndarray) (binary_mask) – A 3D array representing the labels map.

  • (numpy.ndarray) – A 3D binary map used to calculate barycenter binary map.

  • (int) (nb_pts) – Number of points to use for computing barycenters.

  • (bool) (use_manhattan) – If True, use the Manhattan distance instead of the Euclidian distance.

Returns:

numpy.ndarray

Return type:

A 3D array representing the distance map.

scilpy.tractanalysis.distance_to_centroid.compute_labels_map_barycenters(labels_map, is_euclidian=False, nb_pts=False)[source]

Compute the barycenter for each label in a 3D NumPy array by maximizing the distance to the boundary.

Parameters:
  • labels_map ((ndarray)) – The 3D array containing labels from 1-nb_pts. euclidian (bool): If True, the barycenter is the mean of the points in the mask. If False, the barycenter is the medoid of the points in the mask.

  • is_euclidian (bool) – If True, the barycenter is the mean of the points in the mask. If False, the barycenter is the medoid of the points in the mask. This is useful for the hyperplane method.

  • nb_pts (int) – Number of points to use for computing barycenters.

Returns:

An array of size (nb_pts, 3) containing the barycenter for each label.

Return type:

ndarray

scilpy.tractanalysis.distance_to_centroid.correct_labels_jump(labels_map, streamlines, nb_pts)[source]

Correct the labels jump in the labels map by cutting the streamlines where the jump is detected and keeping the longest chunk.

This avoid loops in the labels map and ensure that the labels are consistent along the streamlines.

Parameters:
  • (ndarray) (labels_map) – A 3D array representing the labels map.

  • (ArraySequence) (streamlines) – The streamlines used to compute the labels map.

  • (int) (nb_pts) – Number of points to use for computing barycenters.

Returns:

ndarray

Return type:

A 3D array representing the corrected labels map.

scilpy.tractanalysis.distance_to_centroid.find_medoid(points, max_points=10000)[source]

Find the medoid among a set of points. A medoid is a point that minimizes the sum of the distances to all other points. Unlike a barycenter, the medoid is guaranteed to be one of the points in the set.

Parameters:
  • points (ndarray) – An array of 3D coordinates.

  • max_points (int) – Maximum number of points to use for the computation (will randomly select points if the number of points is greater than max_points).

Returns:

The 3D coordinates of the medoid.

Return type:

np.array

scilpy.tractanalysis.distance_to_centroid.masked_manhattan_distance(mask, target_positions)[source]

Compute the Manhattan distance from every position in a mask to a set of positions, without stepping out of the mask.

Parameters:
  • (ndarray) (mask) – A binary 3D array representing the mask.

  • (list) (target_positions) – A list of target positions within the mask.

Returns:

A 3D array of the same shape as the mask, containing the Manhattan distances.

Return type:

ndarray

scilpy.tractanalysis.distance_to_centroid.subdivide_bundles(sft, sft_centroid, binary_mask, nb_pts, method='centerline', fix_jumps=True)[source]

Function to divide a bundle into multiple section along its length. The resulting labels map is based on the binary_mask, but the streamlines are required for a few internal corrections (for consistency).

The default is to use the euclidian/centerline method, which is fast and works well for most cases.

The hyperplane method allows for more complex shapes and to split the bundles into subsections that follow the geometry of each kind of bundle. However, this method is slower and requires extra quality control to ensure that the labels are correct. This method requires a centroid file that contains multiple streamlines.

Parameters:
  • (StatefulTractogram) (sft_centroid) – Represent the streamlines to be subdivided, streamlines representation is useful fro the fix_jump parameter.

  • (StatefulTractogram) – Centroids used as a reference for subdivision.

  • (ndarray) (binary_mask) – Mask to be converted to a label mask

  • (int) (nb_pts) – Number of subdivision along streamlines’ length

  • (str) (method) – Choice between centerline or hyperplane for subdivision

  • (bool) (fix_jumps) – Run the correction for streamlines to reduce big transition along its length.

Returns:

A 3D array representing the labels map.

Return type:

ndarray

scilpy.tractanalysis.fixel_density module

scilpy.tractanalysis.fixel_density.fixel_density(peaks, bundles, dps_key=None, max_theta=45, nbr_processes=None)[source]

Compute the fixel density map per bundle. Can use parallel processing.

Parameters:
  • peaks (np.ndarray (x, y, z, 15)) – Five principal fiber orientations for each voxel.

  • bundles (list or np.array (N)) – List of (N) paths to bundles.

  • dps_key (string, optional) – Key to the data_per_streamline to use as weight instead of the number of streamlines.

  • max_theta (int, optional) – Maximum angle between streamline and peak to be associated.

  • nbr_processes (int, optional) – The number of subprocesses to use. Default: multiprocessing.cpu_count()

Returns:

fixel_density – Density per fixel per bundle.

Return type:

np.ndarray (x, y, z, 5, N)

scilpy.tractanalysis.fixel_density.fixel_maps_to_masks(maps, abs_thr, rel_thr, norm, nb_bundles)[source]

Compute the fixel density masks from fixel density maps.

Parameters:
  • maps (np.ndarray (x, y, z, 5, N)) – Density per fixel per bundle.

  • abs_thr (float) – Value of density maps threshold to obtain density masks, in number of streamlines or streamline weighting.

  • rel_thr (float) – Value of density maps threshold to obtain density masks, as a ratio of the normalized density. Must be between 0 and 1.

  • norm (string, ["fixel", "voxel", "none"]) – Way of normalizing the density maps. If fixel, will normalize the maps per fixel, in each voxel. If voxel, will normalize the maps per voxel. If none, will not normalize the maps.

  • nb_bundles (int (N)) – Number of bundles (N).

Returns:

  • masks (np.ndarray (x, y, z, 5, N)) – Density masks per fixel per bundle.

  • maps (np.ndarray (x, y, z, 5, N)) – Normalized density maps per fixel per bundle.

scilpy.tractanalysis.json_utils module

scilpy.tractanalysis.json_utils.average_dict(curr_dict)[source]

Compute the mean and std of a metric in a json file

Parameters:

curr_dict (dict)

Return type:

dictionary with mean and std computed

scilpy.tractanalysis.json_utils.merge_dict(dict_1, dict_2, no_list=False, recursive=False)[source]

Merge two dictionary (used for tractometry)

Parameters:
  • dict_1 (dict) – dictionary to merge

  • dict_2 (dict) – dictionary to merge

  • no_list (boolean) – If False it merge list else it will create a list of list

  • recursive (boolean) –

    If true it will try to go as deep as possible in the dict

    to get to the list or value instead of updating dictionaries

Returns:

new_dict – Merged dictionary

Return type:

dict

scilpy.tractanalysis.mrds_along_streamlines module

scilpy.tractanalysis.mrds_along_streamlines.mrds_metric_sums_along_streamlines(sft, mrds_pdds, metrics, max_theta, length_weighting)[source]

Compute a sum map along a bundle for a given fixel-specific metric.

Parameters:
  • sft (StatefulTractogram) – StatefulTractogram containing the streamlines needed.

  • mrds_pdds (ndarray (X, Y, Z, 3*N_TENSORS)) – MRDS principal diffusion directions (PDDs) of the tensors

  • metrics (list of ndarray (X, Y, Z, N_TENSORS)) – Fixel-specific metrics.

  • max_theta (float) – Maximum angle in degrees between the fiber direction and the MRDS principal diffusion direction.

  • length_weighting (bool) – If True, will weight the metric values according to segment lengths.

Returns:

  • metric_sum_map (np.array) – fixel-specific metrics sum map.

  • weight_map (np.array) – Segment lengths.

scilpy.tractanalysis.mrds_along_streamlines.mrds_metrics_along_streamlines(sft, mrds_pdds, metrics, max_theta, length_weighting)[source]

Compute mean map for a given fixel-specific metric along streamlines.

Parameters:
  • sft (StatefulTractogram) – StatefulTractogram containing the streamlines needed.

  • mrds_pdds (ndarray (X, Y, Z, 3*N_TENSORS)) – MRDS principal diffusion directions of the tensors

  • metrics (list of ndarray) – Array of shape (X, Y, Z, N_TENSORS) containing the fixel-specific metric of interest.

  • max_theta (float) – Maximum angle in degrees between the fiber direction and the MRDS principal diffusion direction.

  • length_weighting (bool) – If True, will weigh the metric values according to segment lengths.

scilpy.tractanalysis.multi_bundle_operations module

scilpy.tractanalysis.multi_bundle_operations.filter_by_occurrence(sft_list, vol_dim, ratio_voxels=0.5, ratio_streamlines=0.5)[source]

Use a list of bundles (representing the same bundle in a population, for instance) to keep the voxels representing the average population bundle.

Parameters:
  • sft_list (list[StatefulTractogram]) – The bundles.

  • vol_dim (int, int, int) – The dimension of the volume.

  • ratio_voxels (float, optional) – The threshold on the ratio of bundles that must contain a voxel for it to be kept.

  • ratio_streamlines (float, optional) – The threshold on the ratio of bundles that must contain a streamline for it to be kept.

Returns:

  • volume (np.ndarray) – The volume representing the average bundle.

  • new_sft (StatefulTractogram) – The tractogram reprensenting the average bundle.

scilpy.tractanalysis.reproducibility_measures module

scilpy.tractanalysis.reproducibility_measures.approximate_surface_node(roi)[source]

Compute the number of surface voxels (i.e. nodes connected to at least one zero-valued neighboring voxel)

Parameters:

roi (ndarray) – A ndarray where voxel values represent the density of a bundle and it is binarized.

Returns:

int

Return type:

the number of surface voxels

scilpy.tractanalysis.reproducibility_measures.binary_classification(segmentation_indices, gold_standard_indices, original_count, mask_count=0)[source]

Compute all the binary classification measures using only indices from a dataset and its ground truth in any representation (voxels or streamlines).

Parameters:
  • segmentation_indices (list of int) – Indices of the data that are part of the segmentation.

  • gold_standard_indices (list of int) – Indices of the ground truth.

  • original_count (int) – Total size of the original dataset (before segmentation), e.g len(streamlines) or np.prod(mask.shape).

  • mask_count (int) – Number of non-zeros voxels in the original dataset. Needed in order to get a valid true positive count for the voxel representation.

Returns:

  • A tuple containing

  • - float (Value between 0 and 1 that represent the spatial aggrement) – between both bundles.

  • - list of ndarray (intersection_robust of streamlines in both bundle)

  • - list of ndarray (union_robust of streamlines in both bundle)

scilpy.tractanalysis.reproducibility_measures.compare_bundle_wrapper(density_1, density_2, endpoints_density_1, endpoints_density_2, bundle_1, bundle_2, centroids_1=None, centroids_2=None, voxel_size=1, ratio=False, streamline_dice=False, disable_streamline_distance=False, bundle_adjency_no_overlap=False)[source]

Compute the similarity between two bundles in the voxel representation and streamline representation.

This function conviently computes the similarity between two bundles using different metrics. The function returns a dictionary containing the computed measures.

Density and endpoints density maps, and centroids are pre-computed to speed up the computation (the script using this function does multiple comparisons with the same data).

Similar to compare_volume_wrapper (but for data from streamlines).

Parameters:
  • density_1 (ndarray) – Density map computed from the first bundle.

  • density_2 (ndarray) – Density map computed from the second bundle.

  • endpoints_density_1 (ndarray) – Density map computed from the endpoints of the first bundle.

  • endpoints_density_2 (ndarray) – Density map computed from the endpoints of the second bundle.

  • bundle_1 (list of ndarray) – First set of streamlines.

  • bundle_2 (list of ndarray) – Second set of streamlines.

  • centroids_1 (list of ndarray) – Pre-computed centroids for the first bundle.

  • centroids_2 (list of ndarray) – Pre-computed centroids for the second bundle.

  • voxel_size (float) – Size of the voxel in mm.

  • ratio (bool) – If true, the measures are normalized by the total number of voxels in the first volume.

  • streamline_dice (bool) – If true, compute the dice coefficient between the two sets of streamlines.

  • disable_streamline_distance (bool) – If true, skip the computation of the distance between streamlines

  • bundle_adjency_no_overlap (bool) – If true, exclude overlapping voxels (0mm) from the computation.

Returns:

dict

Return type:

Dictionary containing the computed measures

scilpy.tractanalysis.reproducibility_measures.compare_volume_wrapper(data_1, data_2, voxel_size=1, ratio=False, adjency_no_overlap=False, labels_to_mask=False, one_sided=False)[source]

Compute the similarity between binary mask or labels maps in the voxel representation.

This function conviently computes the similarity between two volumes using different metrics. The function returns a dictionary containing the computed measures.

Similar to compare_bundle_wrapper (but just for Nifti volumes). If coming from scil_volume_pairwise_comparison with ‘–single_compare’, data_1 is a candidate volume and data_2 is the reference.

Parameters:
  • data_1 (ndarray) – First volume to compare.

  • data_2 (ndarray) – Second volume to compare.

  • voxel_size (float) – Size of the voxel in mm.

  • ratio (bool) – If true, the measures are normalized by the total number of voxels in the first volume.

  • adjency_no_overlap (bool) – If true, exclude overlapping voxels (0mm) from the computation.

  • labels_to_mask (bool) – If true, explicitely compare every labels in the first image to the ROI in the second image. Otherwise, the computation presumes both images are either label maps or binary masks.

  • one_sided (bool) – If true, compute the measures presuming the second image is a “ground truth”.

Returns:

dict

Return type:

Dictionary containing the computed measures.

scilpy.tractanalysis.reproducibility_measures.compute_bundle_adjacency_streamlines(bundle_1, bundle_2, non_overlap=False, centroids_1=None, centroids_2=None)[source]

Compute the distance in millimeters between two bundles. Uses centroids to limit computation time. Each centroid of the first bundle is matched to the nearest centroid of the second bundle and vice-versa. Distance between matched paired is averaged for the final results. .. rubric:: References

[1] Garyfallidis et al. Robust and efficient linear

registration of white-matter fascicles in the space of streamlines, Neuroimage, 2015.

Parameters:
  • bundle_1 (list of ndarray) – First set of streamlines.

  • bundle_2 (list of ndarray) – Second set of streamlines.

  • non_overlap (bool) – Exclude overlapping streamlines from the computation.

  • centroids_1 (list of ndarray) – Pre-computed centroids for the first bundle.

  • centroids_2 (list of ndarray) – Pre-computed centroids for the second bundle.

Returns:

float

Return type:

Distance in millimeters between both bundles.

scilpy.tractanalysis.reproducibility_measures.compute_bundle_adjacency_voxel(binary_1, binary_2, non_overlap=False)[source]

Compute the distance in millimeters between two bundles in the voxel representation. Convert the bundles to binary masks. Each voxel of the first bundle is matched to the the nearest voxel of the second bundle and vice-versa. Distance between matched paired is averaged for the final results. :param binary_1: Binary mask computed from the first bundle :type binary_1: ndarray :param binary_2: Binary mask computed from the second bundle :type binary_2: ndarray :param non_overlap: Exclude overlapping voxels from the computation. :type non_overlap: bool

Returns:

float

Return type:

Distance in millimeters between both bundles.

scilpy.tractanalysis.reproducibility_measures.compute_correlation(density_1, density_2)[source]

Compute the overlap (dice coefficient) between two density maps (or binary). Correlation being less robust to extreme case (no overlap, identical array), a lot of check a needed to prevent NaN. :param density_1: Density (or binary) map computed from the first bundle :type density_1: ndarray :param density_2: Density (or binary) map computed from the second bundle :type density_2: ndarray

Returns:

float – between both bundles taking into account density.

Return type:

Value between 0 and 1 that represent the spatial aggrement

scilpy.tractanalysis.reproducibility_measures.compute_dice_streamlines(bundle_1, bundle_2)[source]

Compute the overlap (dice coefficient) between two bundles. Both bundles need to come from the exact same tractogram.

Parameters:
  • bundle_1 (list of ndarray) – First set of streamlines.

  • bundle_2 (list of ndarray) – Second set of streamlines.

Returns:

  • A tuple containing

  • - float (Value between 0 and 1 that represent the spatial aggrement) – between both bundles.

  • - list of ndarray (intersection_robust of streamlines in both bundle)

  • - list of ndarray (union_robust of streamlines in both bundle)

scilpy.tractanalysis.reproducibility_measures.compute_dice_voxel(density_1, density_2)[source]

Compute the overlap (dice coefficient) between two density maps (or binary).

Parameters:
  • density_1 (ndarray) – Density (or binary) map computed from the first bundle

  • density_2 (ndarray) – Density (or binary) map computed from the second bundle

Returns:

  • A tuple containing

  • - float (Value between 0 and 1 that represent the spatial aggrement) – between both bundles.

  • - float (Value between 0 and 1 that represent the spatial aggrement) – between both bundles, weighted by streamlines density.

scilpy.tractanalysis.reproducibility_measures.compute_fractal_dimension(density, n_steps=10, box_size_min=1.0, box_size_max=2.0, threshold=0.0, box_size=None)[source]

Compute the fractal dimension of a bundle to measure the roughness. The code is extracted from https://github.com/FBK-NILab/fractal_dimension Parameters. The result is dependent on voxel size and the number of voxels. If data comparison is performed, the bundles MUST be in same resolution.

Parameters:
  • density (ndarray) – A ndarray where voxel values represent the density of a bundle. This function computes the fractal dimension of the bundle.

  • n_steps (int) – The number of box sizes used to approximate fractal dimension. A larger number of steps will increase the accuracy of the approximation, but will also take more time. The default number of boxes sizes is 10.

  • box_size_min (float) – The minimum size of boxes. This number should be larger than or equal to 1.0 and is defaulted to 1.0.

  • box_size_max (float) – The maximum size of boxes. This number should be larger than the minimum size of boxes.

  • threshold (float) – The threshold to filter the voxels in the density array. The default is set to 0, so only nonzero voxels will be considered.

  • box_size (ndarray) – A ndarray of different sizes of boxes in a linear space in an ascending order.

Returns:

float

Return type:

fractal dimension of a bundle

scilpy.tractanalysis.reproducibility_measures.compute_hausdorff_voxel(map_1, map_2)[source]

Compute the Hausdorff distance between two sets of coordinates, representing density or binary maps.

Parameters:
  • map_1 (ndarray) – Density (or binary) map computed from the first bundle

  • map_2 (ndarray) – Density (or binary) map computed from the second bundle

Returns:

distance – Value representing the Hausdorff distance in voxels. Should be multiplied afterwards by the voxel size to get the distance in mm.

Return type:

float

scilpy.tractanalysis.reproducibility_measures.tractogram_pairwise_comparison(sft_one, sft_two, mask, nbr_cpu=1, skip_streamlines_distance=True)[source]

Compute the difference between two sets of streamlines for each voxel in the mask. This function uses multiprocessing to compute the difference between two sets of streamlines for each voxel.

Parameters:
  • sft_one (StatefulTractogram) – First tractogram to compare.

  • sft_two (StatefulTractogram) – Second tractogram to compare.

  • mask (np.ndarray) – Mask of the data to compare (optional).

  • nbr_cpu (int) – Number of CPU to use (default: 1).

  • skip_streamlines_distance (bool) – If true, skip the computation of the distance between streamlines. (default: True)

Returns:

  • acc_norm (np.ndarray) – Angular correlation coefficient.

  • corr_norm (np.ndarray) – Correlation coefficient of density maps.

  • diff_norm (np.ndarray) – Voxelwise distance between sets of streamlines.

  • heatmap (np.ndarray) – Merged heatmap of the three metrics using harmonic mean.

  • mask (np.ndarray) – Final mask. Intersection of given mask (if any) and density masks of both tractograms.

scilpy.tractanalysis.scoring module

Tractometry

Global connectivity metrics:

  • Computed by default:
    • VS: valid streamlines, belonging to a bundle (i.e. respecting all the

      criteria for that bundle; endpoints, limit_mask, gt_mask.).

    • IS: invalid streamlines. All other streamlines. IS = IC + NC.

  • Optional:
    • WPC: wrong path connections, streamlines connecting correct ROIs but not

      respecting the other criteria for that bundle. Such streamlines always exist but they are only saved separately if specified in the options. Else, they are merged back with the IS. By definition. WPC are only computed if “limits masks” are provided.

    • IC: invalid connections, streamlines joining an incorrect combination of

      ROIs. Use carefully, quality depends on the quality of your ROIs and no analysis is done on the shape of the streamlines.

    • NC: no connections. Invalid streamlines minus invalid connections.

  • Fidelity metrics:
    • OL: Overlap. Percentage of ground truth voxels containing streamline(s)

      for a given bundle.

    • OR: Overreach. Amount of voxels containing streamline(s) when they

      shouldn’t, for a given bundle. We compute two versions : OR_pct_vs = divided by the total number of voxel covered by the bundle. (percentage of the voxels touched by VS). Values range between 0 and 100%. Values are not defined when we recovered no streamline for a bundle, but we set the OR_pct_vs to 0 in that case. OR_pct_gt = divided by the total size of the ground truth bundle mask. Values could be higher than 100%.

    • f1 score: which is the same as the Dice score.

scilpy.tractanalysis.scoring.compute_f1_overlap_overreach(current_vb_voxels, gt_mask, dimensions)[source]

Compute f1, OL and OR/ORn based on a ground truth mask.

Parameters:
  • current_vb_voxels (3D array) – The voxels touched by at least one streamlines for a given bundle.

  • gt_mask (3D array) – The ground truth mask.

  • dimensions (np.array) – The nibabel dimensions of the data (3D).

Returns:

  • f1 (float) – The f1 score.

  • tp_nb_voxels (int) – The TP (true positive) count in number of voxels.

  • fp_nb_voxels (int) – The FP (false positive) count in number of voxels. Hint: Divide it by the ground truth count to get the overreach, or by the recovered bundle count to get the ORn (scores used in the ismrm2015 tractography challenge).

  • fn_nb_voxels (int) – The number of voxels from the gt_mask that have not been recovered; corresponds to the FN count (false negative).

  • overlap (float) – TP divided by the ground truth count (i.e. TP + FN), in percentage.

  • overreach_pct_gt (float) – The overreach, normalized by the ground truth area.

  • overreach_pct_vs (float) – The overreach, normalized by the recovered bundle’s area. (Or 0 if no streamline have been recovered for this bundle).

scilpy.tractanalysis.scoring.compute_f1_score(overlap, overreach)[source]

Compute the F1 score between overlap and overreach (they must be percentages).

Parameters:
  • overlap (float, The overlap value.)

  • overreach (float, The overreach value) – (Version normalized over bundle area, not version normalized over gt).

Returns:

  • f1_score (float, The f1 score.)

  • Ref (https://en.wikipedia.org/wiki/F1_score)

scilpy.tractanalysis.scoring.compute_tractometry(vb_sft_list, wpc_sft_list, ib_sft_list, nc_sft, args, bundles_names, gt_masks, dimensions, ib_names)[source]

Tractometry stats: First in terms of connections (NC, IC, VS, WPC), then in terms of volume (OL, OR, Dice score)

scilpy.tractanalysis.scoring.get_binary_maps(sft)[source]

Extract a mask from a bundle.

Parameters:

sft (StatefulTractogram) – Bundle.

Returns:

  • bundles_voxels (numpy.ndarray) – Mask representing the bundle volume.

  • endpoints_voxels (numpy.ndarray) – Mask representing the bundle’s endpoints.

scilpy.tractanalysis.streamlines_metrics module

scilpy.tractanalysis.streamlines_metrics.compute_tract_counts_map(streamlines, vol_dims)

scilpy.tractanalysis.todi module

class scilpy.tractanalysis.todi.TrackOrientationDensityImaging(img_shape, sphere_type='repulsion724')[source]

Bases: object

compute_average_dir()[source]

Voxel-wise average of TODI orientations.

Average all orientation of each voxel, of the TODI map, into a single direction, warning for crossing.

Returns:

avg_dir – Volume containing a single 3-vector (peak) per voxel.

Return type:

numpy.ndarray (4D)

compute_distance_to_peak(peak_img, normalize_count=True, deg=True, with_avg_dir=True)[source]

Compute distance to peak map.

Compute the distance of the TODI map to peaks at each position, in radian or degree.

Parameters:
  • peak_img (numpy.ndarray (4D)) – Peaks image, as written by most of Scilpy scripts.

  • normalize_count (bool, optional) – Normalize/weight the error map by the density map (default True).

  • deg (bool, optional) – Returned error map as degree instead of radian (default True).

  • with_avg_dir (bool, optional) – Average all orientation of each voxel of the TODI map into a single direction, warning for crossing (default True).

Returns:

error_map – Average angle error map per voxel.

Return type:

numpy.ndarray (3D)

compute_todi(streamlines, length_weights=True, n_steps=1, asymmetric=False)[source]

Compute the TODI map.

At each voxel an histogram distribution of the local streamlines orientations (TODI) is computed.

Parameters:
  • streamlines (list of numpy.ndarray) – List of streamlines.

  • length_weights (bool, optional) – Weights TODI map of each segment’s length (default True).

get_mask()[source]
get_sh(sh_basis, sh_order, smooth=0.006, full_basis=False, is_legacy=True)[source]

Spherical Harmonics (SH) coefficients of the TODI map

Compute the SH representation of the TODI map, converting SF to SH with a smoothing factor.

Parameters:
  • sh_basis ({None, 'tournier07', 'descoteaux07'}) – None for the default DIPY basis, tournier07 for the Tournier 2007 [2] basis, and descoteaux07 for the Descoteaux 2007 [1] basis (None defaults to descoteaux07).

  • sh_order (int) – Maximum SH order in the SH fit. For sh_order, there will be (sh_order + 1) * (sh_order_2) / 2 SH coefficients (default 4).

  • smooth (float, optional) – Smoothing factor for the conversion, Lambda-regularization in the SH fit (default 0.006).

  • is_legacy (bool, optional) – Whether or not the SH basis is in its legacy form.

Returns:

todi_sh – SH representation of the TODI map

Return type:

ndarray

References

get_tdi()[source]

Compute the TDI map.

Compute the Tract Density Image (TDI) from the TODI volume.

Returns:

tdi – Tract Density Image

Return type:

numpy.ndarray (3D)

get_todi()[source]
get_todi_shape()[source]
mask_todi(mask)[source]

Mask the TODI map.

Mask the TODI without having to reshape the whole volume all at once (big in memory).

Parameters:

mask (numpy.ndarray) – Given volume mask for the TODI map.

normalize_todi_per_voxel(p_norm=2)[source]

Normalize TODI map.

Normalize TODI distribution on the sphere for each voxel independently.

Parameters:

p_norm (int, optional) – Chosen Norm to normalize.

Returns:

todi – Normalized TODI map.

Return type:

numpy.ndarray

reshape_to_3d(img_voxelly_masked)[source]

Reshape a complex ravelled image to 3D.

Unravel a given unravel mask (1D), image (1D), SH/SF (2D) to its original 3D shape (with a 4D for SH/SF).

Parameters:

img_voxelly_masked (numpy.ndarray (either in 1D, 2D or 3D)) – That will be reshaped to 3D (input data shape) accordingly. Necessary for future masking operation.

Returns:

unraveled_img – Unravel volume in x, y, z (, c).

Return type:

numpy.ndarray (3D, or 4D)

set_todi(mask, todi)[source]
smooth_todi_dir(order=2)[source]

Smooth orientations on the sphere.

Smooth the orientations / distribution on the sphere. Important for priors construction of BST.

Parameters:

order (int, optional) – Exponent blurring factor, based on the dot product (default 2).

smooth_todi_spatial(sigma=0.5)[source]

Spatial Smoothing of the TODI map.

Blur the TODI map using neighborhood information. Important for priors construction of BST.

Parameters:

sigma (float, optional) – Gaussian blurring factor (default 0.5).

scilpy.tractanalysis.todi.get_sf_from_todi(sft, mask, todi_sigma, sf_threshold, sphere='repulsion724')[source]

Track Orientation Density Imaging (TODI) [1], provides information about the voxel-wise orientation distribution of streamlines. This can be represented by a sphere for each voxel where the radii are weigthed to represent how aligned streamlines in each specific direction.

Here the sphere is represented as spherical harmonics. See also get_sf_from_todi.

[1] Dhollander, Thijs, et al. “Track orientation density imaging (TODI) and track orientation distribution (TOD) based tractography.” NeuroImage 94 (2014): 312-336.

Parameters:
  • sft (StatefulTractogram) – The tractogram

  • mask (np.ndarray) – A binary mask.

  • todi_sigma (int) – Smooth the orientation histogram. A value between 0 and 5.

  • sf_threshold (float) – Relative threshold for sf masking (0.0-1.0).

  • sphere (str) – The name of the sphere.

Returns:

  • sf_data (np.ndarray) – The track orientation distribution per voxel. Shape: n x s, where n in the number of voxels in the final mask (sub_mask_3d), and s is the number of points on the sphere.

  • sub_mask_3d (np.ndarray) – The final mask, from both the input mask and sf_threshold.

scilpy.tractanalysis.todi.get_sh_from_todi(sft, mask)[source]

Track Orientation Density Imaging (TODI) [1], provides information about the voxel-wise orientation distribution of streamlines. This can be represented by a sphere for each voxel where the radii are weigthed to represent how aligned streamlines in each specific direction.

Here the sphere is represented as spherical harmonics. See also get_sf_from_todi. See also dipy.reconst.shm.sf_to_sh, sh_to_sf.

[1] Dhollander, Thijs, et al. “Track orientation density imaging (TODI) and track orientation distribution (TOD) based tractography.” NeuroImage 94 (2014): 312-336.

Parameters:
  • sft (StatefulTractogram) – The tractogram

  • mask (np.ndarray) – A binary mask.

Returns:

sh_data – The track orientation distribution per voxel.

Return type:

np.ndarray

scilpy.tractanalysis.todi_util module

scilpy.tractanalysis.todi_util.compute_vectors_norm(vectors)[source]
scilpy.tractanalysis.todi_util.generate_mask_indices_1d(nb_voxel, indices_1d)[source]
scilpy.tractanalysis.todi_util.get_dir_to_sphere_id(vectors, sphere_vertices)[source]
Find the closest vector on the sphere vertices using a cKDT tree

sphere_vertices must be normed (or all with equal norm).

Parameters:
  • vectors (numpy.ndarray (2D)) – Vectors representing the direction (x,y,z) of segments.

  • sphere_vertices (numpy.ndarray (2D)) – Vertices of a Dipy sphere object.

Returns:

dir_sphere_id – Sphere indices of the closest sphere direction for each vector

Return type:

numpy.ndarray (1D)

scilpy.tractanalysis.todi_util.get_indices_1d(volume_shape, pts)[source]
scilpy.tractanalysis.todi_util.get_segments_dir_and_norm(segments, seg_mid=None, asymmetric=False)[source]
scilpy.tractanalysis.todi_util.get_segments_mid_pts_positions(segments)[source]
scilpy.tractanalysis.todi_util.get_segments_vectors(segments)[source]
scilpy.tractanalysis.todi_util.get_vectors_dir_and_norm(vectors)[source]
scilpy.tractanalysis.todi_util.get_vectors_dir_and_norm_rel_to_center(vectors, seg_mid_pts)[source]

Evaluates vectors direction and norm by taking into account the orientation and position of segments in relation to the center of voxel

scilpy.tractanalysis.todi_util.normalize_vectors(vectors)[source]
scilpy.tractanalysis.todi_util.p_normalize_vectors(vectors, p)[source]
scilpy.tractanalysis.todi_util.psf_from_sphere(sphere_vertices)[source]
scilpy.tractanalysis.todi_util.streamlines_to_endpoints(streamlines)[source]

Equivalent to streamlines resampling to 2 points (first and last).

Parameters:

streamlines (list of numpy.ndarray) – List of streamlines.

Returns:

endpoints – Endpoint array representation with the first and last points.

Return type:

numpy.ndarray (2D)

scilpy.tractanalysis.todi_util.streamlines_to_pts_dir_norm(streamlines, n_steps=1, asymmetric=False)[source]

Evaluate each segment: mid position, direction, length.

Parameters:

streamlines (list of numpy.ndarray) – List of streamlines.

Returns:

  • seg_mid (numpy.ndarray (2D)) – Mid position (x,y,z) of all streamlines’ segments.

  • seg_dir (numpy.ndarray (2D)) – Direction (x,y,z) of all streamlines’ segments.

  • seg_norm (numpy.ndarray (2D)) – Length of all streamlines’ segments.

scilpy.tractanalysis.todi_util.streamlines_to_segments(streamlines, n_steps=1)[source]

Split streamlines into its segments.

Parameters:

streamlines (list of numpy.ndarray) – List of streamlines.

Returns:

segments – Segments array representation with the first and last points.

Return type:

numpy.ndarray (2D)

scilpy.tractanalysis.voxel_boundary_intersection module

scilpy.tractanalysis.voxel_boundary_intersection.subdivide_streamlines_at_voxel_faces(streamlines)

Cut streamlines segments into smaller segments such that a segment covering more than one voxel is split into smaller segments that either end or start at voxel boundaries.

Parameters:

streamlines (list of ndarray) – Streamlines coordinates in voxel space, corner origin.

Returns:

split_coordinates – Updated streamline coordinates with added coordinate points at voxel boundaries.

Return type:

list of ndarray