API Reference

Gradients

These functions deal with gradient pulse sequences.

SpinDoctor.GeneralGradientType
GeneralGradient(gvec, TE)

General gradient sequence $\vec{g}(t) \in \mathbb{R}^3$. The direction and amplitude may vary in time.

source
SpinDoctor.ScalarGradientType
ScalarGradient(dir, profile, amplitude)

Gradient sequence with amplitude amplitude, normal direction dir and scalar time profile profile.

The direction is constant, while the amplitude is controlled by the time profile.

The direction is normalized upon construction.

source
SpinDoctor.PGSEType
f = PGSE(δ, Δ)

Pulsed Gradient Spin Echo sequence with pulse duration δ and time between pulses Δ-δ.

source
SpinDoctor.DoublePGSEType
f = DoublePGSE(δ, Δ)

Double Pulsed Gradient Spin Echo sequence with four pulses of duration δ, separated by pauses of duration Δ-δ, 0, and Δ-δ repsectively.

source
SpinDoctor.CosOGSEType
f = CosOGSE(δ, Δ, nperiod)

Oscillating Gradient Spin Echo sequence with two cos-pulses of duration δ separated by a pause of duration Δ-δ for nperiod periods per pulse.

source
SpinDoctor.SinOGSEType
f = SinOGSE(δ, Δ, nperiod)

Oscillating Gradient Spin Echo sequence with two sin-pulses of duration δ separated by a pause of duration Δ-δ for nperiod periods per pulse.

source
SpinDoctor.integralFunction
integral(f, t = echotime(f))

Integral of time profile f between 0 and t. Unless specified, the echotime is used as the upper integral limit.

For the PGSE, SinOGSE, CosOGSE and DoublePGSE sequences, analytical expressions are available. Otherwise a numerical integral is computed.

source
SpinDoctor.int_F²Function
int_F²(f)

Compute the time profile contribution to the $b$-value given the time profile f. The $b$-value is given by

b = γ^2 * g^2 * int_F²(f)

where γ is the gyromagnetic ratio of the water proton and g is the gradient amplitude.

source
SpinDoctor.intervalsFunction
intervals(f::AbstractTimeProfile)

Get characteristic intervals of the time profile f.

intervals(g::AbstractGradient)

Get characteristic intervals of the time profile f.

source
SpinDoctor.echotimeFunction
echotime(f::AbstractTimeProfile)

Get echo time $TE$ of the time profile f, which is the end of the last characteristic interval.

echotime(grad::AbstractGradient)

Get echo time $TE$ of gradient.

source

Geometry

SpinDoctor.ModelType
Model(; mesh, ρ, D, T₂, κ, γ)

Finite element discretized biological model with initial spin densities ρ, diffusion tensors D, T₂-relaxation times T₂, wall permeabilities κ, and gyromagnetic ratio γ. The vectors ρ, D, T₂, are of length ncompartment, while κ is of length nboundary.

source

Postprocessing

SpinDoctor.savefieldFunction
savefield(mesh, ξ, filename::String; fieldname = "Magnetization")

Save field ξ to a VTK file. It may then be visualized in Paraview.

source
SpinDoctor.fit_adcFunction
fit_adc(bvalues, signals)

Fit the apparent diffusion coefficient (ADC) using a polynomial logfit of the normalized signals signals against the b-values bvalues.

source
SpinDoctor.fit_tensorsFunction
fit_tensors(directions, adcs)

Fit effective diffusion tensors to directionalized ADCs. The six components of the symmetric diffusion tensors are fitted by least squares to the gradient directions and resulting ADCs.

source
SpinDoctor.plot_meshFunction
plot_mesh(femesh, compartments = 1:ncompartment)

Plot finite element mesh, with a subset of the compartments.

source

Problems

SpinDoctor allows for considering a wide range of diffusion MRI problems. These can be solved for using the solve function.

SpinDoctor.MatrixFormalismType
MatrixFormalism(; model, matrices, lap_eig)

Matrix formalism problem. Given a Laplace eigendecomposition lap_eig, this problem consists of computing the MF magnetization.

source
SpinDoctor.AnalyticalMatrixFormalismType
AnalyticalMatrixFormalism(; analytical_laplace, lap_mat, volumes)

Analytical Matrix formalism problem. Given a radial Laplace eigendecomposition analytical_laplace, this problem consists of computing the MF compartment signals.

