KomaMRIBase.ADCType
adc = 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 samples
  • T: (::Float64, [s]) duration to acquire the samples
  • delay: (::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.DelayType
delay = 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)
KomaMRIBase.DiscreteSequenceType
seqd = 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 vector
  • Gy: (::AbstractVector{T<:Real}, [T/m]) y-gradient vector
  • Gz: (::AbstractVector{T<:Real}, [T/m]) z-gradient vector
  • B1: (::AbstractVector{Complex{T<:Real}}, [T]) RF amplitude vector
  • Δf: (::AbstractVector{T<:Real}, [Hz]) RF carrier frequency displacement vector
  • ADC: (::AbstractVector{Bool}) ADC sample vector
  • t: (::AbstractVector{T<:Real}, [s]) time vector
  • Δt: (::AbstractVector{T<:Real}, [s]) delta time vector

Returns

  • seqd: (::DiscreteSequence) DiscreteSequence struct
KomaMRIBase.GradType
gr = 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 waveform
  • T: (::Real, [s]) duration of the gradient waveform
  • N: (::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)
KomaMRIBase.GradType
gr = Grad(A, T)
gr = Grad(A, T, rise)
gr = Grad(A, T, rise, delay)
gr = Grad(A, T, rise, fall, delay)

The Grad struct represents a gradient of a sequence event.

