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. |
peak_values
-
dipy.stats.analysis.peak_values(bundle, peaks, dt, pname, bname, subject, group_id, ind, dir_name)
- 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)
- 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)
Calculates assignment maps of the target bundle with reference to
model bundle centroids.
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. (Default 100)
References
[Chandio2020]
Chandio, B.Q., Risacher, S.L., Pestilli, F., Bullock, D.,
Yeh, FC., Koudoro, S., Rokem, A., Harezlak, J., and Garyfallidis, E.
Bundle analytics, a computational framework for investigating the
shapes and profiles of brain pathways across populations.
Sci Rep 10, 17149 (2020)
gaussian_weights
-
dipy.stats.analysis.gaussian_weights(bundle, n_points=100, return_mahalnobis=False, stat=<function mean>)
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. Default: 100.
- return_mahalanobisbool, optional
Whether to return the Mahalanobis distance instead of the weights.
Default: False.
- 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. Default: np.mean
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)
Calculates a summarized profile of data for a bundle or tract
along its length.
Follows the approach outlined in [Yeatman2012].
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
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
[Yeatman2012]
Yeatman, Jason D., Robert F. Dougherty,
Nathaniel J. Myall, Brian A. Wandell, and Heidi M. Feldman. 2012.
“Tract Profiles of White Matter Properties: Automating Fiber-Tract
Quantification” PloS One 7 (11): e49790.