dmriprep.utils.vectors module¶
Utilities to operate on diffusion gradients.
-
class
dmriprep.utils.vectors.
DiffusionGradientTable
(b0_threshold=50, b_scale=True, bvals=None, bvec_norm_epsilon=0.1, bvecs=None, dwi_file=None, raise_inconsistent=False, rasb_file=None, transforms=None)¶ Bases:
object
Data structure for DWI gradients.
-
property
affine
¶ Get the affine for RAS+/image-coordinates conversions.
-
property
b0mask
¶ Get a mask of low-b frames.
-
property
bvals
¶ Get the N b-values.
-
property
bvecs
¶ Get the N x 3 list of bvecs.
-
generate_rasb
()¶ Compose RAS+B gradient table.
-
generate_vecval
()¶ Compose a bvec/bval pair in image coordinates.
-
property
gradients
¶ Get gradient table (rasb).
-
normalize
()¶ Normalize (l2-norm) b-vectors.
-
property
normalized
¶ Return whether b-vecs have been normalized.
-
property
pole
¶ Check whether the b-vectors cover a full or just half a shell.
If pole is all-zeros then the b-vectors cover a full sphere.
-
reorient_rasb
()¶ Reorient the vectors based o a list of affine transforms.
-
to_filename
(filename, filetype='rasb')¶ Write files (RASB, bvecs/bvals) to a given path.
-
property
-
dmriprep.utils.vectors.
b0mask_from_data
(dwi_file, mask_file, z_thres=3.0)¶ Evaluate B0 locations relative to mean signal variation.
Standardizes (z-score) the average DWI signal within mask and threshold. This is a data-driven way of estimating which volumes in the DWI dataset are really encoding low-b acquisitions.
-
dmriprep.utils.vectors.
bvecs2ras
(affine, bvecs, norm=True, bvec_norm_epsilon=0.2)¶ Convert b-vectors given in image coordinates to RAS+.
Examples
>>> bvecs2ras(2.0 * np.eye(3), [(1.0, 0.0, 0.0), (-1.0, 0.0, 0.0)]).tolist() [[1.0, 0.0, 0.0], [-1.0, 0.0, 0.0]]
>>> affine = np.eye(4) >>> affine[0, 0] *= -1.0 # Make it LAS >>> bvecs2ras(affine, [(1.0, 0.0, 0.0), (-1.0, 0.0, 0.0)]).tolist() [[-1.0, 0.0, 0.0], [1.0, 0.0, 0.0]]
>>> affine = np.eye(3) >>> affine[:2, :2] *= -1.0 # Make it LPS >>> bvecs2ras(affine, [(1.0, 0.0, 0.0), (-1.0, 0.0, 0.0)]).tolist() [[-1.0, 0.0, 0.0], [1.0, 0.0, 0.0]]
>>> affine[:2, :2] = [[0.0, 1.0], [1.0, 0.0]] # Make it ARS >>> bvecs2ras(affine, [(1.0, 0.0, 0.0), (-1.0, 0.0, 0.0)]).tolist() [[0.0, 1.0, 0.0], [0.0, -1.0, 0.0]]
>>> bvecs2ras(affine, [(0.0, 0.0, 0.0)]).tolist() [[0.0, 0.0, 0.0]]
>>> bvecs2ras(2.0 * np.eye(3), [(1.0, 0.0, 0.0), (-1.0, 0.0, 0.0)], ... norm=False).tolist() [[2.0, 0.0, 0.0], [-2.0, 0.0, 0.0]]
-
dmriprep.utils.vectors.
calculate_pole
(bvecs, bvec_norm_epsilon=0.1)¶ Check whether the b-vecs cover a hemisphere, and if so, calculate the pole.
- Parameters
bvecs (numpy.ndarray) – 2D numpy array with shape (N, 3) where N is the number of points. All points must lie on the unit sphere.
- Returns
pole – A zero-vector if
bvecs
covers the full sphere, and the unit vector locating the hemisphere pole othewise.- Return type
numpy.ndarray
Examples
>>> bvecs = [(1.0, 0.0, 0.0), (-1.0, 0.0, 0.0), ... (0.0, 1.0, 0.0), (0.0, -1.0, 0.0), ... (0.0, 0.0, 1.0)] # Just half a sphere >>> calculate_pole(bvecs).tolist() [0.0, 0.0, 1.0]
>>> bvecs += [(0.0, 0.0, -1.0)] # Make it a full-sphere >>> calculate_pole(bvecs).tolist() [0.0, 0.0, 0.0]
References
https://rstudio-pubs-static.s3.amazonaws.com/27121_a22e51b47c544980bad594d5e0bb2d04.html
-
dmriprep.utils.vectors.
normalize_gradients
(bvecs, bvals, b0_threshold=50, bvec_norm_epsilon=0.1, b_scale=True, raise_error=False)¶ Normalize b-vectors and b-values.
The resulting b-vectors will be of unit length for the non-zero b-values. The resultinb b-values will be normalized by the square of the corresponding vector amplitude.
- Parameters
bvecs (m x n 2d array) – Raw b-vectors array.
bvals (1d array) – Raw b-values float array.
b0_threshold (float) – Gradient threshold below which volumes and vectors are considered B0’s.
- Returns
bvecs (m x n 2d array) – Unit-normed b-vectors array.
bvals (1d int array) – Vector amplitude square normed b-values array.
Examples
>>> bvecs = np.vstack((np.zeros(3), 2.0 * np.eye(3), -0.8 * np.eye(3), np.ones(3))) >>> bvals = np.array([1000] * bvecs.shape[0]) >>> normalize_gradients(bvecs, bvals, 50, raise_error=True) Traceback (most recent call last): ValueError:
>>> bvals[0] = 0.0 >>> norm_vecs, norm_vals = normalize_gradients(bvecs, bvals) >>> np.all(norm_vecs[0] == 0) True
>>> norm_vecs[1, ...].tolist() [1.0, 0.0, 0.0]
>>> norm_vals[0] 0 >>> norm_vals[1] 4000 >>> norm_vals[-2] 600 >>> norm_vals[-1] 3000
>>> norm_vecs, norm_vals = normalize_gradients(bvecs, bvals, b_scale=False) >>> norm_vals[0] 0 >>> np.all(norm_vals[1:] == 1000) True
-
dmriprep.utils.vectors.
rasb_dwi_length_check
(dwi_file, rasb_file)¶ Check the number of encoding vectors and number of orientations in the DWI file.