stats#

Module: stats.analysis#

peak_values(bundle, peaks, dt, pname, bname, ...)

Peak_values function finds the generalized fractional anisotropy (gfa)

anatomical_measures(bundle, metric, dt, ...)

Calculates dti measure (eg: FA, MD) per point on streamlines and

assignment_map(target_bundle, model_bundle, ...)

Calculates assignment maps of the target bundle with reference to model bundle centroids.

gaussian_weights(bundle, *[, n_points, ...])

Calculate weights for each streamline/node in a bundle, based on a Mahalanobis distance from the core the bundle, at that node (mean, per default).

afq_profile(data, bundle, affine, *[, ...])

Calculates a summarized profile of data for a bundle or tract along its length.

Module: stats.qc#

find_qspace_neighbors(gtab)

Create a mapping of dwi volume index to its nearest neighbor.

neighboring_dwi_correlation(dwi_data, gtab, *)

Calculate the Neighboring DWI Correlation (NDC) from dMRI data.

Module: stats.sketching#

count_sketch(matrixa_name, matrixa_dtype, ...)

Count Sketching algorithm to reduce the size of the matrix.

peak_values#

dipy.stats.analysis.peak_values(bundle, peaks, dt, pname, bname, subject, group_id, ind, dir_name)[source]#
Peak_values function finds the generalized fractional anisotropy (gfa)

and quantitative anisotropy (qa) values from peaks object (eg: csa) for every point on a streamline used while tracking and saves it in hd5 file.

Parameters:
bundlestring

Name of bundle being analyzed

peakspeaks

contains peak directions and values

dtDataFrame

DataFrame to be populated

pnamestring

Name of the dti metric

bnamestring

Name of bundle being analyzed.

subjectstring

subject number as a string (e.g. 10001)

group_idinteger

which group subject belongs to 1 patient and 0 for control

indinteger list

ind tells which disk number a point belong.

dir_namestring

path of output directory

anatomical_measures#

dipy.stats.analysis.anatomical_measures(bundle, metric, dt, pname, bname, subject, group_id, ind, dir_name)[source]#
Calculates dti measure (eg: FA, MD) per point on streamlines and

save it in hd5 file.

Parameters:
bundlestring

Name of bundle being analyzed

metricmatrix of float values

dti metric e.g. FA, MD

dtDataFrame

DataFrame to be populated

pnamestring

Name of the dti metric

bnamestring

Name of bundle being analyzed.

subjectstring

subject number as a string (e.g. 10001)

group_idinteger

which group subject belongs to 1 for patient and 0 control

indinteger list

ind tells which disk number a point belong.

dir_namestring

path of output directory

assignment_map#

dipy.stats.analysis.assignment_map(target_bundle, model_bundle, no_disks)[source]#

Calculates assignment maps of the target bundle with reference to model bundle centroids.

See [1] for further details about the method.

Parameters:
target_bundlestreamlines

target bundle extracted from subject data in common space

model_bundlestreamlines

atlas bundle used as reference

no_disksinteger, optional

Number of disks used for dividing bundle into disks.

Returns:
indxndarray

Assignment map of the target bundle streamline point indices to the model bundle centroid points.

References

gaussian_weights#

dipy.stats.analysis.gaussian_weights(bundle, *, n_points=100, return_mahalnobis=False, stat=<function mean>)[source]#

Calculate weights for each streamline/node in a bundle, based on a Mahalanobis distance from the core the bundle, at that node (mean, per default).

Parameters:
bundleStreamlines

The streamlines to weight.

n_pointsint, optional

The number of points to resample to. If the `bundle` is an array, this input is ignored.

return_mahalanobisbool, optional

Whether to return the Mahalanobis distance instead of the weights.

statcallable, optional.

