CylinderBundleCompartment {midi}R Documentation

Cylinder bundle compartment class

Description

A class to model restricted diffusion in a bundle of cylinders.

Methods

Public methods


Method new()

Instantiates a new cylinder bundle compartment.

Usage
CylinderBundleCompartment$new(
  axis,
  radius,
  diffusivity,
  cylinder_density,
  axial_diffusivity = NULL,
  radial_diffusivity = NULL,
  n_cylinders = 1L,
  axis_concentration = Inf,
  radius_sd = 0,
  radial_model = c("soderman", "callaghan", "stanisz", "neuman", "vangelderen"),
  seed = 1234
)
Arguments
axis

A numeric vector of length 3 and unit norm specifying the mean axis of the cylinder population.

radius

A positive numeric value specifying the mean radius of the cylinder population in meters.

diffusivity

A positive numeric value specifying the diffusivity within the cylinders in m^2.s^{-1}.

cylinder_density

A numeric value specifying the density of the cylinders in the voxel. Must be between 0 and 1.

axial_diffusivity

A numeric value specifying the axial diffusivity in the space outside the cylinders in m^2.s^{-1}. If not provided, defaults to a tortuosity model reducing the intrinsic diffusivity depending on orientation dispersion. Defaults to NULL.

radial_diffusivity

A numeric value specifying the radial diffusivity in the space outside the cylinders in m^2.s^{-1}. If not provided, defaults to a tortuosity model reducing the axial diffusivity depending on radius heterogeneity. Defaults to NULL.

n_cylinders

An integer value specifying the number of cylinders in the bundle. Defaults to 1L.

axis_concentration

A numeric value specifying the concentration of cylinders along the mean axis. Defaults to Inf.

radius_sd

A numeric value specifying the standard deviation of the radius of the cylinder population in meters. Defaults to 0.

radial_model

A character string specifying the radial model to use. Choices are "soderman", "callaghan", "stanisz", "neuman", and "vangelderen". Defaults to "soderman".

seed

An integer value specifying the seed for the random number generator. Defaults to 1234.

Returns

An instance of the CylinderBundleCompartment class.


Method get_signal()

Computes the signal attenuation predicted by the model.

Usage
CylinderBundleCompartment$get_signal(
  small_delta,
  big_delta,
  G,
  direction,
  echo_time = NULL,
  n_max = 20L,
  m_max = 50L
)
Arguments
small_delta

A numeric value specifying the duration of the gradient pulse in seconds.

big_delta

A numeric value specifying the duration between the gradient pulses in seconds.

G

A numeric value specifying the strength of the gradient in T.m^{-1}.

direction

A length-3 numeric vector specifying the direction of the gradient.

echo_time

A numeric value specifying the echo time in seconds.

n_max

An integer value specifying the maximum order of the Bessel function. Defaults to 20L.

m_max

An integer value specifying the maximum number of extrema for the Bessel function. Defaults to 50L.

Returns

A numeric value storing the predicted signal attenuation.

Examples
cylinderBundleComp <- CylinderBundleCompartment$new(
  axis = c(0, 0, 1),
  radius = 1e-5,
  diffusivity = 2.0e-9,
  cylinder_density = 0.5,
  radial_model = "soderman"
)
cylinderBundleComp$get_signal(
  small_delta = 0.03,
  big_delta = 0.03,
  G = 0.040,
  direction = c(0, 0, 1)
)

Method clone()

The objects of this class are cloneable with this method.

Usage
CylinderBundleCompartment$clone(deep = FALSE)
Arguments
deep

Whether to make a deep clone.

Examples


## ------------------------------------------------
## Method `CylinderBundleCompartment$get_signal`
## ------------------------------------------------

cylinderBundleComp <- CylinderBundleCompartment$new(
  axis = c(0, 0, 1),
  radius = 1e-5,
  diffusivity = 2.0e-9,
  cylinder_density = 0.5,
  radial_model = "soderman"
)
cylinderBundleComp$get_signal(
  small_delta = 0.03,
  big_delta = 0.03,
  G = 0.040,
  direction = c(0, 0, 1)
)

[Package midi version 0.1.0 Index]