sims#
Module: sims._force_core#
FORCE Core Cython Module
Low-level signal generation for multi-compartment diffusion simulations.
|
Compute the angle between two vectors. |
|
Create CSF free water diffusion signal. |
|
Create gray matter isotropic diffusion signal. |
|
Create mixed WM/GM/CSF tissue signal. |
|
Create white matter signal with specified number of fibers. |
|
Generate diffusion signal for a single fiber population. |
|
Generate diffusion signal for three crossing fiber populations. |
|
Generate diffusion signal for two crossing fiber populations. |
|
Check if an angle satisfies minimum separation constraint. |
|
Update the diffusivity ranges used by the simulators. |
Module: sims._multi_tensor_omp#
OpenMP-accelerated Multi-Tensor Signal Generation
High-performance implementation of multi-tensor diffusion signal computation with optional OpenMP parallelization.
|
Compute multi-tensor diffusion signal. |
|
Compute multi-tensor signals for a batch of voxels. |
Module: sims.force#
FORCE: Forward modeling for complex microstructure estimation
This module provides signal simulation for matching-based diffusion MRI reconstruction using multi-compartment tissue models.
The FORCE method generates a library of simulated signals representing various tissue configurations (white matter with 1-3 fibers, gray matter, CSF) and uses this library to reconstruct microstructural parameters from real diffusion MRI data.
References#
|
Generate spherical functions for all directions and ODI values. |
|
Find the smallest non-zero b-value shell. |
|
Initialize a ProcessPoolExecutor worker with a unique RNG state. |
|
Generate FORCE simulations. |
|
Save FORCE simulations to compressed NPZ file. |
|
Load FORCE simulations from NPZ file. |
|
Validate diffusivity configuration dictionary. |
Get default diffusivity configuration. |
Module: sims.phantom#
|
Add noise of specified distribution to a 4D array. |
|
numerical derivatives 2 eigenvectors |
|
Create a phantom based on a 3-D orbit |
Module: sims.voxel#
|
Add noise of specified distribution to the signal from a single voxel. |
|
Simulate the signal for a Sticks & Ball model. |
|
Calculates the perpendicular diffusion signal E(q) in a cylinder of radius R using the Soderman model. |
|
Calculates the parallel Gaussian diffusion signal. |
|
Calculates the three-dimensional signal attenuation E(q) originating from within a cylinder of radius R using the Soderman approximation. |
|
Simulate diffusion-weighted signals with a single tensor. |
|
Simulate a Multi-Tensor signal. |
|
Simulate the diffusion-weight signal, diffusion and kurtosis tensors based on the DKI model |
|
Computes the diffusion kurtosis tensor element (with indexes i, j, k and l) based on the individual diffusion tensor components of a multicompartmental model. |
|
Simulated signal based on the diffusion and diffusion kurtosis tensors of a single voxel. |
|
Simulate a Single-Tensor ODF. |
|
Given the principle tensor axis, return the array of all eigenvectors column-wise (or, the rotation matrix that orientates the tensor). |
|
Simulate a Multi-Tensor ODF. |
|
Simulate a Single-Tensor RTOP. |
|
Simulate a Multi-Tensor RTOP. |
|
Simulate a Single-Tensor PDF. |
|
Simulate a Multi-Tensor PDF. |
|
Simulate a Single-Tensor MSD. |
|
Simulate a Multi-Tensor MSD. |
angle_between#
- dipy.sims._force_core.angle_between(v1, v2)#
Compute the angle between two vectors.
- Parameters:
- v1ndarray (3,)
First vector.
- v2ndarray (3,)
Second vector.
- Returns:
- anglefloat
Angle between the two vectors in radians.
create_csf_signal#
- dipy.sims._force_core.create_csf_signal(bvals, target_sphere)#
Create CSF free water diffusion signal.
- Parameters:
- bvalsndarray
B-values.
- target_spherendarray (N, 3)
Unit sphere vertices.
- Returns:
- tuple
- (signal, labels, num_fibers, placeholder, free_water_fraction,
placeholder, odf, d_par, d_perp)
create_gm_signal#
- dipy.sims._force_core.create_gm_signal(bvals, target_sphere)#
Create gray matter isotropic diffusion signal.
- Parameters:
- bvalsndarray
B-values.
- target_spherendarray (N, 3)
Unit sphere vertices.
- Returns:
- tuple
- (signal, labels, num_fibers, dispersion, placeholder,
neurite_density, odf, d_par, d_perp)
create_mixed_signal#
- dipy.sims._force_core.create_mixed_signal(target_sphere, evecs, bingham_sf, odi_list, bvals, bvecs, multi_tensor_func, wm_threshold, tortuosity)#
Create mixed WM/GM/CSF tissue signal.
This is the main simulation function that generates realistic voxel signals with mixed tissue contributions.
- Parameters:
- target_spherendarray (N, 3)
Unit sphere vertices.
- evecsndarray (N, 3, 3)
Eigenvector matrices for each sphere direction.
- bingham_sfdict
Pre-computed Bingham spherical functions.
- odi_listndarray
List of orientation dispersion index values.
- bvalsndarray
B-values.
- bvecsndarray (M, 3)
Gradient directions.
- multi_tensor_funccallable
Function to compute multi-tensor signal.
- wm_thresholdfloat
Minimum WM fraction to include fiber labels.
- tortuositybool
Whether to use tortuosity constraint.
- Returns:
- tuple
- (signal, labels, num_fibers, dispersion, wm_fraction,
gm_fraction, csf_fraction, neurite_density, odf, ufa_wm, ufa_voxel, fiber_fractions, wm_disp, wm_d_par, wm_d_perp, gm_d_par, csf_d_par, f_ins)
create_wm_signal#
- dipy.sims._force_core.create_wm_signal(num_fib, target_sphere, evecs, bingham_sf, odi_list, bvals, bvecs, multi_tensor_func, tortuosity)#
Create white matter signal with specified number of fibers.
- Parameters:
- num_fibint
Number of fiber populations (1, 2, or 3).
- target_spherendarray (N, 3)
Unit sphere vertices.
- evecsndarray (N, 3, 3)
Eigenvector matrices for each sphere direction.
- bingham_sfdict
Pre-computed Bingham spherical functions.
- odi_listndarray
List of orientation dispersion index values.
- bvalsndarray
B-values.
- bvecsndarray (M, 3)
Gradient directions.
- multi_tensor_funccallable
Function to compute multi-tensor signal.
- tortuositybool
Whether to use tortuosity constraint.
- Returns:
- tuple
Signal and associated parameters.
generate_single_fiber#
- dipy.sims._force_core.generate_single_fiber(target_sphere, evecs, bingham_sf, odi_list, bvals, bvecs, multi_tensor_func, tortuosity)#
Generate diffusion signal for a single fiber population.
- Parameters:
- target_spherendarray (N, 3)
Unit sphere vertices.
- evecsndarray (N, 3, 3)
Eigenvector matrices for each sphere direction.
- bingham_sfdict
Pre-computed Bingham spherical functions.
- odi_listndarray
List of orientation dispersion index values.
- bvalsndarray
B-values.
- bvecsndarray (M, 3)
Gradient directions.
- multi_tensor_funccallable
Function to compute multi-tensor signal.
- tortuositybool
Whether to use tortuosity constraint.
- Returns:
- tuple
- (signal, labels, num_fibers, dispersion, placeholder,
neurite_density, odf, d_par, d_perp, fractions, f_ins)
generate_three_fibers#
- dipy.sims._force_core.generate_three_fibers(target_sphere, evecs, bingham_sf, odi_list, bvals, bvecs, multi_tensor_func, tortuosity)#
Generate diffusion signal for three crossing fiber populations.
- Parameters:
- target_spherendarray (N, 3)
Unit sphere vertices.
- evecsndarray (N, 3, 3)
Eigenvector matrices for each sphere direction.
- bingham_sfdict
Pre-computed Bingham spherical functions.
- odi_listndarray
List of orientation dispersion index values.
- bvalsndarray
B-values.
- bvecsndarray (M, 3)
Gradient directions.
- multi_tensor_funccallable
Function to compute multi-tensor signal.
- tortuositybool
Whether to use tortuosity constraint.
- Returns:
- tuple
Signal and associated parameters.
generate_two_fibers#
- dipy.sims._force_core.generate_two_fibers(target_sphere, evecs, bingham_sf, odi_list, bvals, bvecs, multi_tensor_func, tortuosity)#
Generate diffusion signal for two crossing fiber populations.
- Parameters:
- target_spherendarray (N, 3)
Unit sphere vertices.
- evecsndarray (N, 3, 3)
Eigenvector matrices for each sphere direction.
- bingham_sfdict
Pre-computed Bingham spherical functions.
- odi_listndarray
List of orientation dispersion index values.
- bvalsndarray
B-values.
- bvecsndarray (M, 3)
Gradient directions.
- multi_tensor_funccallable
Function to compute multi-tensor signal.
- tortuositybool
Whether to use tortuosity constraint.
- Returns:
- tuple
Signal and associated parameters.
is_angle_valid#
- dipy.sims._force_core.is_angle_valid(angle, *, threshold=30)#
Check if an angle satisfies minimum separation constraint.
- Parameters:
- anglefloat
Angle to check in radians.
- thresholdfloat, optional
Minimum separation threshold in degrees.
- Returns:
- is_validbool
True if the angle satisfies the separation constraint.
set_diffusivity_ranges#
- dipy.sims._force_core.set_diffusivity_ranges(wm_d_par_range=(0.002, 0.003), wm_d_perp_range=(0.0003, 0.0015), gm_d_iso_range=(0.0007, 0.0012), csf_d=0.003)#
Update the diffusivity ranges used by the simulators.
- Parameters:
- wm_d_par_rangetuple or float
White matter parallel diffusivity range (min, max) in mm^2/s. If a single float, uses fixed value.
- wm_d_perp_rangetuple or float
White matter perpendicular diffusivity range in mm^2/s.
- gm_d_iso_rangetuple or float
Gray matter isotropic diffusivity range in mm^2/s.
- csf_dfloat
CSF isotropic diffusivity in mm^2/s.
multi_tensor#
- dipy.sims._multi_tensor_omp.multi_tensor(mevals, evecs, fractions, bvals, bvecs)#
Compute multi-tensor diffusion signal.
Accumulates signal contributions from multiple diffusion tensors weighted by their volume fractions.
- Parameters:
- mevalsndarray (N, 3)
Eigenvalues for each tensor.
- evecsndarray (N, 3, 3)
Eigenvectors for each tensor.
- fractionsndarray (N,)
Volume fractions (should sum to 100).
- bvalsndarray (M,)
B-values.
- bvecsndarray (M, 3)
Gradient directions.
- Returns:
- signalndarray (M,)
Combined multi-tensor signal.
multi_tensor_batch#
- dipy.sims._multi_tensor_omp.multi_tensor_batch(signals_out, mevals_batch, evecs_batch, fractions_batch, bvals, bvecs, n_threads=1)#
Compute multi-tensor signals for a batch of voxels.
This function is optimized for processing multiple voxels and can utilize OpenMP parallelization.
- Parameters:
- signals_outndarray (B, M)
Output array for computed signals.
- mevals_batchndarray (B, N, 3)
Eigenvalues for each tensor in each voxel.
- evecs_batchndarray (B, N, 3, 3)
Eigenvectors for each tensor in each voxel.
- fractions_batchndarray (B, N)
Volume fractions for each voxel.
- bvalsndarray (M,)
B-values.
- bvecsndarray (M, 3)
Gradient directions.
- n_threadsint, optional
Number of OpenMP threads.
- Returns:
- signals_outndarray (B, M)
Computed signals (same as input, modified in place).
dispersion_lut#
- dipy.sims.force.dispersion_lut(target_sphere, odi_list)[source]#
Generate spherical functions for all directions and ODI values. Currently uses Bingham distribution, but can be extended to other dispersion models in the future.
- Parameters:
- target_spherendarray (N, 3)
Unit sphere vertices.
- odi_listndarray
List of orientation dispersion index values.
- Returns:
- SFdict
Nested dictionary mapping (vertex_index, odi) to spherical function (SF) values.
smallest_shell_bval#
- dipy.sims.force.smallest_shell_bval(bvals, *, b0_threshold=50, shell_tolerance=50, n=1)[source]#
Find the smallest non-zero b-value shell.
- Parameters:
- bvalsndarray
B-values array.
- b0_thresholdfloat, optional
Maximum b-value for b0 volumes.
- shell_tolerancefloat, optional
Tolerance for grouping shells.
- nint, optional
Number of smallest shells to return.
- Returns:
- min_shellfloat
Smallest non-zero shell b-value.
- shell_maskndarray
Boolean mask for volumes in that shell.
init_worker#
- dipy.sims.force.init_worker(base_seed=None)[source]#
Initialize a ProcessPoolExecutor worker with a unique RNG state.
With
base_seed=None, seed = PID + high-resolution time so every worker has a different stream. Withinitargs=(base_seed,)for reproducibility, seed = base_seed + PID so workers differ but the run is reproducible for a fixed worker count.- Parameters:
- base_seedint, optional
Base seed for reproducible worker initialization.
generate_force_simulations#
- dipy.sims.force.generate_force_simulations(gtab, *, num_simulations=100000, output_dir=None, num_cpus=1, batch_size=1000, wm_threshold=0.5, tortuosity=False, odi_range=(0.01, 0.3), num_odi_values=10, diffusivity_config=None, dtype=<class 'numpy.float32'>, compute_dti=True, compute_dki=False, verbose=True)[source]#
Generate FORCE simulations.
Creates a library of simulated diffusion MRI signals with corresponding microstructural parameters for matching-based reconstruction.
- Parameters:
- gtabGradientTable
Gradient table with b-values and b-vectors.
- num_simulationsint, optional
Number of simulated voxels.
- output_dirstr, optional
Directory for output files. If None, uses a temporary directory.
- num_cpusint or None, optional
Number of CPU cores for parallel processing.
1runs in-process (no subprocess, safe everywhere).Noneor-1uses all available cores. Values< -1leave that many cores idle (e.g.-2= all minus one).0raisesValueError.- batch_sizeint, optional
Batch size for processing.
- wm_thresholdfloat, optional
Minimum WM fraction to include fiber labels.
- tortuositybool, optional
Use tortuosity constraint for perpendicular diffusivity.
- odi_rangetuple, optional
(min, max) orientation dispersion index range.
- num_odi_valuesint, optional
Number of ODI values to sample.
- diffusivity_configdict, optional
Custom diffusivity ranges.
- dtypedtype, optional
Data type for outputs.
- compute_dtibool, optional
Compute DTI metrics (FA, MD, RD).
- compute_dkibool, optional
Compute DKI metrics (AK, RK, MK, KFA).
- verbosebool, optional
Enable progress output.
- Returns:
- simulationsdict
Simulations containing signals and parameters.
save_force_simulations#
load_force_simulations#
validate_diffusivity_config#
- dipy.sims.force.validate_diffusivity_config(config)[source]#
Validate diffusivity configuration dictionary.
- Parameters:
- configdict
Diffusivity configuration with keys: wm_d_par_range, wm_d_perp_range, gm_d_iso_range, csf_d.
- Returns:
- validbool
True if configuration is valid.
- Raises:
- ValueError
If configuration is invalid.
get_default_diffusivity_config#
add_noise#
- dipy.sims.phantom.add_noise(vol, *, snr=1.0, S0=None, noise_type='rician', rng=None)[source]#
Add noise of specified distribution to a 4D array.
- Parameters:
- volarray, shape (X,Y,Z,W)
Diffusion measurements in W directions at each
(X, Y, Z)voxel position.- snrfloat, optional
The desired signal-to-noise ratio. (See notes below.)
- S0float, optional
Reference signal for specifying snr.
- noise_typestring, optional
The distribution of noise added. Can be either ‘gaussian’ for Gaussian distributed noise, ‘rician’ for Rice-distributed noise (default) or ‘rayleigh’ for a Rayleigh distribution.
- rngnumpy.random.Generator class, optional
Numpy’s random generator for setting seed values when needed.
- Returns:
- volarray, same shape as vol
Volume with added noise.
Notes
SNR is defined here, following [2], as
S0 / sigma, wheresigmais the standard deviation of the two Gaussian distributions forming the real and imaginary components of the Rician noise distribution (see [3]).References
Examples
>>> signal = np.arange(800).reshape(2, 2, 2, 100) >>> signal_w_noise = add_noise(signal, snr=10, noise_type='rician')
diff2eigenvectors#
orbital_phantom#
- dipy.sims.phantom.orbital_phantom(*, gtab=None, evals=array([0.0015, 0.0004, 0.0004]), func=None, t=None, datashape=(64, 64, 64, 65), origin=(32, 32, 32), scale=(25, 25, 25), angles=None, radii=None, S0=100.0, snr=None, rng=None)[source]#
Create a phantom based on a 3-D orbit
f(t) -> (x,y,z).- Parameters:
- gtabGradientTable
Gradient table of measurement directions.
- evalsarray, shape (3,)
Tensor eigenvalues.
- funcuser defined function f(t)->(x,y,z)
It could be desirable for
-1=<x,y,z <=1. If None creates a circular orbit.- tarray, shape (K,)
Represents time for the orbit. Default is
np.linspace(0, 2 * np.pi, 1000).- datashapearray, shape (X,Y,Z,W)
Size of the output simulated data
- origintuple, shape (3,)
Define the center for the volume
- scaletuple, shape (3,)
Scale the function before applying to the grid
- anglesarray, shape (L,)
Density angle points, always perpendicular to the first eigen vector Default np.linspace(0, 2 * np.pi, 32).
- radiiarray, shape (M,)
Thickness radii. Default
np.linspace(0.2, 2, 6). angles and radii define the total thickness options- S0double, optional
Maximum simulated signal. Default 100.
- snrfloat, optional
The signal to noise ratio set to apply Rician noise to the data. Default is to not add noise at all.
- rngnumpy.random.Generator class, optional
Numpy’s random generator for setting seed values when needed. Default is None.
- Returns:
- dataarray, shape (datashape)
See also
Examples
>>> def f(t): ... x = np.sin(t) ... y = np.cos(t) ... z = np.linspace(-1, 1, len(x)) ... return x, y, z
>>> data = orbital_phantom(func=f)
add_noise#
- dipy.sims.voxel.add_noise(signal, snr, S0, *, noise_type='rician', rng=None)[source]#
Add noise of specified distribution to the signal from a single voxel.
- Parameters:
- signal1-d ndarray
The signal in the voxel.
- snrfloat
The desired signal-to-noise ratio. (See notes below.) If snr is None, return the signal as-is.
- S0float
Reference signal for specifying snr.
- noise_typestring, optional
The distribution of noise added. Can be either ‘gaussian’ for Gaussian distributed noise, ‘rician’ for Rice-distributed noise (default) or ‘rayleigh’ for a Rayleigh distribution.
- rngnumpy.random.Generator, optional
Random number generator for the noise. If
None, uses NumPy’s default random generator.
- Returns:
- signalarray, same shape as the input
Signal with added noise.
Notes
SNR is defined here, following [2], as
S0 / sigma, wheresigmais the standard deviation of the two Gaussian distributions forming the real and imaginary components of the Rician noise distribution (see [3]).References
Examples
>>> signal = np.arange(800).reshape(2, 2, 2, 100) >>> signal_w_noise = add_noise(signal, 10., 100., noise_type='rician')
sticks_and_ball#
- dipy.sims.voxel.sticks_and_ball(gtab, *, d=0.0015, S0=1.0, angles=((0, 0), (90, 0)), fractions=(35, 35), snr=20)[source]#
Simulate the signal for a Sticks & Ball model.
See [4] for a definition of the Sticks & Ball model.
- Parameters:
- gtabGradientTable
Signal measurement directions.
- dfloat, optional
Diffusivity value.
- S0float, optional
Unweighted signal value.
- anglesarray (K, 2) or (K, 3), optional
List of K polar angles (in degrees) for the sticks or array of K sticks as unit vectors.
- fractionsarray-like, optional
Percentage of each stick. Remainder to 100 specifies isotropic component.
- snrfloat, optional
Signal to noise ratio, assuming Rician noise. If set to None, no noise is added.
- Returns:
- S(N,) ndarray
Simulated signal.
- sticks(M,3)
Sticks in cartesian coordinates.
References
[4] Timothy E. J. Behrens, Heidi Johansen-Berg, Saad Jbabdi, Matthew F.S. Rushworth, and Mark W. Woolrich. Probabilistic diffusion tractography with multiple fibre orientations: What can we gain? NeuroImage, 34(1):144–155, 2007. URL: https://doi.org/10.1016/j.neuroimage.2006.09.018, doi:10.1016/j.neuroimage.2006.09.018.
callaghan_perpendicular#
- dipy.sims.voxel.callaghan_perpendicular(q, radius)[source]#
Calculates the perpendicular diffusion signal E(q) in a cylinder of radius R using the Soderman model.
Assumes that the pulse length is infinitely short and the diffusion time is infinitely long.
See [5] for details about the Soderman model.
- Parameters:
- qarray, shape (N,)
q-space value in 1/mm
- radiusfloat
cylinder radius in mm
- Returns:
- Earray, shape (N,)
signal attenuation
References
[5] (1,2) Olle Söderman and Bengt Jöonsson. Restricted Diffusion in Cylindrical Geometry. Journal of Magnetic Resonance, Series A, 117(1):94–97, 1995. URL: https://doi.org/10.1006/jmra.1995.0014, doi:10.1006/jmra.1995.0014.
gaussian_parallel#
cylinders_and_ball_soderman#
- dipy.sims.voxel.cylinders_and_ball_soderman(gtab, tau, *, radii=(0.005, 0.005), D=0.0007, S0=1.0, angles=((0, 0), (90, 0)), fractions=(35, 35), snr=20)[source]#
Calculates the three-dimensional signal attenuation E(q) originating from within a cylinder of radius R using the Soderman approximation.
The diffusion signal is assumed to be separable perpendicular and parallel to the cylinder axis [6].
This function is basically an extension of the ball and stick model. Setting the radius to zero makes them equivalent.
See [5] for details about the Soderman model.
- Parameters:
- gtabGradientTable
Signal measurement directions.
- taufloat
diffusion time in s
- radiiarray-like, optional
cylinder radius in mm
- Dfloat, optional
diffusion constant
- S0float, optional
Unweighted signal value.
- anglesarray (K, 2) or (K, 3), optional
List of K polar angles (in degrees) for the sticks or array of K sticks as unit vectors.
- fractionsarray-like, optional
Percentage of each stick. Remainder to 100 specifies isotropic component.
- snrfloat, optional
Signal to noise ratio, assuming Rician noise. If set to None, no noise is added.
- Returns:
- Earray, shape (N,)
signal attenuation
References
[6] Yaniv Assaf, Raisa Z. Freidlin, Gustavo K. Rohde, and Peter J. Basser. New modeling and experimental framework to characterize hindered and restricted water diffusion in brain white matter. Magnetic Resonance in Medicine, 52(5):965–978, 2004. URL: https://doi.org/10.1002/mrm.20274, doi:10.1002/mrm.20274.
single_tensor#
- dipy.sims.voxel.single_tensor(gtab, S0=1, *, evals=None, evecs=None, snr=None, rng=None)[source]#
Simulate diffusion-weighted signals with a single tensor.
- Parameters:
- gtabGradientTable
Table with information of b-values and gradient directions g. Note that if gtab has a btens attribute, simulations will be performed according to the given b-tensor B information.
- S0double, optional
Strength of signal in the presence of no diffusion gradient (also called the
b=0value).- evals(3,) ndarray, optional
Eigenvalues of the diffusion tensor. By default, values typical for prolate white matter are used.
- evecs(3, 3) ndarray, optional
Eigenvectors of the tensor. You can also think of this as a rotation matrix that transforms the direction of the tensor. The eigenvectors need to be column wise.
- snrfloat, optional
Signal to noise ratio, assuming Rician noise. None implies no noise.
- rngnumpy.random.Generator, optional
Random number generator for the noise. If
None, uses NumPy’s default random generator.
- Returns:
- S(N,) ndarray
- Simulated signal:
S(b, g) = S_0 e^(-b g^T R D R.T g), if gtab.tens=NoneS(B) = S_0 e^(-B:D), if gtab.tens information is given
References
[7] Edward O. Stejskal and John E. Tanner. Spin Diffusion Measurements: Spin Echoes in the Presence of a Time-Dependent Field Gradient. The Journal of Chemical Physics, 42(1):288–292, January 1965. doi:http://dx.doi.org/10.1063/1.1695690.
[8] Maxime Descoteaux. High Angular Resolution Diffusion MRI: from Local Estimation to Segmentation and Tractography. PhD thesis, Université Nice Sophia Antipolis, Valbonne, France, 2007. URL: https://theses.hal.science/tel-00457458, doi:.
multi_tensor#
- dipy.sims.voxel.multi_tensor(gtab, mevals, *, S0=1.0, angles=((0, 0), (90, 0)), fractions=(50, 50), snr=20, rng=None)[source]#
Simulate a Multi-Tensor signal.
- Parameters:
- gtabGradientTable
Table with information of b-values and gradient directions. Note that if gtab has a btens attribute, simulations will be performed according to the given b-tensor information.
- mevalsarray (K, 3)
each tensor’s eigenvalues in each row
- S0float, optional
Unweighted signal value (b0 signal).
- anglesarray (K, 2) or (K, 3), optional
List of K tensor directions in polar angles (in degrees) or unit vectors
- fractionsarray-like, optional
Percentage of the contribution of each tensor. The sum of fractions should be equal to 100%.
- snrfloat, optional
Signal to noise ratio, assuming Rician noise. If set to None, no noise is added.
- rngnumpy.random.Generator, optional
Random number generator for the noise. If
None, uses NumPy’s default random generator.
- Returns:
- S(N,) ndarray
Simulated signal.
- sticks(M,3)
Sticks in cartesian coordinates.
Examples
>>> import numpy as np >>> from dipy.sims.voxel import multi_tensor >>> from dipy.data import get_fnames >>> from dipy.core.gradients import gradient_table >>> from dipy.io.gradients import read_bvals_bvecs >>> fimg, fbvals, fbvecs = get_fnames(name='small_101D') >>> bvals, bvecs = read_bvals_bvecs(fbvals, fbvecs) >>> gtab = gradient_table(bvals, bvecs=bvecs) >>> mevals=np.array(([0.0015, 0.0003, 0.0003],[0.0015, 0.0003, 0.0003])) >>> e0 = np.array([1, 0, 0.]) >>> e1 = np.array([0., 1, 0]) >>> S = multi_tensor(gtab, mevals)
multi_tensor_dki#
- dipy.sims.voxel.multi_tensor_dki(gtab, mevals, *, S0=1.0, angles=((90.0, 0.0), (90.0, 0.0)), fractions=(50, 50), snr=20)[source]#
Simulate the diffusion-weight signal, diffusion and kurtosis tensors based on the DKI model
See [9] for further details.
- Parameters:
- gtabGradientTable
Gradient table.
- mevalsarray (K, 3)
eigenvalues of the diffusion tensor for each individual compartment
- S0float, optional
Unweighted signal value (b0 signal).
- anglesarray (K,2) or (K,3), optional
List of K tensor directions of the diffusion tensor of each compartment in polar angles (in degrees) or unit vectors
- fractionsfloat (K,), optional
Percentage of the contribution of each tensor. The sum of fractions should be equal to 100%.
- snrfloat, optional
Signal to noise ratio, assuming Rician noise. If set to None, no noise is added.
- Returns:
- S(N,) ndarray
Simulated signal based on the DKI model.
- dt(6,)
elements of the diffusion tensor.
- kt(15,)
elements of the kurtosis tensor.
Notes
Simulations are based on multicompartmental models which assumes that tissue is well described by impermeable diffusion compartments characterized by their only diffusion tensor. Since simulations are based on the DKI model, coefficients larger than the fourth order of the signal’s taylor expansion approximation are neglected.
References
[9] (1,2,3) Rafael Neto Henriques, Marta Morgado Correia, Rita Gouveia Nunes, and Hugo Alexandre Ferreira. Exploring the 3D geometry of the diffusion kurtosis tensor—Impact on the development of robust tractography procedures and novel biomarkers. NeuroImage, 111:85–99, 2015. URL: https://doi.org/10.1016/j.neuroimage.2015.02.004, doi:10.1016/j.neuroimage.2015.02.004.
Examples
>>> import numpy as np >>> from dipy.sims.voxel import multi_tensor_dki >>> from dipy.data import get_fnames >>> from dipy.core.gradients import gradient_table >>> from dipy.io.gradients import read_bvals_bvecs >>> fimg, fbvals, fbvecs = get_fnames(name='small_64D') >>> bvals, bvecs = read_bvals_bvecs(fbvals, fbvecs) >>> bvals_2s = np.concatenate((bvals, bvals * 2), axis=0) >>> bvecs_2s = np.concatenate((bvecs, bvecs), axis=0) >>> gtab = gradient_table(bvals_2s, bvecs=bvecs_2s) >>> mevals = np.array([[0.00099, 0, 0],[0.00226, 0.00087, 0.00087]]) >>> S, dt, kt = multi_tensor_dki(gtab, mevals)
kurtosis_element#
- dipy.sims.voxel.kurtosis_element(D_comps, frac, ind_i, ind_j, ind_k, ind_l, *, DT=None, MD=None)[source]#
Computes the diffusion kurtosis tensor element (with indexes i, j, k and l) based on the individual diffusion tensor components of a multicompartmental model.
- Parameters:
- D_comps(K,3,3) ndarray
Diffusion tensors for all K individual compartment of the multicompartmental model.
- frac[float]
Percentage of the contribution of each tensor. The sum of fractions should be equal to 100%.
- ind_iint
Element’s index i (0 for x, 1 for y, 2 for z)
- ind_jint
Element’s index j (0 for x, 1 for y, 2 for z)
- ind_kint
Element’s index k (0 for x, 1 for y, 2 for z)
- ind_l: int
Elements index l (0 for x, 1 for y, 2 for z)
- DT(3,3) ndarray, optional
Voxel’s global diffusion tensor.
- MDfloat, optional
Voxel’s global mean diffusivity.
- Returns:
- wijklfloat
kurtosis tensor element of index i, j, k, l
Notes
wijkl is calculated using equation 8 given in [9].
References
dki_signal#
- dipy.sims.voxel.dki_signal(gtab, dt, kt, *, S0=150, snr=None)[source]#
Simulated signal based on the diffusion and diffusion kurtosis tensors of a single voxel. Simulations are performed assuming the DKI model.
See [9] for further details.
- Parameters:
- gtabGradientTable
Measurement directions.
- dt(6,) ndarray
Elements of the diffusion tensor.
- kt(15, ) ndarray
Elements of the diffusion kurtosis tensor.
- S0float, optional
Strength of signal in the presence of no diffusion gradient.
- snrfloat, optional
Signal to noise ratio, assuming Rician noise. None implies no noise.
- Returns:
- S(N,) ndarray
Simulated signal based on the DKI model:
\[S=S_{0}e^{-bD+\frac{1}{6}b^{2}D^{2}K}\]
References
single_tensor_odf#
- dipy.sims.voxel.single_tensor_odf(r, *, evals=None, evecs=None)[source]#
Simulate a Single-Tensor ODF.
See [10] for further details.
- Parameters:
- r(N,3) or (M,N,3) ndarray
Measurement positions in (x, y, z), either as a list or on a grid.
- evals(3,)
Eigenvalues of diffusion tensor. By default, use values typical for prolate white matter.
- evecs(3, 3) ndarray
Eigenvectors of the tensor, written column-wise. You can also think of these as the rotation matrix that determines the orientation of the diffusion tensor.
- Returns:
- ODF(N,) ndarray
The diffusion probability at
rafter timetau.
References
[10] Iman Aganj, Christophe Lenglet, Guillermo Sapiro, Essa Yacoub, Kamil Ugurbil, and Noam Harel. Reconstruction of the orientation distribution function in single- and multiple-shell q-ball imaging within constant solid angle. Magnetic Resonance in Medicine, 64(2):554–566, 2010. URL: https://onlinelibrary.wiley.com/doi/abs/10.1002/mrm.22365, doi:https://doi.org/10.1002/mrm.22365.
all_tensor_evecs#
- dipy.sims.voxel.all_tensor_evecs(e0)[source]#
Given the principle tensor axis, return the array of all eigenvectors column-wise (or, the rotation matrix that orientates the tensor).
- Parameters:
- e0(3,) ndarray
Principle tensor axis.
- Returns:
- evecs(3,3) ndarray
Tensor eigenvectors, arranged column-wise.
multi_tensor_odf#
- dipy.sims.voxel.multi_tensor_odf(odf_verts, mevals, angles, fractions)[source]#
Simulate a Multi-Tensor ODF.
- Parameters:
- odf_verts(N,3) ndarray
Vertices of the reconstruction sphere.
- mevalssequence of 1D arrays,
Eigen-values for each tensor.
- anglessequence of 2d tuples,
Sequence of principal directions for each tensor in polar angles or cartesian unit coordinates.
- fractionssequence of floats,
Percentages of the fractions for each tensor.
- Returns:
- ODF(N,) ndarray
Orientation distribution function.
Examples
Simulate a MultiTensor ODF with two peaks and calculate its exact ODF.
>>> import numpy as np >>> from dipy.sims.voxel import multi_tensor_odf, all_tensor_evecs >>> from dipy.data import default_sphere >>> vertices, faces = default_sphere.vertices, default_sphere.faces >>> mevals = np.array(([0.0015, 0.0003, 0.0003],[0.0015, 0.0003, 0.0003])) >>> angles = [(0, 0), (90, 0)] >>> odf = multi_tensor_odf(vertices, mevals, angles, [50, 50])
single_tensor_rtop#
- dipy.sims.voxel.single_tensor_rtop(*, evals=None, tau=0.025330295910584444)[source]#
Simulate a Single-Tensor RTOP.
See [11] for further details.
- Parameters:
- evals1D arrays, optional
Eigen-values for the tensor. By default, values typical for prolate white matter are used.
- taufloat, optional
diffusion time. By default the value that makes q=sqrt(b).
- Returns:
- rtopfloat,
Return to origin probability.
References
multi_tensor_rtop#
- dipy.sims.voxel.multi_tensor_rtop(mf, *, mevals=None, tau=0.025330295910584444)[source]#
Simulate a Multi-Tensor RTOP.
See [11] for further details.
- Parameters:
- mfsequence of floats, bounded [0,1]
Percentages of the fractions for each tensor.
- mevalssequence of 1D arrays, optional
Eigen-values for each tensor. By default, values typical for prolate white matter are used.
- taufloat, optional
diffusion time. By default the value that makes q=sqrt(b).
- Returns:
- rtopfloat,
Return to origin probability.
References
single_tensor_pdf#
- dipy.sims.voxel.single_tensor_pdf(r, *, evals=None, evecs=None, tau=0.025330295910584444)[source]#
Simulate a Single-Tensor PDF.
See [11] for further details.
- Parameters:
- r(N,3) or (M,N,3) ndarray
Measurement positions in (x, y, z), either as a list or on a grid.
- evals(3,), optional
Eigenvalues of diffusion tensor. By default, use values typical for prolate white matter.
- evecs(3, 3) ndarray, optional
Eigenvectors of the tensor. You can also think of these as the rotation matrix that determines the orientation of the diffusion tensor.
- taufloat, optional
diffusion time. By default the value that makes q=sqrt(b).
- Returns:
- pdf(N,) ndarray
The diffusion probability at
rafter timetau.
References
multi_tensor_pdf#
- dipy.sims.voxel.multi_tensor_pdf(pdf_points, mevals, angles, fractions, *, tau=0.025330295910584444)[source]#
Simulate a Multi-Tensor PDF.
See [11] for further details.
- Parameters:
- pdf_points(N, 3) ndarray
Points to evaluate the PDF.
- mevalssequence of 1D arrays,
Eigen-values for each tensor. By default, values typical for prolate white matter are used.
- anglessequence,
Sequence of principal directions for each tensor in polar angles or cartesian unit coordinates.
- fractionssequence of floats,
Percentages of the fractions for each tensor.
- taufloat, optional
diffusion time. By default the value that makes q=sqrt(b).
- Returns:
- pdf(N,) ndarray,
Probability density function of the water displacement.
References
single_tensor_msd#
- dipy.sims.voxel.single_tensor_msd(*, evals=None, tau=0.025330295910584444)[source]#
Simulate a Single-Tensor MSD.
See [11] for further details.
- Parameters:
- evals1D arrays, optional
Eigen-values for the tensor. By default, values typical for prolate white matter are used.
- taufloat, optional
diffusion time. By default the value that makes q=sqrt(b).
- Returns:
- msdfloat,
Mean square displacement.
References
multi_tensor_msd#
- dipy.sims.voxel.multi_tensor_msd(mf, *, mevals=None, tau=0.025330295910584444)[source]#
Simulate a Multi-Tensor MSD.
See [11] for further details.
- Parameters:
- mfsequence of floats, bounded [0,1]
Percentages of the fractions for each tensor.
- mevalssequence of 1D arrays, optional
Eigen-values for each tensor. By default, values typical for prolate white matter are used.
- taufloat, optional
diffusion time. By default the value that makes q=sqrt(b).
- Returns:
- msdfloat,
Mean square displacement.
References