source
SpinDoctor.solveFunction
solve(problem::Laplace)

Compute the Laplace eigenvalues, eigenfunctions and first order moments of products of pairs of eigenfunctions.

See Jing-Rebecca Li, Try Tran, Van-Dang Nguyen (2020), Syver Døving Agdestein, Try Nguyen Tran, Jing-Rebecca Li (2022).

source
solve(problem::MatrixFormalism, gradient; ninterval = 500)

Solve for magnetization using Matrix Formalism.

See Jing-Rebecca Li, Try Tran, Van-Dang Nguyen (2020), Syver Døving Agdestein, Try Nguyen Tran, Jing-Rebecca Li (2022).

source
solve(problem::AnalyticalLaplace)

Solve analytical Laplace eigenvalue problem for symmetrical geometries.

source
solve(problem::AnalyticalMatrixFormalism, gradient)

Compute the signal in a multilayered cylinder or sphere using an analytical matrix formalism solution.

This function is based on the following articles and corresponding code: [1] Denis Grebenkov (2007) [2] Denis S. Grebenkov (2010)

source
solve(
    problem::BTPDE,
    gradient,
    odesolver = QNDF();
    abstol = 1e-6,
    reltol = 1e-4,
    callbacks = [].
)

Solve the Bloch-Torrey partial differential equation using P1 finite elements in space and odesolver in time.

source
function solve(
    problem::BTPDE,
    gradient::ScalarGradient,
    odesolver::IntervalConstantSolver;
    callbacks = [],
)

Solve the Bloch-Torrey partial differential equation using P1 finite elements in space and a theta-rule in time. This time stepping scheme requires a degree of implicitness θ and a time step Δt:

  • θ = 0.5: Crank-Nicolson (second order)
  • θ = 1.0: Implicit Euler (first order)

The function only works for interval-wise constant ScalarGradients, and errors otherwise.

source
solve(problem::Karger, gradient, odesolver = MagnusGL6(); timestep)

Solve the finite pulse Karger model (FPK) using precomputed effective diffusion tensors difftensors.

source
solve(
    problem::HADC,
    gradient::ScalarGradient,
    odesolver = QNDF();
    abstol = 1e-6,
    reltol = 1e-4,
)

Compute the ADC using a homogenized ADC model (HADC). This is currently only implemented for scalar gradients.

source

Matrix formalism

Utils

SpinDoctor.unitcircleFunction
unitcircle(ndirection; half = false, normal = [0, 0, 1])

Create ndirection directions unformly distributed on the unit circle defined by normal.

If half is true, the directions lie on a half-circle instead of the whole circle.

source
SpinDoctor.unitsphereFunction
unitsphere(ndirection; half = false, normal = [0, 0, 1])

Create ndirection directions unformly distributed on the unit sphere.

If half is true, the points will be discributed on the hemisphere defined by normal.

source
SpinDoctor.compute_signalFunction
compute_signal(M, ξ)

Compute signal from magnetization ξ, using the mass matrix M for integration.

Given a mesh mesh and a vector of compartment mass matrices M_cmpts, a vector of compartment signals may be obtained by compute_signal.(M_cmpts, split_field(mesh, ξ)).

source

Callbacks

Callbacks are called after every time step when solving the BTPDE.

SpinDoctor.VTKWriterType
VTKWriter(; nupdate = 1, dir = "output", filename = "solution")

Write magnetization field to a VTK file after every nupdate timestep. The files are stored in a ParaView data collection file (PVD). The magnetization time series may be visualized in ParaView by opening the file "$dir/$filename.pvd". The compartments are labeled.

source
SpinDoctor.PlotterType
Plotter(; nupdate = 1)

Plot the evolution of the BTPDE during time stepping. This requires loading the GLMakie plotting backend (]add GLMakie; using GLMakie). The plot is updated every nupdate time step. The resulting figure contains a plot of the time profile, total signal attenuation, and magnetization field (complex magnitude and phase shift).

source

Recipes

SpinDoctor comes with pre-configures recipes for creating finite element meshes of plates, cylinders, spheres, and neurons.

SpinDoctor.create_geometryFunction
create_geometry(setup; savedir = nothing, recreate = true)

Create cells, surfaces and finite element mesh. If savedir is a path, geometry files are saved. If additionally recreate = false, previous geometry files will be reused instead of generating new ones.

source