API
Isochromat Operators
BlochSimulators.Isochromat
— Typestruct 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.decay
— Methoddecay(m::Isochromat{T}, E₁, E₂) where T
Apply T₂ decay to transverse component and T₁ decay to longitudinal component of Isochromat
.
BlochSimulators.initialize_states
— Methodinitialize_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
— Methodinvert(m::Isochromat{T}, p::AbstractTissueProperties) where T
Invert Isochromat
with B₁ insenstive (i.e. adiabatic) inversion pulse
BlochSimulators.invert
— Methodinvert(m::Isochromat{T}, p::AbstractTissueProperties) where T
Invert z-component of Isochromat
(assuming spoiled transverse magnetization so xy-component zero).
BlochSimulators.regrowth
— Methodregrowth(m::Isochromat{T}, E₁) where T
Apply T₁ regrowth to longitudinal component of Isochromat
.
BlochSimulators.rotate
— Methodrotate(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.rotate
— Methodrotate(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!
— Methodsample!(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_xyz!
— Methodsample_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.
EPG Operators
Missing docstring for EPGStates
. Check Documenter's build log for details.
BlochSimulators.F̄₋
— MethodF̄₋(Ω)
View into the second row of the configuration state matrix Ω
, corresponding to the F̄₋
states.
BlochSimulators.F₊
— MethodF₊(Ω)
View into the first row of the configuration state matrix Ω
, corresponding to the F₊
states.
BlochSimulators.Z
— MethodZ(Ω)
View into the third row of the configuration state matrix Ω
, corresponding to the Z
states.
BlochSimulators.decay!
— Methoddecay!(Ω::AbstractConfigurationStates, E₁, E₂)
T₂ decay for F-components, T₁ decay for Z
-component of each state.
BlochSimulators.dephasing!
— Methoddephasing!(Ω::AbstractConfigurationStates)
Shift states around due to dephasing gradient: The F₊
go up one, the F̄₋
go down one and Z
do not change
BlochSimulators.excite!
— Methodexcite!(Ω::AbstractConfigurationStates, 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!
— Methodexcite!(Ω::AbstractConfigurationStates, 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.initial_conditions!
— Methodinitial_conditions!(Ω::AbstractConfigurationStates)
Set all components of all states to 0, except the Z-component of the 0th state which is set to 1.
BlochSimulators.initialize_states
— Methodinitialize_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_states
— Methodinitialize_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.invert!
— Methodinvert!(Ω::AbstractConfigurationStates, p::AbstractTissueProperties)
Invert Z
-component of states of all orders. Assumes fully spoiled transverse magnetization.
BlochSimulators.invert!
— Methodinvert!(Ω::AbstractConfigurationStates)
Invert with B₁ insenstive (i.e. adiabatic) inversion pulse
BlochSimulators.regrowth!
— Methodregrowth!(Ω::AbstractConfigurationStates, E₁)
T₁ regrowth for Z-component of 0th order state.
BlochSimulators.rotate!
— Methodrotate!(Ω::AbstractConfigurationStates, eⁱᶿ::T) where T
Rotate F₊
and F̄₋
states under the influence of eⁱᶿ = exp(i * ΔB₀ * Δt)
BlochSimulators.rotate_decay!
— Methodrotate_decay!(Ω::AbstractConfigurationStates, E₁, E₂, eⁱᶿ)
Rotate and decay combined
BlochSimulators.sample_transverse!
— Methodsample_transverse!(output, index::Union{Integer,CartesianIndex}, Ω::AbstractConfigurationStates)
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_Ω!
— Methodsample_Ω!(output, index::Union{Integer,CartesianIndex}, Ω::AbstractConfigurationStates)
Sample the entire configuration state matrix Ω
. The +=
is needed for 2D sequences where slice profile is taken into account.
BlochSimulators.spoil!
— Methodspoil!(Ω::AbstractConfigurationStates)
Perfectly spoil the transverse components of all states.
BlochSimulators.Ω_eltype
— MethodΩ_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.
Tissue Parameters
BlochSimulators.AbstractTissueProperties
— TypeAbstractTissueProperties{N,T} <: FieldVector{N,T}
Abstract type for custom structs that hold tissue properties used for a simulation within one voxel. For simulations, SimulationParameters
s are used that can be assembled with the @parameters
macro.
Possible fields:
T₁::T
: T₁ relaxation parameters of a voxelT₂::T
: T₂ relaxation parameters of a voxelB₁::T
: Scaling factor for effective B₁ excitation field within a voxelB₀::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:
- FieldVectors/StaticVectors have sizes that are known at compile time. This is beneficial for performance reasons
- The named fields improve readability of the code (e.g.
p.B₁
vsp[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.
Sequences
Abstract Types
BlochSimulators.BlochSimulator
— TypeBlochSimulator{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:
Make a struct that's a subtype of either
IsochromatSimulator
orEPGSimulator
. The struct will hold parameters that are necessary for performing the simulations.Add a method to
simulate_magnetization!
that implements the pulse sequence. For both performance and GPU compatibility, make sure thatsimulate_magnetization!
does not do any heap allocations. Examples forpSSFP
andFISP
sequences are found insrc/sequences
.Add methods to
output_eltype
andoutput_size
that are used to allocate an output array within the simulate function.[Optional] Add a method to show for nicer printing of the sequence in the REPL
[Optional] Add a method to getindex to easily reduce the length of the sequence
[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.jl
and
operators/isochromat.jl` and a properly parametrized sequence struct.
BlochSimulators.IsochromatSimulator
— TypeIsochromatSimulator{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.EPGSimulator
— TypeEPGSimulator{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.
Interface
Examples
BlochSimulators.Generic2D
— TypeGeneric2D{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 intervalGR::M{T}
: Matrix with GRx, GRy and GRz values during each time intervalsample::S
: Vector with Bool's to indicate the sample pointsΔt::V{T}
: Vector with time intervalsz::V{T}
: Vector with different positions along the slice direction
BlochSimulators.Generic3D
— TypeGeneric3D{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 intervalGR::M{T}
: Matrix with GRx, GRy and GRz values during each time intervalsample::S
: Vector with Bool's to indicate the sample pointsΔt::V{T}
: Vector with time intervals
Trajectories
Abstract types
BlochSimulators.AbstractTrajectory
— TypeAbstractTrajectory{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.SpokesTrajectory
— TypeSpokesTrajectory{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.
Interface
Examples
BlochSimulators.CartesianTrajectory2D
— TypeCartesianTrajectory2D{T,I,U,V} <: CartesianTrajectory{T}
Struct that is used to implement a typical Cartesian 2D 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 CartesianTrajectory2D and RadialTrajectory2D are essentially the same 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 trajectorynsamplesperreadout::I
: The total number of samples per readoutΔt::T
: Time between sample pointsk_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 readoutreadout_oversampling::I
: Readout oversampling factor
BlochSimulators.kspace_coordinates
— Methodkspace_coordinates(tr::CartesianTrajectory2D)
Return matrix (nrsamplesperreadout, nrreadouts) with kspace coordinates for the trajectory. Needed for nuFFT reconstructions.
BlochSimulators.magnetization_to_signal
— Methodmagnetization_to_signal(::Union{CPU1,CPUThreads,CUDALibs}, magnetization, parameters, trajectory::CartesianTrajectory2D, 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 a
Matrix{Complex}of size (# samples per readout, # readouts) instead, and do something similar for
m`, 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.sampling_mask
— Methodsampling_mask(tr::CartesianTrajectory2D)
For undersampled Cartesian trajectories, the gradient trajectory can also be described by a sampling mask.
BlochSimulators.RadialTrajectory2D
— TypeRadialTrajectory2D{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 CartesianTrajectory2D and RadialTrajectory2D are essentially the same 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 trajectorynsamplesperreadout::I
: The total number of samples per readoutΔt::T
: Time between sample pointsk_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 readoutreadout_oversampling::I
: Readout oversampling factor
BlochSimulators.add_gradient_delay!
— Methodadd_gradient_delay!(tr::RadialTrajectory2D, 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.kspace_coordinates
— Methodkspace_coordinates(tr::RadialTrajectory2D)
Return matrix (nrsamplesperreadout, nrreadouts) with kspace coordinates for the trajectory. Needed for nuFFT reconstructions.
Dictionary Simulation
Missing docstring for simulate_magnetization(resource, sequence, parameters)
. Check Documenter's build log for details.
BlochSimulators._allocate_array_on_resource
— Method_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_array
— Method_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 specifiedresource
, formatted to store the simulation results for each voxel across the specified echo times.
BlochSimulators.simulate_magnetization!
— Methodsimulate_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_magnetization
— MethodIf the tissue properties for a single voxel are provided only, the simulation is performed on the CPU in a single-threaded fashion.
BlochSimulators.simulate_magnetization
— MethodIf no resource
is provided, the simulation is performed on the CPU in a multi-threaded fashion by default.
BlochSimulators.simulate_magnetization
— MethodHowever, if the parameters
are a CuArray, the simulation is performed on the GPU.
BlochSimulators.simulate_magnetization
— MethodIf the parameters
are a DArray, the simulation is performed on the multiple workers.
BlochSimulators.simulate_magnetization
— Methodfunction 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_magnetization
— Methodsimulate_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
: EitherCPU1()
,CPUThreads()
,CPUProcesses()
orCUDALibs()
sequence::BlochSimulator
: Custom sequence structparameters::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 usinggpu(sequence)
andgpu(parameters)
prior to calling this function. - If
resource == CPUProcesses()
, the parameters must be aDArray
with the first dimension corresponding to the number of workers. The function will distribute the simulation across the workers in the first dimension of theDArray
.
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.
Signal Simulation
Missing docstring for simulate_signal(resource, sequence, parameters, trajectory, coil_sensitivities)
. Check Documenter's build log for details.
BlochSimulators._allocate_signal_array
— Method_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._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.magnetization_to_signal
— Methodmagnetization_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.simulate_signal
— Methodsimulate_signal(sequence, partitioned_parameters::AbstractVector{<:SimulationParameters})
In situations where the number of voxels is too large to store the intermediate magnetization
array, the signal can be calculated in batches: the voxels are divided (by the user) into partitions and the signal is calculated for each partition separately. The final signal is the sum of the signals from all partitions.
BlochSimulators.simulate_signal
— MethodWhen coil sensitivities are not provided, use a single coil with sensitivity = 1 everywhere
BlochSimulators.simulate_signal
— Methodsimulate_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
: EitherCPU1()
,CPUThreads()
,CPUProcesses()
orCUDALibs()
sequence::BlochSimulator
: Custom sequence structparameters::SimulationParameters
: Array with tissue properties for each voxeltrajectory::AbstractTrajectory
: Custom trajectory structcoordinates::StructArray{<:Coordinates}
: Array with spatial coordinates for each voxelcoil_sensitivities::AbstractMatrix
: Sensitivity of coilj
in voxelv
is given bycoil_sensitivities[v,j]
Returns
signal::AbstractArray{<:Complex}
: Simulated MR signal for thesequence
andtrajectory
. The array is of size (# samples per readout, # readouts, # coils).
Utility Functions
BlochSimulators.f32
— Methodf32(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.f64
— Methodf64(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._all_arrays_are_cuarrays
— Method_all_arrays_are_cuarrays(x)
Returns true
if all AbstractArray
fields in x
are CuArray
s and false
otherwise. Will also return
if x
does not have any AbstractArray
fields.
BlochSimulators.gpu
— Methodgpu(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
).
Index
BlochSimulators.AbstractTissueProperties
BlochSimulators.AbstractTrajectory
BlochSimulators.BlochSimulator
BlochSimulators.CartesianTrajectory2D
BlochSimulators.EPGSimulator
BlochSimulators.Generic2D
BlochSimulators.Generic3D
BlochSimulators.Isochromat
BlochSimulators.IsochromatSimulator
BlochSimulators.RadialTrajectory2D
BlochSimulators.SpokesTrajectory
BlochSimulators.F̄₋
BlochSimulators.F₊
BlochSimulators.Z
BlochSimulators._all_arrays_are_cuarrays
BlochSimulators._allocate_array_on_resource
BlochSimulators._allocate_magnetization_array
BlochSimulators._allocate_signal_array
BlochSimulators._signal_per_coil!
BlochSimulators.add_gradient_delay!
BlochSimulators.decay
BlochSimulators.decay!
BlochSimulators.dephasing!
BlochSimulators.excite!
BlochSimulators.excite!
BlochSimulators.f32
BlochSimulators.f64
BlochSimulators.gpu
BlochSimulators.initial_conditions!
BlochSimulators.initialize_states
BlochSimulators.initialize_states
BlochSimulators.initialize_states
BlochSimulators.invert
BlochSimulators.invert
BlochSimulators.invert!
BlochSimulators.invert!
BlochSimulators.kspace_coordinates
BlochSimulators.kspace_coordinates
BlochSimulators.magnetization_to_signal
BlochSimulators.magnetization_to_signal
BlochSimulators.regrowth
BlochSimulators.regrowth!
BlochSimulators.rotate
BlochSimulators.rotate
BlochSimulators.rotate!
BlochSimulators.rotate_decay!
BlochSimulators.sample_transverse!
BlochSimulators.sample_transverse!
BlochSimulators.sample_xyz!
BlochSimulators.sample_Ω!
BlochSimulators.sampling_mask
BlochSimulators.simulate_magnetization
BlochSimulators.simulate_magnetization
BlochSimulators.simulate_magnetization
BlochSimulators.simulate_magnetization
BlochSimulators.simulate_magnetization
BlochSimulators.simulate_magnetization
BlochSimulators.simulate_magnetization!
BlochSimulators.simulate_signal
BlochSimulators.simulate_signal
BlochSimulators.simulate_signal
BlochSimulators.spoil!
BlochSimulators.Ω_eltype