Title: | Microstructure Information from Diffusion Imaging |
---|---|
Description: | An implementation of a taxonomy of models of restricted diffusion in biological tissues parametrized by the tissue geometry (axis, diameter, density, etc.). This is primarily used in the context of diffusion magnetic resonance (MR) imaging to model the MR signal attenuation in the presence of diffusion gradients. The goal is to provide tools to simulate the MR signal attenuation predicted by these models under different experimental conditions. The package feeds a companion 'shiny' app available at <https://midi-pastrami.apps.math.cnrs.fr> that serves as a graphical interface to the models and tools provided by the package. Models currently available are the ones in Neuman (1974) <doi:10.1063/1.1680931>, Van Gelderen et al. (1994) <doi:10.1006/jmrb.1994.1038>, Stanisz et al. (1997) <doi:10.1002/mrm.1910370115>, Soderman & Jonsson (1995) <doi:10.1006/jmra.1995.0014> and Callaghan (1995) <doi:10.1006/jmra.1995.1055>. |
Authors: | Aymeric Stamm [aut, cre] |
Maintainer: | Aymeric Stamm <[email protected]> |
License: | MIT + file LICENSE |
Version: | 0.1.0.9000 |
Built: | 2024-10-28 02:48:48 UTC |
Source: | https://github.com/lmjl-alea/midi |
Plots a cross section of a cylinder bundle from an object of class bundle
as generated by simulate_bundle()
using
ggplot2.
## S3 method for class 'bundle' autoplot(object, grid_size = 100L, ...)
## S3 method for class 'bundle' autoplot(object, grid_size = 100L, ...)
object |
An object of class |
grid_size |
An integer value specifying the number of points on which
the unit circle should be discretized to draw the spheres. Defaults to
|
... |
Additional arguments to be passed to the |
A ggplot2::ggplot()
object.
density <- 0.5 voxel_size <- 5 # micrometers withr::with_seed(1234, { out <- simulate_bundle(density, voxel_size) }) ggplot2::autoplot(out)
density <- 0.5 voxel_size <- 5 # micrometers withr::with_seed(1234, { out <- simulate_bundle(density, voxel_size) }) ggplot2::autoplot(out)
A function to calculate the b-value for a given set of experimental parameters.
bvalue(small_delta, big_delta, G)
bvalue(small_delta, big_delta, G)
small_delta |
A numeric value specifying the duration of the gradient pulse in ms. |
big_delta |
A numeric value specifying the duration between the gradient pulses in ms. |
G |
A numeric value specifying the strength of the gradient in
|
A numeric value storing the predicted b-value in
ms/m
.
bvalue(small_delta = 30, big_delta = 30, G = 0.040)
bvalue(small_delta = 30, big_delta = 30, G = 0.040)
A class to model restricted diffusion in a cylinder using the Callaghan's model.
midi::BaseCompartment
-> midi::CircularlyShapedCompartment
-> CallaghanCompartment
clone()
The objects of this class are cloneable with this method.
CallaghanCompartment$clone(deep = FALSE)
deep
Whether to make a deep clone.
Callaghan, P. T. (1995). Pulsed-gradient spin-echo NMR for planar, cylindrical, and spherical pores under conditions of wall relaxation. Journal of magnetic resonance, Series A, 113(1), 53-59.
A class to model restricted diffusion in a bundle of cylinders.
midi::BaseCompartment
-> CylinderBundleCompartment
new()
Instantiates a new cylinder bundle compartment.
CylinderBundleCompartment$new( cylinder_density, cylinder_compartments, axial_diffusivity = NULL, radial_diffusivity = NULL )
cylinder_density
A numeric value specifying the density of the cylinders in the voxel. Must be between 0 and 1.
cylinder_compartments
A list of instances of the
CylinderCompartment
class specifying the compartments within the
bundle.
axial_diffusivity
A numeric value specifying the axial diffusivity
in the space outside the cylinders in m.s
. If not
provided, defaults to a tortuosity model reducing the intrinsic
diffusivity depending on orientation dispersion. Defaults to
NULL
in
which case the extracellular axial diffusivity is set via a tortuosity
model based on the dispersion in orientation.
radial_diffusivity
A numeric value specifying the radial
diffusivity in the space outside the cylinders in
m.s
. Defaults to
NULL
in which case the
extracellular radial diffusivity is set via a tortuosity model based on
the intracellular density.
An instance of the CylinderBundleCompartment
class.
clone()
The objects of this class are cloneable with this method.
CylinderBundleCompartment$clone(deep = FALSE)
deep
Whether to make a deep clone.
withr::with_seed(1234, { cyls <- rcylinders( n = 100, axis_mean = c(0, 0, 1), radius_mean = 5, diffusivity_mean = 3, axis_concentration = 10, radius_sd = 1, diffusivity_sd = 0 ) }) comp <- CylinderBundleCompartment$new( cylinder_density = 0.5, cylinder_compartments = cyls ) comp$get_signal( small_delta = 30, big_delta = 30, G = 0.040, direction = c(0, 0, 1) ) comp$get_parameters()
withr::with_seed(1234, { cyls <- rcylinders( n = 100, axis_mean = c(0, 0, 1), radius_mean = 5, diffusivity_mean = 3, axis_concentration = 10, radius_sd = 1, diffusivity_sd = 0 ) }) comp <- CylinderBundleCompartment$new( cylinder_density = 0.5, cylinder_compartments = cyls ) comp$get_signal( small_delta = 30, big_delta = 30, G = 0.040, direction = c(0, 0, 1) ) comp$get_parameters()
A class to model restricted diffusion in a cylinder.
midi::BaseCompartment
-> midi::CircularlyShapedCompartment
-> CylinderCompartment
new()
Instantiates a new cylinder compartment.
CylinderCompartment$new( axis = c(0, 0, 1), restricted_compartment = VanGelderenCompartment$new() )
axis
A length-3 numeric vector specifying the axis of the cylinder.
restricted_compartment
An instance of the
CircularlyShapedCompartment
class specifying the restricted
compartment within the sphere. Defaults to a Van Gelderen compartment.
An instance of the CylinderCompartment
class.
clone()
The objects of this class are cloneable with this method.
CylinderCompartment$clone(deep = FALSE)
deep
Whether to make a deep clone.
cylComp <- CylinderCompartment$new() cylComp$get_signal(small_delta = 30, big_delta = 30, G = 0.040) cylComp$get_parameter_names() cylComp$get_parameters()
cylComp <- CylinderCompartment$new() cylComp$get_signal(small_delta = 30, big_delta = 30, G = 0.040) cylComp$get_parameter_names() cylComp$get_parameters()
A class to model free unconstrained diffusion.
midi::BaseCompartment
-> FreeCompartment
new()
Instantiates a new free compartment.
FreeCompartment$new(diffusivity = NULL)
diffusivity
A numeric value specifying the diffusivity within the
sphere in m
.ms
. Defaults to
NULL
, in
which case the default free diffusivity of 3
m
.ms
is used.
An instance of the FreeCompartment
class.
clone()
The objects of this class are cloneable with this method.
FreeCompartment$clone(deep = FALSE)
deep
Whether to make a deep clone.
This class defines the gamma distribution. It provides methods for fitting the distribution to data and generating random samples.
midi::BaseDistribution
-> GammaDistribution
new()
Creates a new gamma distribution.
GammaDistribution$new(shape = 1, scale = 1)
shape
A numeric value specifying the shape parameter of the gamma
distribution. Defaults to 1.0
.
scale
A numeric value specifying the scale parameter of the gamma
distribution. Defaults to 1.0
.
An instance of the gamma distribution as an object of class
GammaDistribution
.
get_shape()
Retrieves the shape parameter of the gamma distribution.
GammaDistribution$get_shape()
A numeric value storing the shape parameter of the gamma distribution.
get_scale()
Retrieves the scale parameter of the gamma distribution.
GammaDistribution$get_scale()
A numeric value storing the scale parameter of the gamma distribution.
clone()
The objects of this class are cloneable with this method.
GammaDistribution$clone(deep = FALSE)
deep
Whether to make a deep clone.
gd <- GammaDistribution$new( shape = 1, scale = 5 ) gd$get_shape() gd$get_scale() gd$get_mean() gd$get_variance() gd$random(10) gd$fit(gd$random(100)) gd$get_shape() gd$get_scale()
gd <- GammaDistribution$new( shape = 1, scale = 5 ) gd$get_shape() gd$get_scale() gd$get_mean() gd$get_variance() gd$random(10) gd$fit(gd$random(100)) gd$get_shape() gd$get_scale()
A class to model restricted diffusion in a cylinder using the Neuman's model.
midi::BaseCompartment
-> midi::CircularlyShapedCompartment
-> NeumanCompartment
clone()
The objects of this class are cloneable with this method.
NeumanCompartment$clone(deep = FALSE)
deep
Whether to make a deep clone.
Neuman, C. H. (1974). Spin echo of spins diffusing in a bounded medium. The Journal of Chemical Physics, 60(11), 4508-4511.
Plots a cross section of a cylinder bundle
## S3 method for class 'bundle' plot(x, grid_size = 100L, ...)
## S3 method for class 'bundle' plot(x, grid_size = 100L, ...)
x |
An object of class |
grid_size |
An integer value specifying the number of points on which
the unit circle should be discretized to draw the spheres. Defaults to
|
... |
Additional arguments to be passed to the |
Nothing.
density <- 0.5 voxel_size <- 5 # micrometers withr::with_seed(1234, { out <- simulate_bundle(density, voxel_size) }) plot(out)
density <- 0.5 voxel_size <- 5 # micrometers withr::with_seed(1234, { out <- simulate_bundle(density, voxel_size) }) plot(out)
Plots a 3D representation of a cylinder bundle from an object of class
bundle
as generated by simulate_bundle()
using
plotly.
plot3d(b, show_linear_mesh = FALSE)
plot3d(b, show_linear_mesh = FALSE)
b |
An object of class |
show_linear_mesh |
A logical value indicating whether the linear mesh of
each cylinder should be displayed. Defaults to |
An HTML widget of class plotly::plotly
storing the 3D
visualization of the cylinder bundle.
density <- 0.5 voxel_size <- 5 # micrometers withr::with_seed(1234, { out <- simulate_bundle(density, voxel_size) }) plot3d(out)
density <- 0.5 voxel_size <- 5 # micrometers withr::with_seed(1234, { out <- simulate_bundle(density, voxel_size) }) plot3d(out)
This function samples cylinder compartments with given axis, radius
and diffusivity distributions.
rcylinders( n, axis_mean, radius_mean, diffusivity_mean, axis_concentration = Inf, radius_sd = 0, diffusivity_sd = 0, restricted_model = c("Callaghan", "Neuman", "Soderman", "Stanisz", "Van Gelderen") )
rcylinders( n, axis_mean, radius_mean, diffusivity_mean, axis_concentration = Inf, radius_sd = 0, diffusivity_sd = 0, restricted_model = c("Callaghan", "Neuman", "Soderman", "Stanisz", "Van Gelderen") )
n |
An integer value specifying the number of compartments to sample. |
axis_mean |
A numeric value specifying the mean of the axis distribution. |
radius_mean |
A numeric value specifying the mean of the radius distribution. |
diffusivity_mean |
A numeric value specifying the mean of the diffusivity distribution. |
axis_concentration |
A numeric value specifying the concentration of the
axis distribution. Defaults to |
radius_sd |
A numeric value specifying the standard deviation of the
radius distribution. Defaults to |
diffusivity_sd |
A numeric value specifying the standard deviation of
the diffusivity distribution. Defaults to |
restricted_model |
A character vector specifying the restricted
diffusion model to use. Defaults to |
The axis distribution is given by a mean and a concentration parameter and
the Dimroth-Watson distribution is used to sample values. The radius and
diffusivity distributions are given by a mean and a standard deviation and
the Gamma distribution is used to sample values. If the concentration
parameter is set to Inf
the axis is fixed to the mean value. If the
standard deviation of the radius and diffusivity distributions are set to 0
the radius and diffusivity are fixed to the mean values. If all parameters
are fixed, only one compartment is sampled.
A list of cylinder compartments of class
CylinderCompartment
.
# Sample 10 cylinder compartments with fixed axis, radius and diffusivity # set.seed(42) cyl_distr <- rcylinders( n = 10L, axis_mean = c(0, 0, 1), radius_mean = 5, diffusivity_mean = 3 )
# Sample 10 cylinder compartments with fixed axis, radius and diffusivity # set.seed(42) cyl_distr <- rcylinders( n = 10L, axis_mean = c(0, 0, 1), radius_mean = 5, diffusivity_mean = 3 )
This is a helper function to run the MIDI Shiny web application in the default web browser.
run_app()
run_app()
Nothing, but launches the Shiny app in the default web browser.
run_app()
run_app()
Generates a cross section of a cylinder bundle with a given density and voxel size. The cross section is generated by simulating a random distribution of cylinders and computing the intersection of the cylinders with a plane. The radii of the cylinders are drawn from a Gamma distribution fitted from data retrieved by histology in the genu of the corpus callosum (Aboitiz et al., 1992). The number of cylinders is determined by the density and the voxel size.
simulate_bundle( density = 0.5, voxel_size = 10, rel_tol = 0.001, verbose = FALSE )
simulate_bundle( density = 0.5, voxel_size = 10, rel_tol = 0.001, verbose = FALSE )
density |
A numeric value between 0 and 1 specifying the density of the
cylinders in the voxel. Defaults to |
voxel_size |
A numeric value specifying the size of the voxel in micro-
meters. Defaults to |
rel_tol |
A numeric value specifying the relative tolerance to reach the
target volume defined as |
verbose |
A logical value specifying whether to print messages. Defaults
to |
A list with the following components:
sections
: A numeric matrix with 3 columns:
x
: The x-coordinates of the centers of the cylinders;
y
: The y-coordinates of the centers of the cylinders;
r
: The radii of the cylinders in micrometers.
voxel_size
: The size of the voxel in micrometers
Aboitiz, F., Scheibel, A. B., Fisher, R. S., & Zaidel, E. (1992). Fiber composition of the human corpus callosum. Brain research, 598(1-2), 143-153.
density <- 0.5 voxel_size <- 5 # micrometers withr::with_seed(1234, { out <- simulate_bundle(density, voxel_size) }) # Actual density in the simulated substrate sum(out$sections[, "r"]^2 * pi) / voxel_size^2
density <- 0.5 voxel_size <- 5 # micrometers withr::with_seed(1234, { out <- simulate_bundle(density, voxel_size) }) # Actual density in the simulated substrate sum(out$sections[, "r"]^2 * pi) / voxel_size^2
A class to model restricted diffusion in a cylinder using the Soderman's model.
midi::BaseCompartment
-> midi::CircularlyShapedCompartment
-> SodermanCompartment
clone()
The objects of this class are cloneable with this method.
SodermanCompartment$clone(deep = FALSE)
deep
Whether to make a deep clone.
Söderman, O., & Jönsson, B. (1995). Restricted diffusion in cylindrical geometry. Journal of Magnetic Resonance, Series A, 117(1), 94-97.
A class to model restricted diffusion in a sphere.
midi::BaseCompartment
-> midi::CircularlyShapedCompartment
-> SphereCompartment
new()
Instantiates a new sphere compartment.
SphereCompartment$new(restricted_compartment = VanGelderenCompartment$new())
restricted_compartment
An instance of the
CircularlyShapedCompartment
class specifying the restricted
compartment within the sphere. Defaults to a Van Gelderen compartment.
An instance of the SphereCompartment
class.
clone()
The objects of this class are cloneable with this method.
SphereCompartment$clone(deep = FALSE)
deep
Whether to make a deep clone.
sphComp <- SphereCompartment$new() sphComp$get_signal(small_delta = 30, big_delta = 30, G = 0.040) sphComp$get_parameter_names() sphComp$get_parameters()
sphComp <- SphereCompartment$new() sphComp$get_signal(small_delta = 30, big_delta = 30, G = 0.040) sphComp$get_parameter_names() sphComp$get_parameters()
A class to model restricted diffusion in a cylinder using the Stanisz's model.
midi::BaseCompartment
-> midi::CircularlyShapedCompartment
-> StaniszCompartment
clone()
The objects of this class are cloneable with this method.
StaniszCompartment$clone(deep = FALSE)
deep
Whether to make a deep clone.
Stanisz, G. J., Wright, G. A., Henkelman, R. M., & Szafer, A. (1997). An analytical model of restricted diffusion in bovine optic nerve. Magnetic Resonance in Medicine, 37(1), 103-111.
A class to model restricted diffusion in a cylinder using the Van Gelderen's model.
midi::BaseCompartment
-> midi::CircularlyShapedCompartment
-> VanGelderenCompartment
clone()
The objects of this class are cloneable with this method.
VanGelderenCompartment$clone(deep = FALSE)
deep
Whether to make a deep clone.
Vangelderen, P., DesPres, D., Vanzijl, P. C. M., & Moonen, C. T. W. (1994). Evaluation of restricted diffusion in cylinders. Phosphocreatine in rabbit leg muscle. Journal of Magnetic Resonance, Series B, 103(3), 255-260.
This class defines the Watson distribution. It provides methods for fitting the distribution to data and generating random samples.
midi::BaseDistribution
-> WatsonDistribution
new()
Creates a new Watson distribution.
WatsonDistribution$new(mu = c(0, 0, 1), kappa = 10)
mu
A numeric vector of length 3 specifying the mean direction of
the Watson distribution. Defaults to (0, 0, 1)
.
kappa
A numeric value specifying the concentration parameter of
the Watson distribution. Defaults to 10.0
.
An instance of the Watson distribution as an object of class
WatsonDistribution
.
get_axis()
Retrieves the mean axis of the Watson distribution.
WatsonDistribution$get_axis()
A numeric vector of length 3 storing the mean axis of the Watson distribution.
get_concentration()
Retrieves the concentration parameter of the Watson distribution.
WatsonDistribution$get_concentration()
A numeric value storing the concentration parameter of the Watson distribution.
clone()
The objects of this class are cloneable with this method.
WatsonDistribution$clone(deep = FALSE)
deep
Whether to make a deep clone.
wd <- WatsonDistribution$new( mu = c(0, 0, 1), kappa = 10 ) wd$get_axis() wd$get_concentration() wd$random(10) wd$fit(wd$random(100)) wd$get_axis() wd$get_concentration()
wd <- WatsonDistribution$new( mu = c(0, 0, 1), kappa = 10 ) wd$get_axis() wd$get_concentration() wd$random(10) wd$fit(wd$random(100)) wd$get_axis() wd$get_concentration()