scilpy.tractanalysis package

scilpy.tractanalysis.afd_along_streamlines module

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

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

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

  • fodf (nibabel.image) – fODF with shape (X, Y, Z, #coeffs). #coeffs depend 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, optional) – Whether or not the SH basis is in its legacy form.

Returns:

  • afd_sum_map (np.array) – AFD map.

  • rd_sum_map (np.array) – fdAFD map.

  • weight_map (np.array) – Segment lengths.

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) and mean Radial fODF (radfODF) maps along a bundle.

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

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.

  • 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.distance_to_centroid module

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

Compute minimal distance to centroids

Parameters:
  • bundles_pts (np.array)

  • centroid_pts (np.array)

  • nb_pts (int)

Return type:

Array

scilpy.tractanalysis.features module

scilpy.tractanalysis.grid_intersections module

scilpy.tractanalysis.grid_intersections.grid_intersections()[source]

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.quick_tools module

scilpy.tractanalysis.quick_tools.get_next_real_point()[source]
scilpy.tractanalysis.quick_tools.get_previous_real_point()[source]

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.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

[Garyfallidis15]

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.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: Angular correlation coefficient. corr_norm: Correlation coefficient of density maps. diff_norm: Voxelwise distance between sets of streamlines. heatmap: Merged heatmap of the three metrics using harmonic mean.

Return type:

List of np.ndarray

scilpy.tractanalysis.scoring module

scilpy.tractanalysis.streamlines_metrics module

scilpy.tractanalysis.streamlines_metrics.compute_tract_counts_map()[source]

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.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.tools module

scilpy.tractanalysis.tools.compute_connectivity(indices, atlas_data, real_labels, segmenting_func)[source]
Parameters:
  • indices (ArraySequence) – The list of 3D indices [i, j, k] of all voxels traversed by all streamlines. This is the output of our uncompress function.

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

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

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

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 ending 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.tools.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]