The statistic used to calculate the central tendency of streamlines in each node. Can be one of {np.mean, np.median} or other functions that have similar API.`

Returns:
warray of shape (n_streamlines, n_points)

Weights for each node in each streamline, calculated as its relative inverse of the Mahalanobis distance, relative to the distribution of coordinates at that node position across streamlines.

afq_profile#

dipy.stats.analysis.afq_profile(data, bundle, affine, *, n_points=100, profile_stat=<function average>, orient_by=None, weights=None, **weights_kwarg)[source]#

Calculates a summarized profile of data for a bundle or tract along its length.

Follows the approach outlined in [2].

Parameters:
data3D volume

The statistic to sample with the streamlines.

bundleStreamLines class instance
The collection of streamlines (possibly already resampled into an array

for each to have the same length) with which we are resampling. See Note below about orienting the streamlines.

affinearray_like (4, 4)

The mapping from voxel coordinates to streamline points. The voxel_to_rasmm matrix, typically from a NIFTI file.

n_points: int, optional

The number of points to sample along the bundle. Default: 100.

orient_by: streamline, optional

A streamline to use as a standard to orient all of the streamlines in the bundle according to.

weights1D array or 2D array or callable, optional

Weight each streamline (1D) or each node (2D) when calculating the tract-profiles. Must sum to 1 across streamlines (in each node if relevant). If callable, this is a function that calculates weights.

profile_statcallable, optional

The statistic used to average the profile across streamlines. If weights is not None, this must take weights as a keyword argument. The default, np.average, is the same as np.mean but takes weights as a keyword argument.

weights_kwargkey-word arguments

Additional key-word arguments to pass to the weight-calculating function. Only to be used if weights is a callable.

Returns:
ndarraya 1D array with the profile of data along the length of

bundle

Notes

Before providing a bundle as input to this function, you will need to make sure that the streamlines in the bundle are all oriented in the same orientation relative to the bundle (use orient_by_streamline()).

References

find_qspace_neighbors#

dipy.stats.qc.find_qspace_neighbors(gtab)[source]#

Create a mapping of dwi volume index to its nearest neighbor.

An approximate q-space is used (the deltas are not included). Note that neighborhood is not necessarily bijective. One neighbor is found per dwi volume.

Parameters:
gtab: dipy.core.gradients.GradientTable

Gradient table.

Returns:
neighbors: list of tuple

A list of 2-tuples indicating the nearest q-space neighbor of each dwi volume.

Examples

>>> from dipy.core.gradients import gradient_table
>>> import numpy as np
>>> gtab = gradient_table(
...     np.array([0, 1000, 1000, 2000]),
...     bvecs=np.array([
...         [1, 0, 0],
...         [1, 0, 0],
...         [0.99, 0.0001, 0.0001],
...         [1, 0, 0]]))
>>> find_qspace_neighbors(gtab)
[(1, 2), (2, 1), (3, 1)]

neighboring_dwi_correlation#

dipy.stats.qc.neighboring_dwi_correlation(dwi_data, gtab, *, mask=None)[source]#

Calculate the Neighboring DWI Correlation (NDC) from dMRI data.

Using a mask is highly recommended, otherwise the FOV will influence the correlations. According to Yeh et al.[3], an NDC less than 0.4 indicates a low quality image.

Parameters:
dwi_data4D ndarray

dwi data on which to calculate NDC

gtabdipy.core.gradients.GradientTable

Gradient table.

mask3D ndarray, optional

Mask of voxels to include in the NDC calculation

Returns:
ndcfloat

The neighboring DWI correlation

References

count_sketch#

dipy.stats.sketching.count_sketch(matrixa_name, matrixa_dtype, matrixa_shape, sketch_rows, tmp_dir)[source]#

Count Sketching algorithm to reduce the size of the matrix.

Parameters:
matrixa_namestr

The name of the memmap file containing the matrix A.

matrixa_dtypedtype

The dtype of the matrix A.

matrixa_shapetuple

The shape of the matrix A.

sketch_rowsint

The number of rows in the sketch matrix.

tmp_dirstr

The directory to save the temporary files.

Returns:
matrixc_file.namestr

The name of the memmap file containing the sketch matrix.

matrixc.dtypedtype

The dtype of the sketch matrix.

matrixc.shapetuple

The shape of the sketch matrix.