scilpy.gradients package
scilpy.gradients.bvec_bval_tools module
- class scilpy.gradients.bvec_bval_tools.B0ExtractionStrategy(value)[source]
Bases:
Enum
An enumeration.
- ALL = 'all'
- FIRST = 'first'
- MEAN = 'mean'
- scilpy.gradients.bvec_bval_tools.check_b0_threshold(min_bval, b0_thr, skip_b0_check, overwrite_with_min=True)[source]
Check if the minimal bvalue is under the threshold. If not, raise an error to ask user to update the b0_thr.
Also verifies if the b0_thr is suspicious (should be included in range [0, DEFAULT_B0_THRESHOLD]).
- Parameters:
min_bval (float) – Minimum bvalue.
b0_thr (float) – Maximum bvalue considered as a b0.
skip_b0_check (bool) – If True, and no b0 is found, only print a warning, do not raise an error.
overwrite_with_min (bool, optional) – If True, no b0 is found, and skip_b0_check is True, the script will proceed with a new b0 threshold equal to the minimal b-value. Else, the script will proceed with the original b0 threshold.
- Returns:
b0_thr – Either the unmodified b0_thr, or, in the case where the minimal b-value is larger than b0_thr, and skip_b0_check is set to True, then returns min_bval.
- Return type:
float
- Raises:
ValueError – If the minimal bvalue is over the threshold (and skip_b0_check is False).
- scilpy.gradients.bvec_bval_tools.flip_gradient_sampling(bvecs, axes, sampling_type)[source]
Flip bvecs on chosen axis.
- Parameters:
bvecs (np.ndarray) – bvecs loaded directly, not re-formatted. Careful! Must respect the format (not verified here).
axes (list of int) – List of axes to flip (e.g. [0, 1]). See str_to_axis_index.
sampling_type (str) –
- Either ‘mrtrix’: bvecs are of shape (N, 4) or
’fsl’: bvecs are of shape (3, N)
- Returns:
bvecs – The final bvecs.
- Return type:
np.array
- scilpy.gradients.bvec_bval_tools.identify_shells(bvals, tol=40.0, round_centroids=False, sort=False)[source]
Guessing the shells from the b-values. Returns the list of shells and, for each b-value, the associated shell.
Starting from the first shell as holding the first b-value in bvals, the next b-value is considered on the same shell if it is closer than threshold, or else we consider that it is on another shell. This is an alternative to K-means considering we don’t already know the number of shells K.
Note. This function should be added in Dipy soon.
- Parameters:
bvals (array (N,)) – Array of bvals
tol (float) – Limit difference to centroid to consider that a b-value is on an existing shell. On or above this limit, the b-value is placed on a new shell.
round_centroids (bool) – If true will round shell values to the nearest 10.
sort (bool) – Sort centroids and shell_indices associated.
- Returns:
centroids (array (K)) – Array of centroids. Each centroid is a b-value representing the shell. K is the number of identified shells.
shell_indices (array (N,)) – For each bval, the associated centroid K.
- scilpy.gradients.bvec_bval_tools.is_normalized_bvecs(bvecs)[source]
Check if b-vectors are normalized.
- Parameters:
bvecs ((N, 3) array) – input b-vectors (N, 3) array
- Return type:
True/False
- scilpy.gradients.bvec_bval_tools.normalize_bvecs(bvecs)[source]
Normalize b-vectors
- Parameters:
bvecs ((N, 3) array) – input b-vectors (N, 3) array
- Returns:
bvecs – normalized b-vectors
- Return type:
(N, 3)
- scilpy.gradients.bvec_bval_tools.round_bvals_to_shell(bvals, shells_to_extract, tol=20)[source]
Return bvals equal to a list of chosen bvals, up to a tolerance.
- Parameters:
bvals (np.array) – All the b-values.
shells_to_extract (list) – The shells of interest.
tol (float, optional) – The tolerance
- scilpy.gradients.bvec_bval_tools.str_to_axis_index(axis)[source]
Convert x y z axis string to 0 1 2 axis index
- Parameters:
axis (str) – Axis value (x, y or z)
- Returns:
index – Axis index
- Return type:
int or None
- scilpy.gradients.bvec_bval_tools.swap_gradient_axis(bvecs, final_order, sampling_type)[source]
Swap bvecs.
- Parameters:
bvecs (np.array) – bvecs loaded directly, not re-formatted. Careful! Must respect the format (not verified here).
final_order (new order) – Final order (ex, 2 1 0)
sampling_type (str) –
- Either ‘mrtrix’: bvecs are of shape (N, 4) or
’fsl’: bvecs are of shape (3, N)
- Returns:
new_bvecs – The final bvecs.
- Return type:
np.array
scilpy.gradients.gen_gradient_sampling module
Most of this code was modified from code by Emmanuel Caruyer <caruyer@gmail.com>.
See his original code on GitHub: https://github.com/ecaruyer/qspace/tree/master
The code was reorganized, but general process is kept the same.
- scilpy.gradients.gen_gradient_sampling.energy_comparison(bvecs1, bvecs2, nb_shells, nb_points_per_shell, alpha=1.0)[source]
Compute the electrostatic-repulsion energy of two sets of b-vectors with the same number of directions per shell.
- Parameters:
bvecs1 (array-like shape (N, 3,)) – First set of b-vectors.
bvecs2 (array-like shape (N, 3,)) – Second set of b-vectors.
nb_shells (int) – Number of shells
nb_points_per_shell (list of ints) – Number of points per shell.
alpha (float) – Controls the power of the repulsion. Default is 1.0
- Returns:
energy1 (float) – Electrostatic-repulsion energy of set bvecs1.
energy2 (float) – Electrostatic-repulsion energy of set bvecs2.
- scilpy.gradients.gen_gradient_sampling.generate_gradient_sampling(nb_samples_per_shell, verbose=1)[source]
Wrapper code to generate gradient sampling from Caruyer’s multiple_shell_energy.py
Code to generate multiple-shell gradient sampling, with optimal angular coverage. This implements the method described in Caruyer et al., MRM 69(6), pp. 1534-1540, 2013.
Generates the bvecs of a multiple shell gradient sampling using generalized Jones electrostatic repulsion.
- Parameters:
nb_samples_per_shell (list[int]) – Number of samples for each shell, starting from lowest.
verbose (int) – 0 = silent, 1 = summary upon completion, 2 = print iterations (To be sent to scipy).
- Returns:
bvecs (numpy.array of shape [n, 3]) – bvecs normalized to 1.
shell_idx (numpy.array) – Shell index for each bvec.
scilpy.gradients.optimize_gradient_sampling module
- scilpy.gradients.optimize_gradient_sampling.add_b0s_to_bvecs(bvecs, shell_idx, start_b0=True, b0_every=None, finish_b0=False)[source]
Add interleaved b0s to gradient sampling.
- Parameters:
bvecs (numpy.array,) – bvecs normalized to 1.
shell_idx (numpy.array) – Shell index for bvecs.
start_b0 (bool) – Option to add a b0 at the beginning.
b0_every (integer or None) – Final gradient sampling will have a b0 every b0_every samples. (start_b0 must be true)
finish_b0 (bool) – Option to add a b0 as last sample.
- Returns:
new_bvecs (numpy.array) – bvecs normalized to 1.
shell_idx (numpy.array) – Shell index for bvecs. Vectors with shells of value -1 are b0 vectors.
nb_new_b0s (int) – The number of b0s interleaved.
- scilpy.gradients.optimize_gradient_sampling.compute_bvalue_lin_b(bmin=0.0, bmax=3000.0, nb_of_b_inside=2, exclude_bmin=True)[source]
Compute bvals linearly distributed in b-value in the interval [bmin, bmax].
- Parameters:
bmin (float) – Minimum b-value, lower b-value bounds.
bmax (float) – Maximum b-value, upper b-value bounds.
nb_of_b_inside (int) – number of b-value excluding bmin and bmax.
exclude_bmin (boolean) – exclude bmin from the interval, useful if bmin = 0.0.
- Returns:
bvals – increasing bvals.
- Return type:
list
- scilpy.gradients.optimize_gradient_sampling.compute_bvalue_lin_q(bmin=0.0, bmax=3000.0, nb_of_b_inside=2, exclude_bmin=True)[source]
Compute bvals linearly distributed in q-value in the interval [bmin, bmax]. This leads to sqrt(b_values) linearly distributed.
- Parameters:
bmin (float) – Minimum b-value, lower b-value bounds.
bmax (float) – Maximum b-value, upper b-value bounds.
nb_of_b_inside (int) – number of b-value excluding bmin and bmax.
exclude_bmin (bool) – exclude bmin from the interval, useful if bmin = 0.0.
- Returns:
bvals – increasing bvals.
- Return type:
list
- scilpy.gradients.optimize_gradient_sampling.compute_min_duty_cycle_bruteforce(bvecs, shell_idx, bvals, ker_size=10, nb_iter=100000, rand_seed=0)[source]
Optimize the ordering of non-b0 samples to optimize gradient duty-cycle.
Philips scanner (and other) will find the peak power requirements with its duty cycle model (this is an approximation) and increase the TR accordingly to the hardware needs. This minimizes this effect by:
Randomly permuting the non-b0s samples
Finding the peak X, Y, and Z amplitude with a sliding-window
Computing the peak power needed as max(peak_x, peak_y, peak_z)
Keeping the permutation yielding the lowest peak power
- Parameters:
bvecs (numpy.array) – bvecs normalized to 1
shell_idx (numpy.array) – Shell index for bvecs.
bvals (list) – increasing bvals, b0 last.
ker_size (int) – kernel size for the sliding window.
nb_iter (int) – number of bruteforce iterations.
rand_seed (int) – seed for the random permutations.
- Returns:
new_bvecs (numpy.array) – bvecs normalized to 1.
shell_idx (numpy.array) – Shell index for bvecs.
- scilpy.gradients.optimize_gradient_sampling.compute_peak_power(q_scheme, ker_size=10)[source]
Function suggested by Guillaume Gilbert.
Optimize the diffusion gradient table by minimizing the maximum gradient load on any of the 3 axes over a preset temporal window (i.e. successive b-vectors).
In short, we want to avoid using the same gradient axis (x, y, or z) intensely for many successive b-vectors.
- Parameters:
q_scheme (nd.array) – Scheme of acquisition.
ker_size (int) – Kernel size (default=10).
- Return type:
Max peak power from q_scheme.
- scilpy.gradients.optimize_gradient_sampling.correct_b0s_philips(bvecs, shell_idx)[source]
Replace the [0.0, 0.0, 0.0] value of b0s bvecs by existing bvecs in the gradient sampling, except possibly the first one.
This is useful because Recon 1.0 of Philips allocates memory proportional to (total nb. of diff. bvals) x (total nb. diff. bvecs) and we can’t leave multiple b0s with b-vector [0.0, 0.0, 0.0] and b-value 0 because (b-vector, b-value) pairs have to be unique.
- Parameters:
bvecs (numpy.array) – bvecs normalized to 1
shell_idx (numpy.array) – Shell index for bvecs. Vectors with shells of value -1 are b0 vectors.
- Returns:
new_bvecs (numpy.array) – bvecs normalized to 1. b0 vectors are now impossible to know as they are replaced by random values from another vector.
shell_idx (numpy.array) – Shell index for bvecs. b0 vectors still have shells of value -1.
- scilpy.gradients.optimize_gradient_sampling.swap_sampling_eddy(bvecs, shell_idx)[source]
Optimize the bvecs of fixed multi-shell gradient sampling for eddy currents correction (fsl EDDY).
Bruteforce approach to maximally spread the bvec, shell per shell.
- For each shell:
- For each vector:
find the closest neighbor,
flips it,
if global system energy is better, keep it flipped
repeat until convergence.
- Parameters:
bvecs (numpy.array) – bvecs normalized to 1.
shell_idx (numpy.array) – Shell index for bvecs.
- Returns:
new_bvecs (numpy.array) – bvecs normalized to 1.
shell_idx (numpy.array) – Shell index for bvecs.
scilpy.gradients.utils module
- scilpy.gradients.utils.get_new_gtab_order(ref_gradient_table, dwi, bvals, bvecs)[source]
Find the sorting order that could be applied to the bval and bvec files to obtain the same order as in the reference gradient table.
This is mostly useful to reorder bval and bvec files in the order they were acquired by the Philips scanner (before version 5.6).
- Parameters:
ref_gradient_table (nd.array) – Gradient table, of shape (N, 4). It will use as reference for the ordering of b-vectors. Ex: Could be the result of scil_gradients_generate_sampling.py
dwi (nibabel image) – dwi of shape (x, y, z, N). Only used to confirm the dwi’s shape.
bvals (array, (N,)) – bvals that need to be reordered.
bvecs (array, (N, 3)) – bvecs that need to be reordered.
- Returns:
new_index – New index to reorder bvals/bvec
- Return type:
nd.array
- scilpy.gradients.utils.random_uniform_on_sphere(nb_vectors)[source]
Creates a set of K pseudo-random unit vectors, following a uniform distribution on the sphere. Reference: Emmanuel Caruyer’s code (https://github.com/ecaruyer).
This is not intended to create a perfect result. It’s usually the initialization step of a repulsion strategy.
- Parameters:
nb_vectors (int) – Number of vectors
- Returns:
bvecs – Pseudo-random unit vectors
- Return type:
nd.array of shape (nb_vectors, 3)