BlochSimulators.EPGStatesType
EPGStates = Union{MMatrix{3}, SizedMatrix{3}}

In the EPG model, the configuration state matrix Ω will be updated inplace. On CPU, we use StaticArrays.MMatrix. The MMatrix will not escape _simulate! and therefore should not result in allocations. On GPU, we use shared memory (CUDA.CuStaticSharedArray). The shared memory is allocated for all threads within a block simultaneously. We then take a @view and wrap it in a SizedArray. Without the SizedArray, the type would be long and unreadable. By wrapping, we can then simply dispatch on a SizedArray instead.

BlochSimulators.AbstractTissuePropertiesType
AbstractTissueProperties{N,T} <: FieldVector{N,T}

Abstract type for custom structs that hold tissue properties used for a simulation within one voxel. For simulations, SimulationParameterss are used that can be assembled with the @parameters macro.

Possible fields:

  • T₁::T: T₁ relaxation parameters of a voxel
  • T₂::T: T₂ relaxation parameters of a voxel
  • B₁::T: Scaling factor for effective B₁ excitation field within a voxel
  • B₀::T: Off-resonance with respect to main magnetic field within a voxel
  • ρˣ::T: Real part of proton density within a voxel
  • ρʸ::T: Imaginary part of proton density within a voxel

Implementation details:

The structs are subtypes of FieldVector, which is a StaticVector with named fields (see the documentation of StaticArrays.jl). There are three reasons for letting the structs be subtypes of FieldVector:

  1. FieldVectors/StaticVectors have sizes that are known at compile time. This is beneficial for performance reasons
  2. The named fields improve readability of the code (e.g. p.B₁ vs p[3])
  3. Linear algebra operations can be performed on instances of the structs. This allows, for example, subtraction (without having to manually define methods) and that is useful for comparing parameter maps.
BlochSimulators.AbstractTrajectoryType
AbstractTrajectory{T}

The abstract type of which all gradient trajectories will be a subtype. The subtypes should contain fields that can describe the full trajectory during a sequence. The type T refers to the precision of the floating point values within the trajectory struct.

BlochSimulators.AdiabaticInversionType
AdiabaticInversion{T<:Real, V<:AbstractVector} <: IsochromatSimulator{T}

This struct is used to simulate an adiabatic inversion pulse. This struct itself could be used as field in other sequence structs.

Fields

  • γΔtA::V: # Time-dependent amplitude modulation.
  • Δω::V # Time-dependent frequency modulation
  • Δt::T: Time discretization step, assumed constant
BlochSimulators.BlochSimulatorType
BlochSimulator{T}

The abstract type of which all sequence simulators will be a subtype. The parameter T should be a number type (e.g. Float64, Float32) and the tissueparameters that are used as input to the simulator should have the same number type. By convention, a BlochSimulator will be used to simulate magnetization at echo times only without taking into account spatial encoding gradients (i.e. readout or phase encoding gradients). To simulate the magnetization at other readout times, including phase from spatial encoding gradients, an AbstractTrajectory will be needed as well.

To make a simulator for a particular pulse sequence:

  1. Make a struct that's a subtype of either IsochromatSimulator or EPGSimulator. The struct will hold parameters that are necessary for performing the simulations.

  2. Add a method to simulate_magnetization! that implements the pulse sequence. For both performance and GPU compatibility, make sure that simulate_magnetization! does not do any heap allocations. Examples for pSSFP and FISP sequences are found in src/sequences.

  3. Add methods to output_eltype and output_size that are used to allocate an output array within the simulate function.

  4. [Optional] Add a method to show for nicer printing of the sequence in the REPL

  5. [Optional] Add a method to getindex to easily reduce the length of the sequence

  6. [Optional] Add a constructor for the struct that takes in data from Matlab or something else and assembles the struct

IMPORTANT