Arguments

  • A: (::Real or ::Vector, [T/m]) amplitude of the gradient
  • T: (::Real or ::Vector, [s]) duration of the flat-top
  • rise: (::Real, [s]) duration of the rise
  • fall: (::Real, [s]) duration of the fall
  • delay: (::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.PhantomType
obj = Phantom(name, x, y, z, ρ, T1, T2, T2s, Δw, Dλ1, Dλ2, Dθ, ux, uy, uz)

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 name
  • x: (::AbstractVector{T<:Real}, [m]) spin x-position vector
  • y: (::AbstractVector{T<:Real}, [m]) spin y-position vector
  • z: (::AbstractVector{T<:Real}, [m]) spin z-position vector
  • ρ: (::AbstractVector{T<:Real}) spin proton density vector
  • T1: (::AbstractVector{T<:Real}, [s]) spin T1 parameter vector
  • T2: (::AbstractVector{T<:Real}, [s]) spin T2 parameter vector
  • T2s: (::AbstractVector{T<:Real}, [s]) spin T2s parameter vector
  • Δw: (::AbstractVector{T<:Real}, [rad/s]) spin off-resonance parameter vector
  • Dλ1: (::AbstractVector{T<:Real}) spin Dλ1 (diffusion) parameter vector
  • Dλ2: (::AbstractVector{T<:Real}) spin Dλ2 (diffusion) parameter vector
  • : (::AbstractVector{T<:Real}) spin Dθ (diffusion) parameter vector
  • ux: (::Function) displacement field in the x-axis
  • uy: (::Function) displacement field in the y-axis
  • uz: (::Function) displacement field in the z-axis

Returns

  • obj: (::Phantom) Phantom struct

Examples

julia> obj = Phantom(x=[0.0])

julia> obj.ρ
KomaMRIBase.RFType
rf = RF_fun(f::Function, T::Real, N::Int64)

Generate an RF sequence with amplitudes sampled from a function waveform.

Note

This function is not being used in this KomaMRI version.

Arguments

  • f: (::Function, [T]) function for the RF amplitud waveform
  • T: (::Real, [s]) duration of the RF pulse
  • N: (::Int64) number of samples of the RF pulse

Returns

  • rf:(::RF) RF struct with amplitud defined by the function f
KomaMRIBase.RFType
rf = 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: (::Real or ::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.ScannerType
sys = 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 strength
  • B1: (::Real, =10e-6, [T]) maximum RF amplitude
  • Gmax: (::Real, =60e-3, [T/m]) maximum gradient amplitude
  • Smax: (::Real, =500, [mT/m/ms]) gradient's maximum slew-rate
  • ADC_Δt: (::Real, =2e-6, [s]) ADC raster time
  • seq_Δt: (::Real, =1e-5, [s]) sequence-block raster time
  • GR_Δt: (::Real, =1e-5, [s]) gradient raster time
  • RF_Δt: (::Real, =1e-6, [s]) RF raster time
  • RF_ring_down_T: (::Real, =20e-6, [s]) RF ring down time
  • RF_dead_time_T: (::Real, =100e-6, [s]) RF dead time
  • ADC_dead_time_T: (::Real, =10e-6, [s]) ADC dead time

Returns

  • sys: (::Scanner) Scanner struct

Examples

julia> sys = Scanner()

julia> sys.B0
KomaMRIBase.SequenceType
seq = 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 blocks
  • RF: (::Matrix{RF}) RF matrix. The 1 row is for the coil and columns are for blocks
  • ADC: (::Array{ADC,1}) ADC block vector
  • DUR: (::Vector, [s]) duration block vector
  • DEF: (::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
Base.:*Method

Scalar multiplication of a phantom

Base.:+Method

Addition of phantoms

Base.:+Method
seq = +(s::Sequence, d::Delay)
seq = +(d::Delay, s::Sequence)

Add a delay to sequence struct. It ultimately affects to the duration of the gradients of a sequence.

Arguments

  • s: (::Sequence) sequence struct
  • d: (::Delay) delay struct

Returns

  • seq: (::Sequence) delayed sequence
Base.getpropertyMethod
y = getproperty(x::Vector{ADC}, f::Symbol)

Overloads Base.getproperty(). It is meant to access properties of the ADC vector x directly without the need to iterate elementwise.

Arguments

  • x: (::Vector{ADC}) vector of ADC structs
  • f: (::Symbol, opts: [:N, :T, :delay, :Δf, , :dur]) input symbol that represents a property of the ADC structs

Returns

  • y: (::Vector{Any}) vector with the property defined by the f for all elements of the ADC vector x
Base.getpropertyMethod
y = getproperty(x::Vector{Grad}, f::Symbol)
y = getproperty(x::Matrix{Grad}, f::Symbol)

Overloads Base.getproperty(). It is meant to access properties of the Grad vector x directly without the need to iterate elementwise.

Arguments

  • x: (::Vector{Grad} or ::Matrix{Grad}) vector or matrix of Grad structs
  • f: (::Symbol, opts: [:x, :y, :z, :T, :delay, :rise, :delay, :dur, :A, f]) input symbol that represents a property of the vector or matrix of Grad structs

Returns

  • y: (::Vector{Any} or ::Matrix{Any}) vector or matrix with the property defined by the symbol f for all elements of the Grad vector or matrix x
Base.getpropertyMethod
y = getproperty(x::Vector{RF}, f::Symbol)
y = getproperty(x::Matrix{RF}, f::Symbol)

Overloads Base.getproperty(). It is meant to access properties of the RF vector x directly without the need to iterate elementwise.

Arguments

  • x: (::Vector{RF} or ::Matrix{RF}) vector or matrix of RF structs
  • f: (::Symbol, opts: [:A, :Bx, :By, :T, :Δf, :delay and :dur]) input symbol that represents a property of the vector or matrix of RF structs

Returns

  • y: (::Vector{Any} or ::Matrix{Any}) vector with the property defined by the symbol f for all elements of the RF vector or matrix x
Base.showMethod
str = show(io::IO, s::Delay)

Displays the delay time in m[s] of the delay struct s in the julia REPL.

Arguments

  • s: (::Delay) delay struct

Returns

  • str: (::String) output string message
Base.showMethod
str = show(io::IO, x::Grad)

Displays information about the Grad struct x in the julia REPL.

Arguments

  • x: (::Grad) Grad struct

Returns

  • str (::String) output string message
Base.showMethod
str = show(io::IO, x::RF)

Displays information about the RF struct x in the julia REPL.

Arguments

  • x: (::RF) RF struct

Returns

  • str: (::String) output string message
Base.showMethod
str = show(io::IO, s::Sequence)

Displays information about the Sequence struct s in the julia REPL.

Arguments

  • s: (::Sequence) Sequence struct

Returns

  • str (::String) output string message
Base.sizeMethod

Size and length of a phantom

Base.viewMethod

Separate object spins in a sub-group (lightweigth).

KomaMRIBase.brain_phantom2DMethod
obj = brain_phantom2D(; axis="axial", ss=4)

Creates a two-dimensional brain Phantom struct.

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 phantom
  • ss: (::Integer, =4) subsampling parameter in all axis

Returns

  • obj: (::Phantom) Phantom struct

Examples

julia> obj = brain_phantom2D(; axis="sagittal", ss=1)

julia> plot_phantom_map(obj, :ρ)
KomaMRIBase.brain_phantom3DMethod
obj = brain_phantom3D(; ss=4)

Creates a three-dimentional brain Phantom struct.

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, =4) subsampling parameter in all axes

Returns

  • obj: (::Phantom) 3D Phantom struct

Examples

julia> obj = brain_phantom3D(; ss=5)

julia> plot_phantom_map(obj, :ρ)
KomaMRIBase.cumtrapzMethod
y = 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 array
  • x: (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.discretizeMethod
seqd = 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.durMethod
y = dur(x::Grad)
y = dur(x::Vector{Grad})

Duration time in [s] of Grad struct or Grad array. When the input is a gradient vector, then the duration is the maximum duration of all the elements of the gradient vector.

Arguments

  • x: (::Grad or ::Vector{Grad}) RF struct or RF array

Returns

  • y: (::Float64, [s]) duration of the RF struct or RF array
KomaMRIBase.durMethod
y = dur(x::RF)
y = dur(x::Array{RF,1})
y = dur(x::Array{RF,2})

Duration time in [s] of RF struct or RF array.

Arguments

  • x: (::RF or ::Array{RF,1} or ::Array{RF,2}) RF struct or RF array

Returns

  • y: (::Float64, [s]) duration of the RF struct or RF array
KomaMRIBase.durMethod
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.dursMethod
ΔT = durs(x::Sequence)

Returns the array of durations of sequence's blocks in [s].

Arguments

  • x: (::Sequence) Sequence struct

Returns

  • ΔT: (::Vector{Real}, [s]) array of durations of sequence's blocks
KomaMRIBase.get_M1Method
M1, M1_adc = get_M1(seq::Sequence; Δt=1)

Outputs the designed M1 of the Sequence seq.

Arguments

  • seq: (::Sequence) Sequence struct
  • Δt: (::Real, =1, [s]) nominal delta time separation between two time samples for ADC acquisition and Gradients

Returns

  • M1: (3-column ::Matrix{Float64}) First moment
  • M1_adc: (3-column ::Matrix{Float64}) First moment sampled at ADC points
KomaMRIBase.get_M2Method
M2, M2_adc = get_M2(seq::Sequence; Δt=1)

Outputs the designed M2 of the Sequence seq.

Arguments

  • seq: (::Sequence) Sequence struct
  • Δt: (::Real, =1, [s]) nominal delta time separation between two time samples for ADC acquisition and Gradients

Returns

  • M2: (3-column ::Matrix{Float64}) Second moment
  • M2_adc: (3-column ::Matrix{Float64}) Second moment sampled at ADC points
KomaMRIBase.get_RF_centerMethod
t = get_RF_center(x::RF)

Calculates the time where is the center of the RF pulse x. This calculation includes the RF delay.

Arguments

  • x: (::RF) RF struct

Returns

  • t: (::Int64, [s]) time where is the center of the RF pulse x
KomaMRIBase.get_RF_typesMethod
rf_idx, rf_type = get_RF_types(seq, t)

Get RF centers and types (excitation or precession). Useful for k-space calculations.

Arguments

  • seq: (::Sequence) Sequence struct
  • t: (::Vector{Float64}, [s]) time values

Returns

  • rf_idx: (::Vector{Int64}) indices of the RF centers
  • rf_type: (::Vector{Int64}, opts: [0, 1]) RF type (0: excitation, 1: precession)
KomaMRIBase.get_adc_phase_compensationMethod
phase = get_adc_phase_compensation(seq)

Returns the array of phases for every acquired sample in the sequence seq.

Note

This function is useful to compensate the phase when the RF pulse has a phase too. Refer to the end of the run_sim_time_iter function to see its usage.

Arguments

  • seq: (::Sequence) sequence struct

Returns

  • phase: (::Vector{Complex{Int64}}, [rad]) array of phases for every acquired sample
KomaMRIBase.get_adc_sampling_timesMethod
times = 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_block_start_timesMethod
T0 = 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_eddy_currentsMethod
EC, EC_adc = get_eddy_currents(seq::Sequence; Δt=1, λ=80e-3)

Outputs the designed eddy currents of the Sequence seq.

Arguments

  • seq: (::Sequence) Sequence struct
  • Δt: (::Real, =1, [s]) nominal delta time separation between two time samples for ADC acquisition and Gradients
  • λ: (::Float64, =80e-3, [s]) eddy current decay constant time

Returns

  • EC: (3-column ::Matrix{Float64}) Eddy currents
  • EC_adc: (3-column ::Matrix{Float64}) Eddy currents sampled at ADC points
KomaMRIBase.get_flip_angleMethod
α = 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 struct x
KomaMRIBase.get_flip_anglesMethod
y = 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
KomaMRIBase.get_gradsMethod
Gx, Gy, Gz = get_grads(seq, t::Vector)
Gx, Gy, Gz = get_grads(seq, t::Matrix)

Get the gradient array from sequence seq evaluated in time points t.

Arguments

  • seq: (::Sequence) Sequence struct
  • t: (::Vector{Float64} or 1-row ::Matrix{Float64}, [s]) times to evaluate

Returns

  • Gx: (Vector{Float64} or 1-row ::Matrix{Float64}, [T]) gradient vector values in the x direction
  • Gy: (Vector{Float64} or 1-row ::Matrix{Float64}, [T]) gradient vector values in the y direction
  • Gz: (Vector{Float64} or 1-row ::Matrix{Float64}, [T]) gradient vector values in the z direction
KomaMRIBase.get_kspaceMethod
kspace, kspace_adc = get_kspace(seq::Sequence; Δt=1)

Outputs the designed k-space trajectory of the Sequence seq.

Arguments

  • seq: (::Sequence) Sequence struct
  • Δt: (::Real, =1, [s]) nominal delta time separation between two time samples for ADC acquisition and Gradients

Returns

  • kspace: (3-column ::Matrix{Float64}) kspace
  • kspace_adc: (3-column ::Matrix{Float64}) adc kspace
KomaMRIBase.get_rfsMethod
B1, Δf_rf  = get_rfs(seq::Sequence, t)

Returns the RF pulses and the delta frequency.

Arguments

  • seq: (::Sequence) Sequence struct
  • t: (1-row ::Matrix{Float64}, [s]) time points

Returns

  • B1: (1-row ::Matrix{ComplexF64}, [T]) vector of RF pulses
  • Δf_rf: (1-row ::Matrix{Float64}, [Hz]) delta frequency vector
KomaMRIBase.get_samplesMethod
samples = 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 to Inf
  • max_rf_samples: (::Integer, =Inf) maximum number of samples for the RF struct

Returns

  • samples: (::NamedTuple) contains samples for gx, gy, gz, rf, and adc events. Each event, represented by e::NamedTuple, includes time samples (e.t) and amplitude samples (e.A)
KomaMRIBase.get_slew_rateMethod
SR, SR_adc = get_slew_rate(seq::Sequence; Δt=1)

Outputs the designed slew rate of the Sequence seq.

Arguments

  • seq: (::Sequence) Sequence struct
  • Δt: (::Real, =1, [s]) nominal delta time separation between two time samples for ADC acquisition and Gradients

Returns

  • SR: (3-column ::Matrix{Float64}) Slew rate
  • SR_adc: (3-column ::Matrix{Float64}) Slew rate sampled at ADC points
KomaMRIBase.get_theo_AMethod
A = get_theo_A(g::Grad; off_val=0)
A = get_theo_A(r::RF; off_val=0, max_rf_samples=Inf)
A = get_theo_A(d::ADC; off_val=0)

Get the theoretical amplitudes of a rectangle waveform for Grad, RF or ADC structs. This are 5 points: delay, start, rise, stop and fall.

Note

In some cases the array result can have duplicated points, so it is necessary to remove them whenever necessary.

Arguments

  • g: (::Grad) Gradient struct
  • r: (::RF) RF struct
  • d: (::ADC) ADC truct

Keywords

  • off_val: (::Float64, =0) offset value for amplitude. In general, it is used for not showing some points in plots by giving an Inf value
  • max_rf_samples: (::Float64, =Inf) number of maximum samples for the RF struct. In general, this parameter is not necessary to set

Returns

  • A: (::Vector{Float64}) vector with the amplitude key points of the rectangle waveform
KomaMRIBase.get_theo_GiMethod
t, g = get_theo_Gi(seq, idx)

Get the theoretical gradient for a sequence in a defined axis.

Arguments

  • seq: (::Sequence) Sequence struct
  • idx: (::Int64, opts=[1, 2, 3]) axis x, y or z for the gradient

Returns

  • t: (::Vector{Float64}) time key points
  • g: (::Vector{Float64}) amplitude key points
KomaMRIBase.get_theo_tMethod
t = get_theo_t(g::Grad)
t = get_theo_t(r::RF; max_rf_samples=Inf)
t = get_theo_t(d::ADC)

Get the theoretical times of a rectangle waveform for Grad, RF or ADC structs. This are 5 points: delay, start, rise, stop and fall.

Note

In some cases the array result can have duplicated points, so it is necessary to remove them whenever necessary.

Arguments

  • g: (::Grad) Gradient struct
  • r: (::RF) RF struct
  • d: (::ADC) ADC truct

Keywords

  • max_rf_samples: (::Float64, =Inf) number of maximum samples for the RF struct. In general, this parameter is not necessary to set

Returns

  • t: (::Vector{Float64}) vector with the time key points of the rectangle waveform
KomaMRIBase.get_variable_timesMethod
t, Δt = get_variable_times(seq; Δt=1e-3, Δt_rf=1e-5)

This function returns non-uniform time points that are relevant in the sequence seq.

Arguments

  • seq: (::Sequence) Sequence struct

Keywords

  • Δt: (::Real, =1e-3, [s]) nominal delta time separation between two time samples for ADC acquisition and Gradients (by nominal we mean that the time separation should be at most Δt when the samples are regarded by KomaMRI.is_ADC_on or KomaMRI.is_GR_on), otherwise the time points are not necessary and the separation will be bigger)
  • Δt_rf: (::Real, =1e-5, [s]) nominal delta time separation between two time samples for RF excitation (by nominal we mean that the time separation should be at most Δt_rf when the samples are regarded by KomaMRI.is_RF_on, otherwise the time points are not necessary and the separation will be bigger)

Returns

  • t: (::Vector{Float64}, [s]) time array with non-uniform time values
  • Δt: (::Vector{Float64}, [s]) delta time array with the separation between two adjacent time points of the t time array
KomaMRIBase.is_ADC_onMethod
y = 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 struct
  • t: (::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
KomaMRIBase.is_DelayMethod
y = is_Delay(x::Sequence)

Tells if the sequence seq is a delay.

Arguments

  • x: (::Sequence) Sequence struct

Returns

  • y::Bool: boolean that tells whether or not the sequence is a delay
KomaMRIBase.is_GR_onMethod
y = 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_onMethod
y = 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_onMethod
y = 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_onMethod
y = 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_RF_onMethod
y = 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 struct
  • t: (::Vector{Float64}, [s]) time to check

Returns

  • y: (::Bool) boolean that tells whether or not the RF in the sequence is active
KomaMRIBase.kfoldpermMethod
array_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 ordered
  • k: (::Integer) number of groups to divide the N elements.

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 is k groups, but this could increase by adding elements to the breaks input array
KomaMRIBase.pelvis_phantom2DMethod
obj = pelvis_phantom2D(; ss=4)

Creates a two-dimensional pelvis Phantom struct.

Keywords

  • ss: (::Integer, =4) subsampling parameter

Returns

  • obj: (::Phantom) Phantom struct

Examples

julia> obj = pelvis_phantom2D(; ss=1)

julia> pelvis_phantom2D(obj, :ρ)
KomaMRIBase.points_from_key_timesMethod
t = points_from_key_times(times; dt)

Returns a vector which contains the same points as times but with additional points that have a separation of at most dt.

Arguments

  • times: (::Vector{Float64}, [s]) time array with key points you want to keep

Keywords

  • dt: (::Float64, [s]) maximum delta time separation between two time samples

Returns

  • t: (::Vector{Float64}, [s]) time array with the same points as the input array but with additional points that have a separation of at most dt.
KomaMRIBase.rotxMethod
Rx = rotx(θ::Real)

Rotates vector counter-clockwise with respect to the x-axis.

Arguments

  • θ: (::Real, [rad]) rotation angle

Returns

  • Rx: (::Matrix{Int64}) rotation matrix
KomaMRIBase.rotyMethod
Ry = roty(θ::Real)

Rotates vector counter-clockwise with respect to the y-axis.

Arguments

  • θ: (::Real, [rad]) rotation angle

Returns

  • Ry: (::Matrix{Int64}) rotation matrix
KomaMRIBase.rotzMethod
Rz = rotz(θ::Real)

Rotates vector counter-clockwise with respect to the z-axis.

Arguments

  • θ: (::Real, [rad]) rotation angle

Returns

  • Rz: (::Matrix{Int64}) rotation matrix
KomaMRIBase.trapzMethod
y = trapz(Δt, x)

Trapezoidal integration for every spin of a phantom.

Note

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 array
  • x: (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.⏢Method
y = ⏢(A, t, ΔT, ζ1, ζ2, delay)

Generates a trapezoidal waveform vector.

Arguments

  • A: (::Real) amplitude
  • t: (::Vector{Float64}, [s]) times to evaluate (actually it's a 1-row ::Matrix{Float64})
  • ΔT: (::Real, [s]) time duration of the top-flat
  • ζ1: (::Real, [s]) rise time duration
  • ζ2: (::Real, [s]) fall time duration
  • delay: (::Real, [s]) delay time

Returns

  • y: (::Vector{Float64}) trapezoidal waveform (actually it's a 1-row ::Matrix{Float64})
KomaMRIBase.PulseDesigner.EPIMethod
seq = EPI(FOV::Real, N::Integer, sys::Scanner)

Returns a sequence with EPI gradients.

Arguments

  • FOV: (::Real, [m]) field of view
  • N: (::Integer) number of pixels in the x and y axis
  • sys: (::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.EPI_exampleMethod
seq = 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)
KomaMRIBase.PulseDesigner.RF_hardMethod
seq = 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 amplitude
  • T: (::Real, [s]) RF pulse duration
  • sys: (::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_sincMethod
seq = 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 amplitude
  • T: (::Real, [s]) RF sinc duration
  • sys: (::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
  • a: (::Real, =0.46) height appodization window parameter
  • TBP: (::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.radial_baseMethod
seq = radial_base(FOV::Real, Nr::Integer, sys::Scanner)

Returns a sequence with radial gradients for a single trajectory.

Arguments

  • FOV: (::Real, [m]) field of view
  • N: (::Integer) number of pixels along the diameter
  • sys: (::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_baseMethod
spiral = 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 view
  • N: (::Integer) number of pixels along the radious
  • sys: (::Scanner) Scanner struct

Keywords

  • S0: (::Vector{Real}, =sys.Smax*2/3, [T/m/s]) slew rate reference
  • Nint: (::Integer, =8) number of interleaves
  • λ: (::Real, =Nint/FOV, [1/m]) kspace spiral parameter
  • BW: (::Real, =60e3, [Hz]) adquisition parameter

Returns

  • spiral: (::Function) function that returns a Sequence struct 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)