KomaMRIBase
Scanner-related functions
KomaMRIBase.Scanner — Typesys = Scanner(B0, B1, Gmax, Smax, ADC_Δt, seq_Δt, GR_Δt, RF_Δt,
RF_ring_down_T, RF_dead_time_T, ADC_dead_time_T)The Scanner struct. It contains hardware limitations of the MRI resonator. It is an input for the simulation.
Arguments
B0: (::Real,=1.5,[T]) main magnetic field strengthB1: (::Real,=10e-6,[T]) maximum RF amplitudeGmax: (::Real,=60e-3,[T/m]) maximum gradient amplitudeSmax: (::Real,=500,[mT/m/ms]) gradient's maximum slew-rateADC_Δt: (::Real,=2e-6,[s]) ADC raster timeseq_Δt: (::Real,=1e-5,[s]) sequence-block raster timeGR_Δt: (::Real,=1e-5,[s]) gradient raster timeRF_Δt: (::Real,=1e-6,[s]) RF raster timeRF_ring_down_T: (::Real,=20e-6,[s]) RF ring down timeRF_dead_time_T: (::Real,=100e-6,[s]) RF dead timeADC_dead_time_T: (::Real,=10e-6,[s]) ADC dead time
Returns
sys: (::Scanner) Scanner struct
Examples
julia> sys = Scanner()
julia> sys.B0Phantom-related functions
KomaMRIBase.Phantom — Typeobj = Phantom(name, x, y, z, ρ, T1, T2, T2s, Δw, Dλ1, Dλ2, Dθ, motion)The Phantom struct. Most of its field names are vectors, with each element associated with a property value representing a spin. This struct serves as an input for the simulation.
Arguments
name: (::String) phantom namex: (::AbstractVector{T<:Real},[m]) spin x-position vectory: (::AbstractVector{T<:Real},[m]) spin y-position vectorz: (::AbstractVector{T<:Real},[m]) spin z-position vectorρ: (::AbstractVector{T<:Real}) spin proton density vectorT1: (::AbstractVector{T<:Real},[s]) spin T1 parameter vectorT2: (::AbstractVector{T<:Real},[s]) spin T2 parameter vectorT2s: (::AbstractVector{T<:Real},[s]) spin T2s parameter vectorΔw: (::AbstractVector{T<:Real},[rad/s]) spin off-resonance parameter vectorDλ1: (::AbstractVector{T<:Real}) spin Dλ1 (diffusion) parameter vectorDλ2: (::AbstractVector{T<:Real}) spin Dλ2 (diffusion) parameter vectorDθ: (::AbstractVector{T<:Real}) spin Dθ (diffusion) parameter vectormotion: (::MotionModel{T<:Real}) motion model
Returns
obj: (::Phantom) Phantom struct
Examples
julia> obj = Phantom(x=[0.0])
julia> obj.ρKomaMRIBase.brain_phantom2D — Functionphantom = brain_phantom2D(;axis="axial", ss=4)Creates a two-dimensional brain Phantom struct. Default ss=4 sample spacing is 2 mm. Original file (ss=1) sample spacing is .5 mm.
References
- B. Aubert-Broche, D.L. Collins, A.C. Evans: "A new improved version of the realistic digital brain phantom" NeuroImage, in review - 2006
- B. Aubert-Broche, M. Griffin, G.B. Pike, A.C. Evans and D.L. Collins: "20 new digital brain phantoms for creation of validation image data bases" IEEE TMI, in review - 2006
- https://brainweb.bic.mni.mcgill.ca/brainweb
Keywords
axis: (::String,="axial", opts=["axial","coronal","sagittal"]) orientation of the phantomss: (::Integer or ::Vector{Integer},=4) subsampling parameter for all axes if scaler, per axis if 2 element vector [ssx, ssy]us: (::Integer or ::Vector{Integer},=1) upsampling parameter for all axes if scaler, per axis if 2 element vector [usx, usy], if used ss is set to ss=1
Returns
obj: (::Phantom) Phantom struct
Examples
julia> obj = brain_phantom2D(; axis="sagittal", ss=1)
julia> obj = brain_phantom2D(; axis="axial", us=[1, 2])
julia> plot_phantom_map(obj, :ρ)KomaMRIBase.brain_phantom3D — Functionobj = brain_phantom3D(; ss=4, us=1)Creates a three-dimentional brain Phantom struct. Default ss=4 sample spacing is 2 mm. Original file (ss=1) sample spacing is .5 mm.
References
- B. Aubert-Broche, D.L. Collins, A.C. Evans: "A new improved version of the realistic digital brain phantom" NeuroImage, in review - 2006
- B. Aubert-Broche, M. Griffin, G.B. Pike, A.C. Evans and D.L. Collins: "20 new digital brain phantoms for creation of validation image data bases" IEEE TMI, in review - 2006
- https://brainweb.bic.mni.mcgill.ca/brainweb
Keywords
ss: (::Integer or ::Vector{Integer},=4) subsampling parameter for all axes if scaler, per axis if 3 element vector [ssx, ssy, ssz]us: (::Integer or ::Vector{Integer},=1) upsampling parameter for all axes if scaler, per axis if 3 element vector [usx, usy, usz]start_end: (::Vector{Integer},=[160,200]) z index range of presampled phantom, 180 is center
Returns
obj: (::Phantom) 3D Phantom struct
Examples
julia> obj = brain_phantom3D(; ss=5)
julia> obj = brain_phantom3D(; us=[2, 2, 1])
julia> plot_phantom_map(obj, :ρ)KomaMRIBase.pelvis_phantom2D — Functionobj = pelvis_phantom2D(; ss=4, us=1)Creates a two-dimensional pelvis Phantom struct. Default ss=4 sample spacing is 2 mm. Original file (ss=1) sample spacing is .5 mm.
Keywords
ss: (::Integer or ::Vector{Integer},=4) subsampling parameter for all axes if scaler, per axis if 2 element vector [ssx, ssy]us: (::Integer or ::Vector{Integer},=1) upsampling parameter for all axes if scaler, per axis if 2 element vector [usx, usy]
Returns
obj: (::Phantom) Phantom struct
Examples
julia> obj = pelvis_phantom2D(; ss=2])
julia> obj = pelvis_phantom2D(; us=[1, 2])
julia> pelvis_phantom2D(obj, :ρ)KomaMRIBase.heart_phantom — Functionobj = heart_phantom(...)Heart-like LV 2D phantom. The variable circumferential_strain and radial_strain are for streching (if positive) or contraction (if negative). rotation_angle is for rotation.
Arguments
circumferential_strain: (::Real,=-0.3) contraction parameterradial_strain: (::Real,=-0.3) contraction parameterrotation_angle: (::Real,=1) rotation parameter
Returns
phantom: (::Phantom) Heart-like LV phantom struct
MotionModel-related functions
KomaMRIBase.get_spin_coords — Functionx, y, z = get_spin_coords(motion, x, y, z, t)Calculates the position of each spin at a set of arbitrary time instants, i.e. the time steps of the simulation. For each dimension (x, y, z), the output matrix has $N_{ ext{spins}}$ rows and length(t) columns.
Arguments
motion: (::MotionModel) phantom motionx: (::AbstractVector{T<:Real},[m]) spin x-position vectory: (::AbstractVector{T<:Real},[m]) spin y-position vectorz: (::AbstractVector{T<:Real},[m]) spin z-position vectort: (::AbstractArray{T<:Real}) horizontal array of time instants
Returns
x, y, z: (::Tuple{AbstractArray, AbstractArray, AbstractArray}) spin positions over time
SimpleMotion <: MotionModel
KomaMRIBase.SimpleMotion — Typemotion = SimpleMotion(types)SimpleMotion model. It allows for the definition of motion by means of simple parameters. The SimpleMotion struct is composed by only one field, called types, which is a tuple of simple motion types. This tuple will contain as many elements as simple motions we want to combine.
Arguments
types: (::Tuple{Vararg{<:SimpleMotionType{T}}}) tuple of simple motion types
Returns
motion: (::SimpleMotion) SimpleMotion struct
Examples
julia> motion = SimpleMotion(
Translation(dx=0.01, dy=0.02, dz=0.0, t_start=0.0, t_end=0.5),
Rotation(pitch=15.0, roll=0.0, yaw=20.0, t_start=0.1, t_end=0.5),
HeartBeat(circumferential_strain=-0.3, radial_strain=-0.2, longitudinal_strain=0.0, t_start=0.2, t_end=0.5)
)SimpleMotion types
KomaMRIBase.Translation — Typetranslation = Translation(dx, dy, dz, t_start, t_end)Translation motion struct. It produces a linear translation of the phantom. Its fields are the final displacements in the three axes (dx, dy, dz) and the start and end times of the translation.
Arguments
dx: (::Real,[m]) translation in xdy: (::Real,[m]) translation in ydz: (::Real,[m]) translation in zt_start: (::Real,[s]) initial timet_end: (::Real,[s]) final time
Returns
translation: (::Translation) Translation struct
Examples
julia> tr = Translation(dx=0.01, dy=0.02, dz=0.03, t_start=0.0, t_end=0.5)KomaMRIBase.Rotation — Typerotation = Rotation(pitch, roll, yaw, t_start, t_end)Rotation motion struct. It produces a rotation of the phantom in the three axes: x (pitch), y (roll), and z (yaw). We follow the RAS (Right-Anterior-Superior) orientation, and the rotations are applied following the right-hand rule (counter-clockwise):
The applied rotation matrix is obtained as follows:
\[\begin{equation} \begin{aligned} R &= R_z(\alpha) R_y(\beta) R_x(\gamma) \\ &= \begin{bmatrix} \cos \alpha & -\sin \alpha & 0 \\ \sin \alpha & \cos \alpha & 0 \\ 0 & 0 & 1 \end{bmatrix} \begin{bmatrix} \cos \beta & 0 & \sin \beta \\ 0 & 1 & 0 \\ -\sin \beta & 0 & \cos \beta \end{bmatrix} \begin{bmatrix} 1 & 0 & 0 \\ 0 & \cos \gamma & -\sin \gamma \\ 0 & \sin \gamma & \cos \gamma \end{bmatrix} \\ &= \begin{bmatrix} \cos \alpha \cos \beta & \cos \alpha \sin \beta \sin \gamma - \sin \alpha \cos \gamma & \cos \alpha \sin \beta \cos \gamma + \sin \alpha \sin \gamma \\ \sin \alpha \cos \beta & \sin \alpha \sin \beta \sin \gamma + \cos \alpha \cos \gamma & \sin \alpha \sin \beta \cos \gamma - \cos \alpha \sin \gamma \\ -\sin \beta & \cos \beta \sin \gamma & \cos \beta \cos \gamma \end{bmatrix} \end{aligned} \end{equation}\]
Arguments
pitch: (::Real,[º]) rotation in xroll: (::Real,[º]) rotation in yyaw: (::Real,[º]) rotation in zt_start: (::Real,[s]) initial timet_end: (::Real,[s]) final time
Returns
rotation: (::Rotation) Rotation struct
Examples
julia> rt = Rotation(pitch=15.0, roll=0.0, yaw=20.0, t_start=0.1, t_end=0.5)KomaMRIBase.HeartBeat — Typeheartbeat = HeartBeat(circumferential_strain, radial_strain, longitudinal_strain, t_start, t_end)HeartBeat struct. It produces a heartbeat-like motion, characterised by three types of strain: Circumferential, Radial and Longitudinal
Arguments
circumferential_strain: (::Real) contraction parameterradial_strain: (::Real) contraction parameterlongitudinal_strain: (::Real) contraction parametert_start: (::Real,[s]) initial timet_end: (::Real,[s]) final time
Returns
heartbeat: (::HeartBeat) HeartBeat struct
Examples
julia> hb = HeartBeat(circumferential_strain=-0.3, radial_strain=-0.2, longitudinal_strain=0.0, t_start=0.2, t_end=0.5)KomaMRIBase.PeriodicTranslation — Typeperiodic_translation = PeriodicTranslation(dx, dy, dz, period, asymmetry)PeriodicTranslation motion struct. It produces a periodic translation of the phantom in the three directions x, y and z. The amplitude of the oscillation will be defined by dx, dy and dz
Arguments
dx: (::Real,[m]) translation in xdy: (::Real,[m]) translation in ydz: (::Real,[m]) translation in zperiod: (::Real,[s]) periodasymmetry: (::Real) asymmetry factor, between 0 and 1
Returns
periodic_translation: (::PeriodicTranslation) PeriodicTranslation struct
KomaMRIBase.PeriodicRotation — Typeperiodic_rotation = PeriodicRotation(pitch, roll, yaw, period, asymmetry)PeriodicRotation motion struct. It produces a rotation of the phantom in the three axes: x (pitch), y (roll), and z (yaw)
Arguments
pitch: (::Real,[º]) rotation in xroll: (::Real,[º]) rotation in yyaw: (::Real,[º]) rotation in zperiod: (::Real,[s]) periodasymmetry: (::Real) asymmetry factor, between 0 and 1
Returns
periodic_rotation: (::PeriodicRotation) PeriodicRotation struct
KomaMRIBase.PeriodicHeartBeat — Typeperiodic_heartbeat = PeriodicHeartBeat(circumferential_strain, radial_strain, longitudinal_strain, period, asymmetry)HeartBeat struct. It produces a heartbeat-like motion, characterised by three types of strain: Circumferential, Radial and Longitudinal
Arguments
circumferential_strain: (::Real) contraction parameterradial_strain: (::Real) contraction parameterlongitudinal_strain: (::Real) contraction parameterperiod: (::Real,[s]) periodasymmetry: (::Real) asymmetry factor, between 0 and 1
Returns
periodic_heartbeat: (::PeriodicHeartBeat) PeriodicHeartBeat struct
ArbitraryMotion <: MotionModel
KomaMRIBase.ArbitraryMotion — Typemotion = ArbitraryMotion(t_start, t_end, dx, dy, dz)ArbitraryMotion model. For this motion model, it is necessary to define motion for each spin independently, in x (dx), y (dy) and z (dz). dx, dy and dz are three matrixes, of ($N_{spins}$ x $N_{discrete\,times}$) each. This means that each row corresponds to a spin trajectory over a set of discrete time instants. The motion duration is determined by t_start and t_end. The discrete time instants are evenly spaced, ($dt = rac{t_{end} - t_{start}}{N_{discrete\,times}}$).
This motion model is useful for defining arbitrarly complex motion, specially for importing the spin trajectories from another source, like XCAT or a CFD.
Arguments
t_start: (::T,[s])t_end: (::T,[s])dx: (::Array{T,2},[m]) matrix for displacements in xdy: (::Array{T,2},[m]) matrix for displacements in ydz: (::Array{T,2},[m]) matrix for displacements in z
Returns
motion: (::ArbitraryMotion) ArbitraryMotion struct
Examples
julia> motion = ArbitraryMotion(
0.0,
1.0,
0.01.*rand(1000, 10),
0.01.*rand(1000, 10),
0.01.*rand(1000, 10)
)Sequence-related functions
KomaMRIBase.Sequence — Typeseq = Sequence()
seq = Sequence(GR)
seq = Sequence(GR, RF)
seq = Sequence(GR, RF, ADC)
seq = Sequence(GR, RF, ADC, DUR)
seq = Sequence(GR::Array{Grad,1})
seq = Sequence(GR::Array{Grad,1}, RF::Array{RF,1})
seq = Sequence(GR::Array{Grad,1}, RF::Array{RF,1}, A::ADC, DUR, DEF)The Sequence struct. It contains events of an MRI sequence. Most field names (except for the DEF field) consist of matrices or vectors, where each column index represents a sequence block. This struct serves as an input for the simulation.
Arguments
GR: (::Matrix{Grad}) gradient matrix. Rows for x-y-z amplitudes and columns are for blocksRF: (::Matrix{RF}) RF matrix. The 1 row is for the coil and columns are for blocksADC: (::Array{ADC,1}) ADC block vectorDUR: (::Vector,[s]) duration block vectorDEF: (::Dict{String, Any}) dictionary with relevant information of the sequence. Possible keys could be ["AdcRasterTime","GradientRasterTime","Name","Nz","Num_Blocks","Nx","Ny","PulseqVersion","BlockDurationRaster","FileName","RadiofrequencyRasterTime"]
Returns
seq: (::Sequence) Sequence struct
KomaMRIBase.dur — Functiony = dur(x::Grad)
y = dur(x::Vector{Grad})
y = dur(x::Matrix{Grad})Duration time in [s] of Grad struct or Grad Array.
Arguments
x: (::Grador::Vector{Grad}or::Matrix{Grad}) Grad struct or Grad Array
Returns
y: (::Float64,[s]) duration of the Grad struct or Grad Array
y = dur(x::RF)
y = dur(x::Vector{RF})
y = dur(x::Matrix{RF})Duration time in [s] of RF struct or RF Array.
Arguments
x: (::RFor::Vector{RF}or::Matrix{RF}) RF struct or RF array
Returns
y: (::Float64, [s]) duration of the RF struct or RF array
T = dur(x::Sequence)The total duration of the sequence in [s].
Arguments
x: (::Sequence) Sequence struct
Returns
T: (::Real,[s]) total duration of the sequence
KomaMRIBase.get_block_start_times — FunctionT0 = get_block_start_times(seq::Sequence)Returns a vector containing the start times of blocks in a sequence. The initial time is always zero, and the final time corresponds to the duration of the sequence.
Arguments
seq: (::Sequence) Sequence struct
Returns
T0: (::Vector,[s]) start times of the blocks in a sequence
KomaMRIBase.get_flip_angles — Functiony = get_flip_angles(x::Sequence)Returns all the flip angles of the RF pulses in the sequence x.
Arguments
x: (::Sequence) Sequence struct
Returns
y: (::Vector{Float64},[deg]) flip angles
Grad
KomaMRIBase.Grad — Typegr = Grad(A, T)
gr = Grad(A, T, rise)
gr = Grad(A, T, rise, delay)
gr = Grad(A, T, rise, fall, delay)
gr = Grad(A, T, rise, fall, delay, first, last)The Grad struct represents a gradient of a sequence event.
Arguments
A: (::Realor::Vector,[T/m]) amplitude of the gradientT: (::Realor::Vector,[s]) duration of the flat-toprise: (::Real,[s]) duration of the risefall: (::Real,[s]) duration of the falldelay: (::Real,[s]) duration of the delay
Returns
gr: (::Grad) gradient struct
Examples
julia> gr = Grad(1, 1, 0.1, 0.1, 0.2)
julia> seq = Sequence([gr]); plot_seq(seq)KomaMRIBase.Grad — Methodgr = Grad(f::Function, T::Real, N::Integer; delay::Real)Generates an arbitrary gradient waveform defined by the function f in the interval t ∈ [0,T]. The time separation between two consecutive samples is given by T/(N-1).
Arguments
f: (::Function) function that describes the gradient waveformT: (::Real,[s]) duration of the gradient waveformN: (::Integer,=300) number of samples of the gradient waveform
Keywords
delay: (::Real,=0,[s]) delay time of the waveform
Returns
gr: (::Grad) gradient struct
Examples
julia> gx = Grad(t -> sin(π*t / 0.8), 0.8)
julia> seq = Sequence([gx]); plot_seq(seq)RF
KomaMRIBase.RF — Typerf = RF(A, T)
rf = RF(A, T, Δf)
rf = RF(A, T, Δf, delay)The RF struct represents a Radio Frequency excitation of a sequence event.
Arguments
A: (::Complex,[T]) RF complex amplitud modulation (AM), $B_1(t) = |B_1(t)| e^{i\phi(t)} = B_{1}(t) + iB_{1,y}(t)$T: (::Real, [s]) RF durationΔf: (::Realor::Vector, [Hz]) RF frequency difference with respect to the Larmor frequency. This can be a number but also a vector to represent frequency modulated signals (FM).delay: (::Real, [s]) RF delay time
Returns
rf: (::RF) the RF struct
Examples
julia> rf = RF(1, 1, 0, 0.2)
julia> seq = Sequence(); seq += rf; plot_seq(seq)KomaMRIBase.RF — Methodrf = RF_fun(f::Function, T::Real, N::Int64)Generate an RF sequence with amplitudes sampled from a function waveform.
This function is not being used in this KomaMRI version.
Arguments
f: (::Function, [T]) function for the RF amplitud waveformT: (::Real, [s]) duration of the RF pulseN: (::Int64) number of samples of the RF pulse
Returns
rf:(::RF) RF struct with amplitud defined by the functionf
KomaMRIBase.get_flip_angle — Functionα = get_flip_angle(x::RF)Calculates the flip angle α [deg] of an RF struct. α = γ ∫ B1(τ) dτ
Arguments
x: (::RF) RF struct
Returns
α: (::Int64,[deg]) flip angle RF structx
ADC
KomaMRIBase.ADC — Typeadc = ADC(N, T)
adc = ADC(N, T, delay)
adc = ADC(N, T, delay, Δf, ϕ)The ADC struct represents the Analog to Digital Converter (ADC) of a sequence event.
Arguments
N: (::Int64) number of acquired samplesT: (::Float64, [s]) duration to acquire the samplesdelay: (::Float64, [s]) delay time to start the acquisitionΔf: (::Float64, [Hz]) delta frequency. It is meant to compensate RF pulse phasesϕ: (::Float64,[rad]) phase. It is meant to compensate RF pulse phases
Returns
adc: (::ADC) ADC struct
Examples
julia> adc = ADC(16, 1, 0.1)
julia> seq = Sequence(); seq += adc; plot_seq(seq)KomaMRIBase.get_adc_sampling_times — Functiontimes = get_adc_sampling_times(seq)Returns an array of times when the samples of the sequence seq are acquired.
Arguments
seq: (::Sequence) sequence struct
Returns
times: (::Vector{Float64},[s]) time array when samples are acquired
KomaMRIBase.get_adc_phase_compensation — Functioncomp = get_adc_phase_compensation(seq)Returns an array of phase compensation factors, $\exp(-\mathrm{i}\varphi)$, which are used to compensate the acquired signal $S$ by applying the operation $S_{\mathrm{comp}} = S \exp(-\mathrm{i}\varphi)$ after the simulation. This compensation is necessary because the signal typically exhibits a phase offset of $\varphi$ following RF excitation with a phase of $\varphi$. Such pulses are commonly employed in sequences involving RF spoiling.
Arguments
seq: (::Sequence) sequence struct
Returns
comp: (::Vector{Complex},[rad]) array of phase compensations for every acquired sample
Delay
KomaMRIBase.Delay — Typedelay = Delay(T)The Delay struct is meant to add a delay to a sequence by using a sum operator.
Arguments
T: (::Real,[s]) time delay value
Returns
delay: (::Delay) delay struct
Examples
julia> delay = Delay(0.5)
julia> s = Sequence([Grad(1, 1, 0.1)])
julia> seq = delay + s; plot_seq(seq)Rotation matrices
KomaMRIBase.rotx — FunctionRx = rotx(θ::Real)Rotates vector counter-clockwise with respect to the x-axis.
Arguments
θ: (::Real,[rad]) rotation angle
Returns
Rx: (::Matrix{Int64}) rotation matrix
KomaMRIBase.roty — FunctionRy = roty(θ::Real)Rotates vector counter-clockwise with respect to the y-axis.
Arguments
θ: (::Real,[rad]) rotation angle
Returns
Ry: (::Matrix{Int64}) rotation matrix
KomaMRIBase.rotz — FunctionRz = rotz(θ::Real)Rotates vector counter-clockwise with respect to the z-axis.
Arguments
θ: (::Real,[rad]) rotation angle
Returns
Rz: (::Matrix{Int64}) rotation matrix
Moments
KomaMRIBase.get_Mk — FunctionMk, Mk_adc = get_Mk(seq::Sequence, k; Δt=1, skip_rf=zeros(Bool, sum(is_RF_on.(seq))))Computes the $k$th-order moment of the Sequence seq given by the formula $\int_0^T t^k G(t) dt$.
Arguments
seq: (::Sequence) Sequence structk: (::Integer) order of the moment to be computedΔt: (::Real,=1,[s]) nominal delta time separation between two time samples for ADC acquisition and Gradientsskip_rf: (::Vector{Bool},=zeros(Bool, sum(is_RF_on.(seq)))) boolean vector which indicates whether to skip the computation for resetting the integral for excitation or refocusing RF type
Returns
Mk: (3-column ::Matrix{Real}) $k$th-order momentMk_adc: (3-column ::Matrix{Real}) $k$th-order moment sampled at ADC times
KomaMRIBase.get_kspace — FunctionKomaMRIBase.get_M0 — FunctionComputes the zero-order moment of the Sequence seq. Refer to get_Mk and get_kspace
KomaMRIBase.get_M1 — FunctionComputes the 1st-order moment of the Sequence seq. Refer to get_Mk
KomaMRIBase.get_M2 — FunctionComputes the 2nd-order moment of the Sequence seq. Refer to get_Mk
Event checks
KomaMRIBase.is_RF_on — Functiony = is_RF_on(x::Sequence)
y = is_RF_on(x::Sequence, t::Vector{Float64})Tells if the sequence seq has elements with RF active, or active during time t.
Arguments
x: (::Sequence) Sequence structt: (::Vector{Float64},[s]) time to check
Returns
y: (::Bool) boolean that tells whether or not the RF in the sequence is active
KomaMRIBase.is_GR_on — Functiony = is_GR_on(x::Sequence)Tells if the sequence seq has elements with GR active.
Arguments
x: (::Sequence) Sequence struct
Returns
y: (::Bool) boolean that tells whether or not the GR in the sequence is active
KomaMRIBase.is_Gx_on — Functiony = is_Gx_on(x::Sequence)Tells if the sequence seq has elements with GR active in x direction.
Arguments
x: (::Sequence) Sequence struct
Returns
y: (::Bool) boolean that tells whether or not the GRx in the sequence is active
KomaMRIBase.is_Gy_on — Functiony = is_Gy_on(x::Sequence)Tells if the sequence seq has elements with GR active in y direction.
Arguments
x: (::Sequence) Sequence struct
Returns
y: (::Bool) boolean that tells whether or not the GRy in the sequence is active
KomaMRIBase.is_Gz_on — Functiony = is_Gz_on(x::Sequence)Tells if the sequence seq has elements with GR active in z direction.
Arguments
x: (::Sequence) Sequence struct
Returns
y: (::Bool) boolean that tells whether or not the GRz in the sequence is active
KomaMRIBase.is_ADC_on — Functiony = is_ADC_on(x::Sequence)
y = is_ADC_on(x::Sequence, t::Union{Array{Float64,1}, Array{Float64,2}})Tells if the sequence seq has elements with ADC active, or active during time t.
Arguments
x: (::Sequence) sequence structt: (::Union{Array{Float64,1}, Array{Float64,2}},[s]) time to check
Returns
y: (::Bool) boolean that tells whether or not the ADC in the sequence is active
DiscreteSequence
KomaMRIBase.DiscreteSequence — Typeseqd = DiscreteSequence(Gx, Gy, Gz, B1, Δf, ADC, t, Δt)A sampled version of a Sequence struct, containing vectors for event amplitudes at specified times. DiscreteSequence is the struct used for simulation.
Arguments
Gx: (::AbstractVector{T<:Real},[T/m]) x-gradient vectorGy: (::AbstractVector{T<:Real},[T/m]) y-gradient vectorGz: (::AbstractVector{T<:Real},[T/m]) z-gradient vectorB1: (::AbstractVector{Complex{T<:Real}},[T]) RF amplitude vectorΔf: (::AbstractVector{T<:Real},[Hz]) RF carrier frequency displacement vectorADC: (::AbstractVector{Bool}) ADC sample vectort: (::AbstractVector{T<:Real},[s]) time vectorΔt: (::AbstractVector{T<:Real},[s]) delta time vector
Returns
seqd: (::DiscreteSequence) DiscreteSequence struct
KomaMRIBase.discretize — Functionseqd = discretize(seq::Sequence; sampling_params=default_sampling_params())This function returns a sampled Sequence struct with RF and gradient time refinements based on simulation parameters.
Arguments
seq: (::Sequence) sequence
Keywords
sampling_params: (::Dict{String, Any},=default_sampling_params()) sampling parameter dictionary
Returns
seqd: (::DiscreteSequence) DiscreteSequence struct
KomaMRIBase.get_samples — Functionsamples = get_samples(seq::Sequence; off_val=0, max_rf_samples=Inf)Returns the samples of the events in seq.
Arguments
seq: (::Sequence) Sequence struct
Keywords
off_val: (::Number,=0) offset value for amplitude. Typically used to hide points in plots by setting it toInfmax_rf_samples: (::Integer,=Inf) maximum number of samples for the RF struct
Returns
samples: (::NamedTuple) contains samples forgx,gy,gz,rf, andadcevents. Each event, represented bye::NamedTuple, includes time samples (e.t) and amplitude samples (e.A)
KomaMRIBase.times — Functiont = times(gr::Grad)
t = times(rf::RF)
t = times(adc::ADC)Get time samples of MRI sequence event.
Arguments
gr: (::Grad) Gradient structrf: (::RF) RF structadc: (::ADC) ADC struct
Returns
t: (::Vector{Number}) vector with time samples
limits = times(obj.motion)KomaMRIBase.ampls — FunctionA = ampls(g::Grad)
A = ampls(r::RF)
A = ampls(d::ADC)Get amplitude samples of MRI sequence event.
Arguments
gr: (::Grad) Gradient structrf: (::RF) RF structadc: (::ADC) ADC struct
Returns
A: (::Vector{Number}) vector with amplitude samples
KomaMRIBase.freqs — Functionf = freqs(r::RF)Get frequency samples of MRI sequence event.
Arguments
rf: (::RF) RF struct
Returns
f: (::Vector{Number}) vector with frequency samples
Other functions
KomaMRIBase.trapz — Functiony = trapz(Δt, x)Trapezoidal integration for every spin of a phantom.
In practice, this function is used to integrate (Gx * x + Gy * y + Gz * z) * Δt for all the spins. NΔt is the length of Δt. Ns stands for the number of spins of a phantom. x is a matrix which rows represents different spins and columns are different times and the elements are the field Gx * x + Gy * y + Gz * z values.
Arguments
Δt: (1 x NΔt ::Matrix{Float64},[s]) delta time 1-row arrayx: (Ns x (NΔt+1) ::Matrix{Float64},[T]) magnitude of the field Gx * x + Gy * y + Gz * z
Returns
y: (Ns x 1 ::Matrix{Float64},[T*s]) vector where every element is the integral of (Gx * x + Gy * y + Gz * z) * Δt for every spin of a phantom
KomaMRIBase.cumtrapz — Functiony = cumtrapz(Δt, x)Trapezoidal cumulative integration over time for every spin of a phantom.
Arguments
Δt: (1 x NΔt ::Matrix{Float64},[s]) delta time 1-row arrayx: (Ns x (NΔt+1) ::Matrix{Float64},[T]) magnitude of the field Gx * x + Gy * y + Gz * z
Returns
y: (Ns x NΔt ::Matrix{Float64},[T*s]) matrix where every column is the cumulative integration over time of (Gx * x + Gy * y + Gz * z) * Δt for every spin of a phantom
cumtrapzA more efficient GPU implementation of the cumtrapz method defined in TrapezoidalIntegration.jl. Uses a kernel to compute cumsum along the second dimension.
Arguments
Δt: (1 x NΔt ::Matrix{Float64},[s]) delta time 1-row arrayx: (Ns x (NΔt+1) ::Matrix{Float64},[T]) magnitude of the field Gx * x + Gy * y + Gz * z
Returns
y: (Ns x NΔt ::Matrix{Float64},[T*s]) matrix where every column is the cumulative integration over time of (Gx * x + Gy * y + Gz * z) * Δt for every spin of a phantom
KomaMRIBase.kfoldperm — Functionarray_of_ranges = kfoldperm(N, k; breaks=[])Divides a list of indices from 1 to N into k groups.
Arguments
N: (::Integer) number of elements to be orderedk: (::Integer) number of groups to divide theNelements.
Keywords
breaks: (::Vector{<:Integer},=[]) array of indices where predefined breakpoints are placed.
Returns
array_of_ranges: (::Vector{UnitRange{<:Integer}}) array containing ranges of different groups. The target iskgroups, but this could increase by adding elements to thebreaksinput array
Sequence Building Blocks (SBB)
KomaMRIBase.PulseDesigner — ModulePulseDesignerA module to define different pulse sequences.
KomaMRIBase.PulseDesigner.RF_hard — Functionseq = RF_hard(B1, T, sys; G=[0, 0, 0], Δf=0)Returns a sequence with a RF excitation pulse.
Arguments
B1: (::Number,[T]) RF pulse amplitudeT: (::Real,[s]) RF pulse durationsys: (::Scanner) Scanner struct
Keywords
G: (::Vector{Real},=[0, 0, 0],[T/m]) gradient amplitudes for x, y, zΔf: (::Real,=0,[Hz]) RF pulse carrier frequency displacement
Returns
seq: (::Sequence) Sequence struct with a RF pulse
Examples
julia> sys = Scanner(); durRF = π / 2 / (2π * γ * sys.B1);
julia> seq = PulseDesigner.RF_hard(sys.B1, durRF, sys);
julia> plot_seq(seq)KomaMRIBase.PulseDesigner.RF_sinc — Functionseq = RF_sinc(B1, T, sys; G=[0, 0, 0], Δf=0, a=0.46, TBP=4)Returns a sequence with a RF sinc waveform.
References
- Matt A. Bernstein, Kevin F. King, Xiaohong Joe Zhou, Chapter 2 - Radiofrequency Pulse
Shapes, Handbook of MRI Pulse Sequences, 2004, Pages 35-66, https://doi.org/10.1016/B978-012092861-3/50006-6.
Arguments
B1: (::Number,[T]) RF sinc amplitudeT: (::Real,[s]) RF sinc durationsys: (::Scanner) Scanner struct
Keywords
G: (::Vector{Real},=[0, 0, 0],[T/m]) gradient amplitudes for x, y, zΔf: (::Real,=0,[Hz]) RF pulse carrier frequency displacementa: (::Real,=0.46) height appodization window parameterTBP: (::Real,=4) width appodization window parameter
Returns
seq: (::Sequence) Sequence struct with a RF pulse
Examples
julia> sys = Scanner(); durRF = π / 2 / (2π * γ * sys.B1);
julia> seq = PulseDesigner.RF_sinc(sys.B1, durRF, sys);
julia> plot_seq(seq)KomaMRIBase.PulseDesigner.EPI — Functionseq = EPI(FOV::Real, N::Integer, sys::Scanner)Returns a sequence with EPI gradients.
Arguments
FOV: (::Real,[m]) field of viewN: (::Integer) number of pixels in the x and y axissys: (::Scanner) Scanner struct
Returns
seq: (::Sequence) Sequence struct with EPI gradients
Examples
julia> sys, FOV, N = Scanner(), 23e-2, 101
julia> seq = PulseDesigner.EPI(FOV, N, sys)
julia> plot_seq(seq)
julia> plot_kspace(seq)KomaMRIBase.PulseDesigner.radial_base — Functionseq = radial_base(FOV::Real, Nr::Integer, sys::Scanner)Returns a sequence with radial gradients for a single trajectory.
Arguments
FOV: (::Real,[m]) field of viewN: (::Integer) number of pixels along the diametersys: (::Scanner) Scanner struct
Returns
seq: (::Sequence) Sequence struct of a single radial trajectory
Examples
julia> sys, FOV, N = Scanner(), 23e-2, 101
julia> seq = PulseDesigner.radial_base(FOV, N, sys)
julia> plot_seq(seq)
julia> plot_kspace(seq)KomaMRIBase.PulseDesigner.spiral_base — Functionspiral = spiral_base(FOV, N, sys; S0=sys.Smax*2/3, Nint=8, λ=Nint/FOV, BW=60e3)Definition of a spiral base sequence.
References
- Glover, G.H. (1999), Simple analytic spiral K-space algorithm. Magn. Reson. Med.,
42: 412-415. https://doi.org/10.1002/(SICI)1522-2594(199908)42:2<412::AID-MRM25>3.0.CO;2-U
Arguments
FOV: (::Real,[m]) field of viewN: (::Integer) number of pixels along the radioussys: (::Scanner) Scanner struct
Keywords
S0: (::Vector{Real},=sys.Smax*2/3,[T/m/s]) slew rate referenceNint: (::Integer,=8) number of interleavesλ: (::Real,=Nint/FOV,[1/m]) kspace spiral parameterBW: (::Real,=60e3,[Hz]) adquisition parameter
Returns
spiral: (::Function) function that returns aSequencestruct when evaluated
Examples
julia> sys, FOV, N = Scanner(), 23e-2, 101
julia> spiral = PulseDesigner.spiral_base(FOV, N, sys)
julia> seq = spiral(0)
julia> plot_seq(seq)KomaMRIBase.PulseDesigner.EPI_example — Functionseq = EPI_example(; sys=Scanner())Returns a sequence suitable for acquiring the 2D brain example in the provided examples.
Keywords
sys: (::Scanner) Scanner struct
Returns
seq: (::Sequence) EPI example Sequence struct
Examples
julia> seq = PulseDesigner.EPI_example();
julia> plot_seq(seq)