The simulate_magnetization! functions (which dispatch on the provided sequence) are assumed to be type-stable and non-allocating Should be possible to achieve when using functions from operators/epg.jlandoperators/isochromat.jl` and a properly parametrized sequence struct.

BlochSimulators.CartesianTrajectoryType
CartesianTrajectory{T,I,U,V} <: SpokesTrajectory{T}

Struct that is used to implement a typical Cartesian gradient trajectory. The trajectory is described in a compact fashion by only storing the starting position in k-space (k_start_readout) for each readout as well as the step in k-space per readout point Δk_adc.

Note that CartesianTrajectory and RadialTrajectory are essentially the same in when using when using this compact description. A SpokesTrajectory struct is therefore defined as a supertype of both and methods are defined for SpokesTrajectory instead to avoid code repetition.

The type parameters are intentionally left vague. The J, for example, may be an integer for sequences where each readout has the same number of samples, but for sequences with different numbers of samples per readout it may be a vector of integers.

Fields

  • nreadouts::I: The total number of readouts for this trajectory
  • nsamplesperreadout::I: The total number of samples per readout
  • Δt::T: Time between sample points
  • k_start_readout::U: Starting position in k-space for each readout
  • Δk_adc::U: k-space step Δkₓ per sample point (same for all readouts)
  • py::V: Phase encoding index for each readout
  • readout_oversampling::I: Readout oversampling factor
BlochSimulators.CoordinatesType
Coordinates{T<:Real}

Basic type that holds the spatial coordinates of a voxel. Note that when performing signal simulations, StructArray{<:Coordinates}s are used to store the coordinates of all voxels. Such arrays are created using the make_coordinates.

Fields

  • x::T: Position of voxel along the x direction
  • y::T: Position of voxel along the y direction
  • z::T: Position of voxel along the z direction (not used for 2D simulations)
BlochSimulators.EPGSimulatorType
EPGSimulator{T,Ns} <: BlochSimulator{T}

Abstract type of which all sequence simulators that are based on the EPG model will be a subtype. The parameter T should be a number type (e.g. Float64, Float32) and the tissueparameters that are used as input to the simulator should have the same number type. The parameter Ns corresponds to the maximum order of configuration states that are tracked in the simulations.

BlochSimulators.FISP2DType
FISP2D{T, Ns, U<:AbstractVector, V<:AbstractMatrix} <: EPGSimulator{T,Ns}

This struct is used to simulate gradient-spoiled sequence with varying flip angle scheme and adiabatic inversion prepulse using the EPG model. The TR and TE are fixed throughout the sequence. Instantenous RF excitations are assumed. Slice profile correction is done using the partitioned EPG model, where for each flip angle a vector of RF scaling factors are determined prior to the simulation (using, for example, the small tip-angle approximation or Shinnar LeRoux forward model).

Within each TR, a single time steps is used to simulate the RF excitation. Then, in one time step we go from the end of the RF excitation to the echo time (applying T₁ and T₂ decay, T₁ regrowth and B₀ rotation), and again in one time step from the echo time to the start of the next RF excitation.

Fields

  • RF_train::U Vector with flip angle for each TR with abs.(RFtrain) the RF flip angles in degrees and angle.(RFtrain) should be the RF phases in degrees.
  • sliceprofiles::V # Matrix with RF scaling factors (a.u.) to simulate slice profile effects. Each column represents the (flip angle dependent) scaling factors for one position along the slice direction.
  • TR::T: Repetition time in seconds, assumed constant during the sequence
  • TE::T: Echo time in seconds, assumed constant during the sequence
  • max_state::Val{Ns}: Maximum number of states to keep track of in EPG simulation
  • TI::T: Inversion delay after the inversion prepulse in seconds
BlochSimulators.Generic2DType
Generic2D{T,V,M,S} where {T<:AbstractFloat, V<:AbstractVector, M<:AbstractMatrix, S} <: IsochromatSimulator{T}

Simulate a generic 2D sequence defined by arrays containing RF and gradient waveforms. Contains a loop over z locations to take into account slice profile effects. The Δt vector stores the time intervals for the waveforms.

Fields

  • RF::V{Complex{T}}: Vector with (complex) RF values during each time interval
  • GR::M{T}: Matrix with GRx, GRy and GRz values during each time interval
  • sample::S: Vector with Bool's to indicate the sample points
  • Δt::V{T}: Vector with time intervals
  • z::V{T}: Vector with different positions along the slice direction
BlochSimulators.Generic3DType
Generic3D{T,V<:AbstractVector{Complex{T}},W<:AbstractVector{T},M<:AbstractMatrix{T},S} <: IsochromatSimulator{T}

Simulate a generic sequence defined by arrays containing RF and gradient waveforms. Unlike the Generic2D sequence, it is assumed that the excitation is homogenous over the voxel and therefore no summation over a slice direction is applied. The Δt vector stores the time intervals for the waveforms.

Fields

  • RF::V{Complex{T}}: Vector with (complex) RF values during each time interval
  • GR::M{T}: Matrix with GRx, GRy and GRz values during each time interval
  • sample::S: Vector with Bool's to indicate the sample points
  • Δt::V{T}: Vector with time intervals
BlochSimulators.IsochromatType
struct Isochromat{T<:Real} <: FieldVector{3,T}
    x::T
    y::T
    z::T
end

Holds the x,y,z components of a spin isochromat in a FieldVector, which is a StaticVector (from the package StaticArrays) with custom fieldnames.

BlochSimulators.IsochromatSimulatorType
IsochromatSimulator{T} <: BlochSimulator{T}

Abstract type of which all sequence simulators that are based on the isochromat model will be a subtype. The parameter T should be a number type (e.g. Float64, Float32) and the tissueparameters that are used as input to the simulator should have the same number type.

BlochSimulators.RadialTrajectoryType
RadialTrajectory{T,I,U,V} <: SpokesTrajectory{T}

Struct that is used to implement a typical radial gradient trajectory. The trajectory can is described in a compact fashion by only storing the starting position in k-space (k_start_readout) for each readout as well as the step in k-space per readout point Δk_adc.

Note that CartesianTrajectory and RadialTrajectory are essentially the same in when using when using this compact description. A SpokesTrajectory struct is therefore defined as a supertype of both and methods are defined for SpokesTrajectory instead to avoid code repetition.

The type parameters are intentionally left vague. The J, for example, may be an integer for sequences where each readout has the same number of samples, but for sequences with different numbers of samples per readout it may be a vector of integers.

Fields

  • nreadouts::I: The total number of readouts for this trajectory
  • nsamplesperreadout::I: The total number of samples per readout
  • Δt::T: Time between sample points
  • k_start_readout::U: Starting position in k-space for each readout
  • Δk_adc::U: k-space step Δk between each readout
  • φ::V: Radial angle for each readout
  • readout_oversampling::I: Readout oversampling factor
BlochSimulators.SpokesTrajectoryType
SpokesTrajectory{T} <: AbstractTrajectory{T}

Typical Cartesian and radial trajectories have a lot in common: a readout can be described by a starting point in k-space and a Δk per sample point. To avoid code repetition, both type of trajectories are made a subtype of SpokesTrajectory such that some methods that would be the same for both trajectories otherwise are written for SpokesTrajectory instead.

BlochSimulators.pSSFP2DType
pSSFP2D{T<:AbstractFloat,N,M,U<:AbstractVector{Complex{T}},V<:Number} <: IsochromatSimulator{T}

This struct is used to simulate a inversion-recovery, gradient-balanced transient-state sequence with varying flip angle scheme based on the isochromat model. The TR and TE are fixed throughout the sequence. The TR and TE are fixed throughout the sequence. Slice profile correction is done by discretizing the RF excitation waveform in time and using multiple Isochromats with different positions along the slice direction (z) per voxel. The sequence also uses an 'α/2' prepulse after the inversion.

Within each TR, multiple time steps are used to simulate the RF excitation. Then, in one time step we go from the end of the RF excitation to the echo time (applying slice refocussing gradient, T₂ decay and B₀ rotation), and again in one time step from the echo time to the start of the next RF excitation.

Fields

  • RF_train::U Vector with flip angle for each TR with abs.(RFtrain) the RF flip angles in degrees and angle.(RFtrain) should be the RF phases in degrees.
  • TR::T: Repetition time in seconds, assumed constant during the sequence
  • γΔtRF::SVector{N}{V}: Time-discretized RF waveform, normalized to flip angle of 1 degree
  • Δt::NamedTuple{(:ex, :inv, :pr),NTuple{3,T}}: Time interval for each sample of excitation pulse (ex), inversion delay (inv) and time between RF and TE (pr)
  • γΔtGRz::NamedTuple{(:ex, :inv, :pr),NTuple{3,T}}: Slice select gradients for ex, inv and pr
  • z::SVector{M}{T} # Vector with different positions along the slice direction.
BlochSimulators.pSSFP3DType
pSSFP3D{T<:AbstractFloat,N,U<:AbstractVector{Complex{T}},V<:Number} <: IsochromatSimulator{T}

This struct is used to simulate an inversion-recovery, gradient-balanced, transient-state sequence with varying flip angle scheme based on the isochromat model. The TR and TE are fixed throughout the sequence. The RF excitation waveform can be discretized in time but no slice profile mechanism is provided. The sequence also uses an 'α/2' prepulse after the inversion.

Within each TR, multiple time steps are used to simulate the RF excitation. Then, in one time step we go from the end of the RF excitation to the echo time (applying slice refocussing gradient, T₂ decay and B₀ rotation), and again in one time step from the echo time to the start of the next RF excitation.

Fields

  • RF_train::U Vector with flip angle for each TR with abs.(RFtrain) the RF flip angles in degrees and angle.(RFtrain) should be the RF phases in degrees.
  • TR::T: Repetition time in seconds, assumed constant during the sequence
  • γΔtRF::SVector{N}{V}: Time-discretized RF waveform, normalized to flip angle of 1 degree
  • Δt::NamedTuple{(:ex, :inv, :pr),NTuple{3,T}}: Time interval for each sample of excitation pulse (ex), inversion delay (inv) and time between RF and TE (pr)
BlochSimulators.F̄₋Method
F̄₋(Ω)

View into the second row of the configuration state matrix Ω, corresponding to the F̄₋ states.

BlochSimulators.F₊Method
F₊(Ω)

View into the first row of the configuration state matrix Ω, corresponding to the F₊ states.

BlochSimulators.ZMethod
Z(Ω)

View into the third row of the configuration state matrix Ω, corresponding to the Z states.

BlochSimulators._all_arrays_are_cuarraysMethod
_all_arrays_are_cuarrays(x)

Returns true if all AbstractArray fields in x are CuArrays and false otherwise. Will also return if x does not have any AbstractArray fields.

BlochSimulators._allocate_array_on_resourceMethod
_allocate_array_on_resource(resource, _eltype, _size)

Allocate an array on the specified resource with the given element type _eltype and size _size. If resource is CPU1() or CPUThreads(), the array is allocated on the CPU. If resource is CUDALibs(), the array is allocated on the GPU. For CPUProcesses(), the array is distributed in the "voxel"-dimension over multiple CPU workers.

This function is called by _allocate_magnetization_array and is not intended considered part of te public API.

BlochSimulators._allocate_magnetization_arrayMethod
_allocate_magnetization_array(resource, sequence, parameters)

Allocate an array to store the output of the Bloch simulations (per voxel, echo times only) to be performed with the sequence. For each BlochSimulator, methods should have been added to output_eltype and output_size for this function to work properly.

This function is called by simulate_magnetization and is not intended considered part of te public API.

Returns

  • magnetization_array: An array allocated on the specified resource, formatted to store the simulation results for each voxel across the specified echo times.
BlochSimulators._allocate_signal_arrayMethod
_allocate_signal_array(resource, trajectory::AbstractTrajectory, coil_sensitivities)

Allocate an array to store the output of the signal simulation (all readout points, integrated over all voxels).

BlochSimulators._get_readout_and_sample_idxMethod
_get_readout_and_sample_idx(trajectory, t)

Given time index t, compute the associated readout and sample indices (r,s). For trajectories where each readout has the same length, this is equivalent to r,s = fld1(t,ns), mod1(t,ns) with ns being the (constant) number of samples per readout.

BlochSimulators._signal_per_coil!Method
_signal_per_coil!(signal, resource, magnetization, parameters, trajectory, coordinates, coil_sensitivities)

Compute the signal for a given coil by calculating a volume integral of the transverse magnetization in each voxel for each time point separately (using the signal_at_time_point! function). Each time point is computed in parallel for multi-threaded CPU computation and on the GPU for CUDA computation.

BlochSimulators.add_gradient_delay!Method
add_gradient_delay!(tr::RadialTrajectory, S)

Apply gradient delay to radial trajectory in in-place fashion. The delay is described by the 2x2 matrix S and is assumed to influence the start of the readout only, not the readout direction.

BlochSimulators.decay!Method
decay!(Ω::EPGStates, E₁, E₂)

T₂ decay for F-components, T₁ decay for Z-component of each state.

BlochSimulators.decayMethod
decay(m::Isochromat{T}, E₁, E₂) where T

Apply T₂ decay to transverse component and T₁ decay to longitudinal component of Isochromat.

BlochSimulators.dephasing!Method
dephasing!(Ω::EPGStates)

Shift states around due to dephasing gradient: The F₊ go up one, the F̄₋ go down one and Z do not change

BlochSimulators.excite!Method
excite!(Ω::EPGStates, RF::Complex, p::AbstractTissueProperties)

Mixing of states due to RF pulse. Magnitude of RF is the flip angle in degrees. Phase of RF is the phase of the pulse. If RF is real, the computations simplify a little bit.

BlochSimulators.excite!Method
excite!(Ω::EPGStates, RF::T, p::AbstractTissueProperties) where T<:Union{Real, Quantity{<:Real}}

If RF is real, the calculations simplify (and probably Ω is real too, reducing memory (access) requirements).

BlochSimulators.f32Method
f32(x)

Change precision of x to Float32. It uses Functors.fmap to recursively traverse the fields of the struct x. For custom structs (e.g. <:BlochSimulator or <:AbstractTrajectory), it is required that typeof(x) be made a Functors.@functor (e.g. @functor FISP).

It may be necessary to add new adapt rules (by adding new methods to adapt_storage) if new structs with complicated nested fields are introduced.

BlochSimulators.f64Method
f64(x)

Change precision of x to Float64. It uses Functors.fmap to recursively traverse the fields of the struct x. For custom structs (e.g. <:BlochSimulator or <:AbstractTrajectory), it is required that typeof(x) be made a Functors.@functor (e.g. @functor FISP).

It may be necessary to add new adapt rules (by adding new methods to adapt_storage) if new structs with complicated nested fields are introduced.

BlochSimulators.gpuMethod
gpu(x)

Move x to CUDA device. It uses Functors.fmap to recursively traverse the fields of the struct x, converting <:AbstractArrays to CuArrays, and ignoring isbitsarrays. For custom structs (e.g. <:BlochSimulator or <:AbstractTrajectory), it is required that typeof(x) be made a Functors.@functor (e.g. @functor FISP).

BlochSimulators.initial_conditions!Method
initial_conditions!(Ω::EPGStates)

Set all components of all states to 0, except the Z-component of the 0th state which is set to 1.

BlochSimulators.initialize_statesMethod
initialize_states(::AbstractResource, sequence::EPGSimulator{T,Ns}) where {T,Ns}

Initialize an MMatrix of EPG states on CPU to be used throughout the simulation.

BlochSimulators.initialize_statesMethod
initialize_states(::CUDALibs, sequence::EPGSimulator{T,Ns}) where {T,Ns}

Initialize an array of EPG states on a CUDA GPU to be used throughout the simulation.

BlochSimulators.initialize_statesMethod
initialize_states(::AbstractResource, ::IsochromatSimulator{T}) where T

Initialize a spin isochromat to be used throughout a simulation of the sequence.

This may seem redundant but to is necessary to share the same programming interface with EPGSimulators.

BlochSimulators.invert!Method
invert!(Ω::EPGStates, p::AbstractTissueProperties)

Invert Z-component of states of all orders. Assumes fully spoiled transverse magnetization.

BlochSimulators.invert!Method
invert!(Ω::EPGStates)

Invert with B₁ insenstive (i.e. adiabatic) inversion pulse

BlochSimulators.invertMethod
invert(m::Isochromat{T}, p::AbstractTissueProperties) where T

Invert Isochromat with B₁ insenstive (i.e. adiabatic) inversion pulse

BlochSimulators.invertMethod
invert(m::Isochromat{T}, p::AbstractTissueProperties) where T

Invert z-component of Isochromat (assuming spoiled transverse magnetization so xy-component zero).

BlochSimulators.kspace_coordinatesMethod
kspace_coordinates(tr::CartesianTrajectory)

Return matrix (nrsamplesperreadout, nrreadouts) with kspace coordinates for the trajectory. Needed for nuFFT reconstructions.

BlochSimulators.kspace_coordinatesMethod
kspace_coordinates(tr::RadialTrajectory)

Return matrix (nrsamplesperreadout, nrreadouts) with kspace coordinates for the trajectory. Needed for nuFFT reconstructions.

BlochSimulators.magnetization_to_signalMethod
magnetization_to_signal(resource, magnetization, parameters, trajectory, coordinates, coil_sensitivities)

Allocates memory for the signal and computes the signal for each coil separately using the _signal_per_coil! function.

Implementation details

The _signal_per_coil! function has different implementations depending on the computational resources (i.e. the type of resource). The default implementations loop over all time points and compute the volume integral of the transverse magnetization in each voxel for each time point separately. This loop order is not necessarily optimal (and performance may be) across all trajectories and computational resources. If a better implementation is available, add new methods to this function for those specific combinations of resources and trajectories.

The "voxels" are assumed to be distributed over the workers. Each worker computes performs a volume integral over the voxels that it owns only (for all time points) using the CPU1() code. The results are then summed up across all workers.

Note

When using multiple CPU's, the "voxels" are distributed over the workers. Each worker computes the signal for its own voxels in parallel and the results are summed up across all workers.

BlochSimulators.magnetization_to_signalMethod
magnetization_to_signal(::Union{CPU1,CPUThreads,CUDALibs}, magnetization, parameters, trajectory::CartesianTrajectory, coordinates, coil_sensitivities)

Arguments

  • magnetization: Matrix{Complex} of size (# readouts, # voxels) with phase-encoded magnetization at echo times.
  • parameters: Tissue parameters of all voxels, including spatial coordinates.
  • trajectory: Cartesian trajectory struct.
  • coordinates: Vector{Coordinates} with spatial coordinates for each voxel.
  • coil_sensitivities: Matrix{Complex} of size (# voxels, # coils) with coil sensitivities.

Returns

  • signal: Vector of length (# coils) with each element a Matrix{Complex} of size (# readouts, # samples per readout)

Extended help

As noted in the description of the simulate_signal function (see src/simulate/signal.jl), we simulate the MR signal at timepoint t from coil i as: signalᵢ[t] = sum(m[t,v] * cᵢ[v] * ρ[v] for v in 1:(# voxels)), where cᵢis the coil sensitivity profile of coil i, ρ is the proton density map and m the matrix with the magnetization at all timepoints for each voxel obtained through Bloch simulations.

The output (signalᵢ) for each coil is in principle a Vector{Complex}of length (# samples per readout) * (# readouts). If we reshape the output into aMatrix{Complex}of size (# samples per readout, # readouts) instead, and do something similar form`, then the signal value associated with the s-th sample point of the r-th readout can be expressed as signalᵢ[r,s] = sum( m[r,s,v]] * cᵢ[v] * ρ[v] for v in 1:(# voxels)).

The problem here is that we typically cannot store the full m. Instead, we compute the magnetization at echo times only. The reason is that, if mᵣ is the magnetization at the r-th echo time in some voxel, and E = exp(-ΔtR₂[v]) * exp(im(Δkₓ*x[v])) is the change per sample point (WHICH FOR CARTESIAN SEQUENCES IS THE SAME FOR ALL READOUTS AND SAMPLES), then the magnetization at the s-th sample relative the the echo time can can be computed as mₛ = mᵣ * E[v]^s

Therefore we can write

signalⱼ[r,s] = sum( magnetization[r,v] * E[v]^s * ρ[v] * cⱼ[v] for v in 1:(# voxels)) signalⱼ[r,s] = magnetization[r,:] * (E.^s .* ρ .* cⱼ)

Because the (E.^s .* ρ .* cⱼ)-part is the same for all readouts, we can simply perform this computation for all readouts simultaneously as signalⱼ[:,s] = magnetization * (E.^s .* ρ .* cⱼ)

If we define the matrix Eˢ as E .^ (-(ns÷2):(ns÷2)-1), then we can do the computation for all different sample points at the same time as well using a single matrix-matrix multiplication: signalⱼ = magnetization * (Eˢ .* (ρ .* cⱼ))

The signalⱼ array is of size (# readouts, # samples per readout). We prefer to have it transposed, therefore we compute signalⱼ = transpose(Eˢ .* (ρ .* cⱼ)) * transpose(magnetization) instead.

For the final output, we do this calculation for each coil j and get a vector of signal matrices (one matrix for each coil) as a result.

Note that this implementation relies entirely on vectorized code and works on both CPU and GPU. The matrix-matrix multiplications are - I think - already multi-threaded so a separate multi-threaded implementation is not needed.

BlochSimulators.make_coordinatesMethod

CoordinatesCollection(x::T, y::T, z::T) where {T<:AbstractArray{<:Real}}

Create a 3D meshgrid of Coordinates from arrays x, y, and z and return it as a StructArray.

Arguments

  • x::T: Array of x-coordinates per voxel.
  • y::T: Array of y-coordinates per voxel.
  • z::T: Array of z-coordinates per voxel.
BlochSimulators.nreadoutsMethod
nreadouts(::AbstractTrajectory)

For each ::AbstractTrajectory, a method should be added to this function that specifies how many readouts the trajectory consists of.

BlochSimulators.nsamplesMethod
nsamples(trajectory::AbstractTrajectory)

Determines the total number of samples acquired with the trajectory. Requires nreadouts and nsamplesperreadout to be implemented.

BlochSimulators.nsamplesperreadoutMethod
nsamplesperreadout(::AbstractTrajectory, readout_idx)

For each ::AbstractTrajectory, a method should be added to this function that specifies how many samples in total are acquired during the trajectory.

BlochSimulators.output_eltypeMethod
output_eltype(::BlochSimulator)

For each <:BlochSimulator, a method should be added to this function that specifies the output type of the simulation. For MR signal simulation, this is typically a complex number representing the transverse magnetization. For other types of simulations, one may want to retrieve the x,y and z components of an isochromat as output (implemented as a FieldVector perhaps) or state configuration matrices Ω.

BlochSimulators.output_sizeMethod
output_size(::BlochSimulator)

For each <:BlochSimulator, a method should be added to this function that specifies the output size of the simulation for a single ::AbstractTissueProperties.

BlochSimulators.phase_encoding!Method
phase_encoding!(magnetization, trajectory::AbstractTrajectory, parameters)

For each ::AbstractTrajectory, a method should be added to this function if it does any kind of phase encoding (so far Cartesian only).

BlochSimulators.regrowthMethod
regrowth(m::Isochromat{T}, E₁) where T

Apply T₁ regrowth to longitudinal component of Isochromat.

BlochSimulators.rotate!Method
rotate!(Ω::EPGStates, eⁱᶿ::T) where T

Rotate F₊ and F̄₋ states under the influence of eⁱᶿ = exp(i * ΔB₀ * Δt)

BlochSimulators.rotateMethod
rotate(m::Isochromat, γΔtGRz, z, Δt, p::AbstractTissueProperties)

Rotation of Isochromat without RF (so around z-axis only) due to gradients and B0 (i.e. refocussing slice select gradient).

BlochSimulators.rotateMethod
rotate(m::Isochromat{T}, γΔtRF::Complex, γΔtGR::Tuple, (x,y,z), Δt, p::AbstractTissueProperties, Δω = zero(T)) where T

RF, gradient and/or ΔB₀ induced rotation of Isochromat computed using Rodrigues rotation formula (https://en.wikipedia.org/wiki/Rodrigues%27rotationformula).

BlochSimulators.sample_transverse!Method
sample!(output, index::Union{Integer,CartesianIndex}, m::Isochromat)

Sample transverse magnetization from Isochromat. The "+=" is needed for 2D sequences where slice profile is taken into account.

BlochSimulators.sample_transverse!Method
sample_transverse!(output, index::Union{Integer,CartesianIndex}, Ω::EPGStates)

Sample the measurable transverse magnetization, that is, the F₊ component of the 0th state. The += is needed for 2D sequences where slice profile is taken into account.

BlochSimulators.sample_xyz!Method
sample_xyz!(output, index::Union{Integer,CartesianIndex}, m::Isochromat)

Sample m.x, m.y and m.z components from Isochromat. The "+=" is needed for 2D sequences where slice profile is taken into account.

BlochSimulators.sample_Ω!Method
sample_Ω!(output, index::Union{Integer,CartesianIndex}, Ω::EPGStates)

Sample the entire configuration state matrix Ω. The += is needed for 2D sequences where slice profile is taken into account.

BlochSimulators.sampling_maskMethod
sampling_mask(tr::CartesianTrajectory)

For undersampled Cartesian trajectories, the gradient trajectory can also be described by a sampling mask.

BlochSimulators.simulate_magnetization!Method
simulate_magnetization!(magnetization, sequence::BlochSimulator, state, p::AbstractTissueProperties) end

For each <:BlochSimulator, a method should be added to this function that implements the actual pulse sequence using information contained in the sequence struct together with the operators from src/operators/{isochromat,epg}.jl. For performance reasons as well as GPU compatibility it is important that the implementation is type-stable and non-allocating.

Arguments

  • magnetization: Pre-allocated array with size(magnetization) = output_size(sequence) and eltype(magnetization) = output_eltype(sequence) to store the output of the simulation.
  • sequence: Sequence struct containing fields that are used to implement the actual pulse sequence.
  • state: Either an Isochromat or EPGStates, depending on which model is used.
  • p: Custom struct (<:AbstractTissueProperties) containing input parameters to the simulation (e.g. T₁T₂)
BlochSimulators.simulate_magnetization!Method
simulate_magnetization!(magnetization, resource, sequence, parameters)

Simulate the magnetization response for all combinations of tissue properties contained in parameters and stores the results in the pre-allocated magnetization array. The actual implementation depends on the computational resource specified in resource.

This function is called by simulate_magnetization and is not intended considered part of te public API.

BlochSimulators.simulate_magnetizationMethod
function simulate_magnetization(sequence, parameters)

Convenience function to simulate magnetization without specifying the computational resource. The function automatically selects the appropriate resource based on the type of the sequence and parameters. The fallback case is to use multi-threaded CPU computations.

BlochSimulators.simulate_magnetizationMethod
simulate_magnetization(resource, sequence, parameters)

Simulate the magnetization response (typically the transverse magnetization at echo times without any spatial encoding gradients applied) for all combinations of tissue properties contained in parameters.

This function can also be used to generate dictionaries for MR Fingerprinting purposes.

Arguments

  • resource::AbstractResource: Either CPU1(), CPUThreads(), CPUProcesses() or CUDALibs()
  • sequence::BlochSimulator: Custom sequence struct
  • parameters::SimulationParameters: Array with different combinations of tissue properties for each voxel.

Note

  • If resource == CUDALibs(), the sequence and parameters must have been moved to the GPU using gpu(sequence) and gpu(parameters) prior to calling this function.
  • If resource == CPUProcesses(), the parameters must be a DArray with the first dimension corresponding to the number of workers. The function will distribute the simulation across the workers in the first dimension of the DArray.

Returns

  • magnetization::AbstractArray: Array of size (output_size(sequence), length(parameters)) containing the magnetization response of the sequence for all combinations of input tissue properties.

Note

If no resource is provided, the simulation is performed on the CPU in a multi-threaded fashion by default. If the parameters are a CuArray, the simulation is performed on the GPU. If the parameters are a DArray, the simulation is performed on the multiple workers.

BlochSimulators.simulate_signalMethod
simulate_signal(resource, sequence, parameters, trajectory, coil_sensitivities)

Simulate the MR signal at timepoint t from coil i as: sᵢ(t) = ∑ⱼ cᵢⱼρⱼmⱼ(t), where cᵢⱼis the coil sensitivity of coil i at position of voxel j, ρⱼ is the proton density of voxel j and mⱼ(t) the (normalized) transverse magnetization in voxel j obtained through Bloch simulations.

Arguments

  • resource::AbstractResource: Either CPU1(), CPUThreads(), CPUProcesses() or CUDALibs()
  • sequence::BlochSimulator: Custom sequence struct
  • parameters::SimulationParameters: Array with tissue properties for each voxel
  • trajectory::AbstractTrajectory: Custom trajectory struct
  • coordinates::StructArray{<:Coordinates}: Array with spatial coordinates for each voxel
  • coil_sensitivities::AbstractMatrix: Sensitivity of coil j in voxel v is given by coil_sensitivities[v,j]

Returns

  • signal::AbstractArray{<:Complex}: Simulated MR signal for the sequence and trajectory. The array is of size (# samples per readout, # readouts, # coils).
BlochSimulators.spoil!Method
spoil!(Ω::EPGStates)

Perfectly spoil the transverse components of all states.

BlochSimulators.to_sample_pointMethod
to_sample_point(m, trajectory, readout_idx, sample_idx, parameters)

For each ::AbstractTrajectory, a method should be added to this function that, given the magnetization m at the readout with index readout_idx, it computes the magnetization at the readout point with index sample_idx (by applying spatial encoding gradients, T₂ decay, B₀ rotation, etc...) based on the trajectory and parameters.

Arguments

  • m: Magnetization at the echo with index readout_idx.
  • trajectory: Trajectory struct containing fields that are used to compute the magnetization at other readout times, including the effects of spatial encoding gradients.
  • readout_idx: Index that corresponds to the current readout.
  • sample_idx: Index for the desired sample during this readout.
  • parameters: Tissue parameters of current voxel, including spatial coordinates.

Output:

  • mₛ: Magnetization at sample with index sample_idx
BlochSimulators.Ω_eltypeMethod
Ω_eltype(sequence::EPGSimulator{T,Ns}) where {T,Ns} = Complex{T}

By default, configuration states are complex. For some sequences, they will only ever be real (no RF phase, no complex slice profile correction) and for these sequences a method needs to be added to this function.

BlochSimulators.@parametersMacro
macro parameters(args...)

Create a SimulationParameters with the actual struct type being determined by the arguments passed to the macro.

Examples

# Create a StructArray{T₁T₂} with T₁ and T₂ values
T₁, T₂ = rand(100), 0.1*rand(100)
parameters = @parameters T₁ T₂

# Create a StructArray{T₁T₂B₁} with T₁, T₂ and B₁ values
T₁, T₂, B₁ = rand(100), 0.1*rand(100), ones(100)
parameters = @parameters T₁ T₂ B₁

# Create a StructArray{T₁T₂B₀} with T₁, T₂ and B₀ values
# This time use the aliases that don't use unicode characters
T1, T2, B0 = rand(100), 0.1*rand(100), ones(100)
parameters = @parameters T1 T2 B0