KomaMRIBase

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
source
KomaMRIBase.PhantomType
obj = 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 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
  • motion: (::AbstractMotion{T<:Real}) motion

Returns

  • obj: (::Phantom) Phantom struct

Examples

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

julia> obj.ρ
source
KomaMRIBase.brain_phantom2DFunction
phantom = 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/tissuemrparameters.txt

Keywords

  • axis: (::String, ="axial", opts=["axial", "coronal", "sagittal"]) orientation of the phantom
  • 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], 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, :ρ)
source
KomaMRIBase.brain_phantom3DFunction
obj = brain_phantom3D(; ss=4, us=1, start_end=[160,200])

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/tissuemrparameters.txt

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, :ρ)
source
KomaMRIBase.pelvis_phantom2DFunction
obj = 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, :ρ)
source
KomaMRIBase.heart_phantomFunction
obj = heart_phantom(
    circumferential_strain, radial_strain, rotation_angle; 
    heart_rate, asymmetry
)

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.

Keywords

  • circumferential_strain: (::Real, =-0.3) contraction parameter. Between -1 and 1
  • radial_strain: (::Real, =-0.3) contraction parameter. Between -1 and 1
  • rotation_angle: (::Real, =15.0, [º]) maximum rotation angle
  • heart_rate: (::Real, =60, [bpm]) heartbeat frequency
  • temporal_asymmetry: (::Real, =0.2) time fraction of the period in which the systole occurs. Therefore, diastole lasts for period * (1 - temporal_asymmetry)

Returns

  • obj: (::Phantom) Heart-like LV phantom struct
source
KomaMRIBase.NoMotionType
nomotion = NoMotion{T<:Real}()

NoMotion struct. It is used to create static phantoms.

Returns

  • nomotion: (::NoMotion) NoMotion struct

Examples

julia> nomotion = NoMotion{Float64}()
source
KomaMRIBase.MotionListType
motionlist = MotionList(motions...)

MotionList struct. The other option, instead of NoMotion, is to define a dynamic phantom by means of the MotionList struct. It is composed by one or more Motion instances.

Arguments

  • motions: (::Vector{Motion{T<:Real}}) vector of Motion instances

Returns

  • motionlist: (::MotionList) MotionList struct

Examples

julia>  motionlist = MotionList(
            Motion(
                action = Translate(0.01, 0.0, 0.02),
                time = TimeRange(0.0, 1.0),
                spins = AllSpins()
            ),
            Motion(
                action = Rotate(0.0, 0.0, 45.0),
                time = Periodic(1.0),
                spins = SpinRange(1:10)
            )
        )
source
KomaMRIBase.get_spin_coordsFunction
x, y, z = get_spin_coords(motionset, 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_{ {spins}}$ rows and length(t) columns.

Arguments

  • motionset: (::AbstractMotion{T<:Real}) phantom motion
  • 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
  • t: horizontal array of time instants

Returns

  • x, y, z: (::Tuple{AbstractArray, AbstractArray, AbstractArray}) spin positions over time
source

Motion

KomaMRIBase.MotionType
motion = Motion(action)
motion = Motion(action, time)
motion = Motion(action, time, spins)

Motion struct. It defines the motion, during a certain time interval, of a given group of spins. It is composed by three fields: action, which defines the motion itself, time, which accounts for the time during which the motion takes place, and spins, which indicates the spins that are affected by that motion.

Arguments

  • action: (::AbstractAction{T<:Real}) action, such as Translate or Rotate
  • time: (::AbstractTimeSpan{T<:Real}, =TimeRange(0.0)) time information about the motion
  • spins: (::AbstractSpinSpan, =AllSpins()) spin indexes affected by the motion

Returns

  • motion: (::Motion) Motion struct

Examples

julia> motion =  Motion(
            action = Translate(0.01, 0.0, 0.02),
            time = TimeRange(0.0, 1.0),
            spins = SpinRange(1:10)
       )
source

AbstractAction types

KomaMRIBase.TranslateType
translate = Translate(dx, dy, dz)

Translate struct. It produces a linear translation. Its fields are the final displacements in the three axes (dx, dy, dz).

Arguments

  • dx: (::Real, [m]) translation in x
  • dy: (::Real, [m]) translation in y
  • dz: (::Real, [m]) translation in z

Returns

  • translate: (::Translate) Translate struct

Examples

julia> translate = Translate(dx=0.01, dy=0.02, dz=0.03)
source
KomaMRIBase.TranslateMethod
translate = Translate(dx, dy, dz, time, spins)

Arguments

  • dx: (::Real, [m]) translation in x
  • dy: (::Real, [m]) translation in y
  • dz: (::Real, [m]) translation in z
  • time: (::AbstractTimeSpan{T<:Real}) time information about the motion
  • spins: (::AbstractSpinSpan) spin indexes affected by the motion

Returns

  • translate: (::Motion) Motion struct

Examples

julia> translate = Translate(0.01, 0.02, 0.03, TimeRange(0.0, 1.0), SpinRange(1:10))
source
KomaMRIBase.RotateType
rotate = Rotate(pitch, roll, yaw)

Rotate struct. It produces a rotation 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):

Head Rotation Axis

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 x
  • roll: (::Real, [º]) rotation in y
  • yaw: (::Real, [º]) rotation in z

Returns

  • rotate: (::Rotate) Rotate struct

Examples

julia> rotate = Rotate(pitch=15.0, roll=0.0, yaw=20.0)
source
KomaMRIBase.RotateMethod
rotate = Rotate(pitch, roll, yaw, spins)

Arguments

  • pitch: (::Real, [º]) rotation in x
  • roll: (::Real, [º]) rotation in y
  • yaw: (::Real, [º]) rotation in z
  • time: (::AbstractTimeSpan{T<:Real}) time information about the motion
  • spins: (::AbstractSpinSpan) spin indexes affected by the motion

Returns

  • rotate: (::Motion) Motion struct with Rotate action

Examples

julia> rotate = Rotate(15.0, 0.0, 20.0, TimeRange(0.0, 1.0), SpinRange(1:10))
source
KomaMRIBase.HeartBeatType
heartbeat = HeartBeat(circumferential_strain, radial_strain, longitudinal_strain)

HeartBeat struct. It produces a heartbeat-like motion, characterised by three types of strain: circumferential, radial and longitudinal

Arguments

  • circumferential_strain: (::Real) contraction parameter
  • radial_strain: (::Real) contraction parameter
  • longitudinal_strain: (::Real) contraction parameter

Returns

  • heartbeat: (::HeartBeat) HeartBeat struct

Examples

julia> heartbeat = HeartBeat(circumferential_strain=-0.3, radial_strain=-0.2, longitudinal_strain=0.0)
source
KomaMRIBase.HeartBeatMethod
heartbeat = HeartBeat(circumferential_strain, radial_strain, longitudinal_strainl, time, spins)

Arguments

  • circumferential_strain: (::Real) contraction parameter
  • radial_strain: (::Real) contraction parameter
  • longitudinal_strain: (::Real) contraction parameter
  • time: (::AbstractTimeSpan{T<:Real}) time information about the motion
  • spins: (::AbstractSpinSpan) spin indexes affected by the motion

Returns

  • heartbeat: (::Motion) Motion struct with HeartBeat action

Examples

julia> heartbeat = HeartBeat(-0.3, -0.2, 0.0, TimeRange(0.0, 1.0), SpinRange(1:10))
source
KomaMRIBase.PathType
path = Path(dx, dy, dz)

Path struct. For this action (and for FlowPath), motion is not defined solely on the basis of three numerical parameters, one for each spatial direction, as occurs for the Translate, Rotate and HeartBeat actions.

For this action, it is necessary to define motion for each spin independently, in x (dx), y (dy) and z (dz). dx, dy and dz are now three matrixes, of ($N_{spins}* \times \; N_{discrete\,times}$) each. This means that each row corresponds to a spin trajectory over a set of discrete time instants.

Note

*When creating a motion with Flow or FlowPath, you must make sure that the number of rows of the matrices dx, dy and dz matches the number of spins that are affected by the motion.

Remember that the range of spins affected by a motion is defined by the spins field of the Motion struct

example:

julia> motion = Motion(
    action = Path(
        dx=[0.01 0.02; 0.02 0.03],  # 2 rows
        dy=[0.02 0.03; 0.03 0.04], 
        dz=[0.03 0.04; 0.04 0.05]),
    time = TimeRange(0.0, 1.0),
    spins = SpinRange(1:2)          # 2 spins
)

Arguments

  • dx: (::AbstractArray{T<:Real}, [m]) displacements in x
  • dy: (::AbstractArray{T<:Real}, [m]) displacements in y
  • dz: (::AbstractArray{T<:Real}, [m]) displacements in z

Returns

  • path: (::Path) Path struct

Examples

julia> path = Path(
           dx=[0.01 0.02; 0.02 0.03], 
           dy=[0.02 0.03; 0.03 0.04], 
           dz=[0.03 0.04; 0.04 0.05]
       )
source
KomaMRIBase.PathMethod
path = Path(dx, dy, dz, time, spins)

Arguments

  • dx: (::AbstractArray{T<:Real}, [m]) displacements in x
  • dy: (::AbstractArray{T<:Real}, [m]) displacements in y
  • dz: (::AbstractArray{T<:Real}, [m]) displacements in z
  • time: (::AbstractTimeSpan{T<:Real}) time information about the motion
  • spins: (::AbstractSpinSpan) spin indexes affected by the motion

Returns

  • path: (::Motion) Motion struct with Path action

Examples

julia> path = Path(
          [0.01 0.02; 0.02 0.03], 
          [0.02 0.03; 0.03 0.04], 
          [0.03 0.04; 0.04 0.05], 
          TimeRange(0.0, 1.0), 
          SpinRange(1:10)
       )
source
KomaMRIBase.FlowPathType
flowpath = FlowPath(dx, dy, dz, spin_reset)

FlowPath struct. This action is the same as Path, except that it includes an additional field, called spin_reset, which accounts for spins leaving the volume and being remapped to another input position. When this happens, the magnetization state of these spins must be reset during the simulation.

As with the dx, dy and dz matrices, spin_reset has a size of ($N_{spins} \times \; N_{discrete\,times}$).

Arguments

  • dx: (::AbstractArray{T<:Real}, [m]) displacements in x
  • dy: (::AbstractArray{T<:Real}, [m]) displacements in y
  • dz: (::AbstractArray{T<:Real}, [m]) displacements in z
  • spin_reset: (::AbstractArray{Bool}) reset spin state flags

Returns

  • flowpath: (::FlowPath) FlowPath struct

Examples

julia> flowpath = FlowPath(
           dx=[0.01 0.02; 0.02 0.03], 
           dy=[0.02 0.03; 0.03 0.04], 
           dz=[0.03 0.04; 0.04 -0.04],
           spin_reset=[false false; false true]
       )
source
KomaMRIBase.FlowPathMethod
flowpath = FlowPath(dx, dy, dz, spin_reset, time, spins)

Arguments

  • dx: (::AbstractArray{T<:Real}, [m]) displacements in x
  • dy: (::AbstractArray{T<:Real}, [m]) displacements in y
  • dz: (::AbstractArray{T<:Real}, [m]) displacements in z
  • spin_reset: (::AbstractArray{Bool}) reset spin state flags
  • time: (::AbstractTimeSpan{T<:Real}) time information about the motion
  • spins: (::AbstractSpinSpan) spin indexes affected by the motion

Returns

  • flowpath: (::Motion) Motion struct with FlowPath action

Examples

julia> flowpath = FlowPath(
          [0.01 0.02; 0.02 0.03], 
          [0.02 0.03; 0.03 0.04], 
          [0.03 0.04; 0.04 0.05], 
          [false false; false true],
          TimeRange(0.0, 1.0), 
          SpinRange(1:10)
       )
source
KomaMRIBase.TimeRangeType
timerange = TimeRange(t_start, t_end)

TimeRange struct. It is a specialized type that inherits from AbstractTimeSpan and defines a time interval, with start and end times.

Arguments

  • t_start: (::Real, [s]) start time
  • t_end: (::Real, [s]) end time

Returns

  • timerange: (::TimeRange) TimeRange struct

Examples

julia> timerange = TimeRange(0.0, 1.0)
source
KomaMRIBase.PeriodicType
periodic = Periodic(period, asymmetry)

Periodic struct. It is a specialized type that inherits from AbstractTimeSpan, designed to work with time intervals that repeat periodically. It includes a measure of asymmetry in order to recreate a asymmetric period.

Arguments

  • period: (::Real, [s]) period duration
  • asymmetry: (::Real, =0.5) temporal asymmetry factor. Between 0 and 1.

Returns

  • periodic: (::Periodic) Periodic struct

Examples

julia> periodic = Periodic(1.0, 0.2)
source
KomaMRIBase.unit_timeFunction
t_unit = unit_time(t, time_range)

The unit_time function normalizes a given array of time values t to a unit interval [0, 1] based on a specified start time t_start and end time t_end. This function is used for non-periodic motions, where each element of t is transformed to fit within the range [0, 1] based on the provided start and end times.

Unit Time

Arguments

  • t: (::AbstractArray{T<:Real}, [s]) array of time values to be normalized
  • time_range: (::TimeRange{T<:Real}, [s]) time interval (defined by t_start and t_end) over which we want to normalise

Returns

  • t_unit: (::AbstractArray{T<:Real}, [s]) array of normalized time values

Examples

julia> t_unit = KomaMRIBase.unit_time([0.0, 1.0, 2.0, 3.0, 4.0, 5.0], TimeRange(1.0, 4.0))
6-element Vector{Float64}:
 0.0
 0.0
 0.333
 0.666
 1.0
 1.0
source
t_unit = unit_time(t, periodic)

The unit_time function normalizes a given array of time values t to a unit interval [0, 1] for periodic motions, based on a specified period and an asymmetry factor. This function is useful for creating triangular waveforms or normalizing time values in periodic processes.

Unit Time Triangular

Arguments

  • t: (::AbstractArray{T<:Real}, [s]) array of time values to be normalized
  • periodic: (::Periodic{T<:Real}, [s]) information about the period and the temporal asymmetry

Returns

  • t_unit: (::AbstractArray{T<:Real}, [s]) array of normalized time values

Examples

julia> t_unit = KomaMRIBase.unit_time([0.0, 1.0, 2.0, 3.0, 4.0, 5.0], Periodic(4.0, 0.5))
6-element Vector{Float64}:
 0.0
 0.5
 1.0
 0.5
 0.0
 0.5
source

AbstractSpinSapn types

KomaMRIBase.AllSpinsType
allspins = AllSpins()

AllSpins struct. It is a specialized type that inherits from AbstractSpinSpan and is used to cover all the spins of a phantom.

Returns

  • allspins: (::AllSpins) AllSpins struct

Examples

julia> allspins = AllSpins()
source
KomaMRIBase.SpinRangeType
spinrange = SpinRange(range)

SpinRange struct. It is a specialized type that inherits from AbstractSpinSpan and is used to select a certain range and number of spins.

Arguments

  • range: (::AbstractVector) spin id's. This argument can be a Range, a Vector or a BitVector

Returns

  • spinrange: (::SpinRange) SpinRange struct

Examples

julia> spinrange = SpinRange(1:10)
julia> spinrange = SpinRange([1, 3, 5, 7])
julia> spinrange = SpinRange(obj.x .> 0)
source
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
source
KomaMRIBase.durFunction
y = 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: (::Grad or ::Vector{Grad} or ::Matrix{Grad}) Grad struct or Grad Array

Returns

  • y: (::Float64, [s]) duration of the Grad struct or Grad Array
source
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: (::RF or ::Vector{RF} or ::Matrix{RF}) RF struct or RF array

Returns

  • y: (::Float64, [s]) duration of the RF struct or RF array
source
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
source
KomaMRIBase.get_block_start_timesFunction
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
source
KomaMRIBase.get_flip_anglesFunction
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
source

Grad

KomaMRIBase.GradType
gr = 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: (::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)
source
KomaMRIBase.GradMethod
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)
source

RF

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)
source
KomaMRIBase.RFMethod
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
source
KomaMRIBase.get_flip_angleFunction
α = 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
source

ADC

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)
source
KomaMRIBase.get_adc_sampling_timesFunction
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
source
KomaMRIBase.get_adc_phase_compensationFunction
comp = 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
source

Delay

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)
source

Rotation matrices

KomaMRIBase.rotxFunction
Rx = rotx(θ::Real)

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

Arguments

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

Returns

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

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

Arguments

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

Returns

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

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

Arguments

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

Returns

  • Rz: (::Matrix{Int64}) rotation matrix
source

Moments

KomaMRIBase.get_MkFunction
Mk, 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 struct
  • k: (::Integer) order of the moment to be computed
  • Δt: (::Real, =1, [s]) nominal delta time separation between two time samples for ADC acquisition and Gradients
  • skip_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 moment
  • Mk_adc: (3-column ::Matrix{Real}) $k$th-order moment sampled at ADC times
source

Event checks

KomaMRIBase.is_RF_onFunction
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
source
KomaMRIBase.is_GR_onFunction
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
source
KomaMRIBase.is_Gx_onFunction
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
source
KomaMRIBase.is_Gy_onFunction
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
source
KomaMRIBase.is_Gz_onFunction
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
source
KomaMRIBase.is_ADC_onFunction
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
source

DiscreteSequence

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
source
KomaMRIBase.discretizeFunction
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
source
KomaMRIBase.get_samplesFunction
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)
source
KomaMRIBase.timesFunction
t = times(gr::Grad)
t = times(rf::RF)
t = times(adc::ADC)

Get time samples of MRI sequence event.

Arguments

  • gr: (::Grad) Gradient struct
  • rf: (::RF) RF struct
  • adc: (::ADC) ADC struct

Returns

  • t: (::Vector{Number}) vector with time samples
source
times = times(motion)
source

times

source

times

source
times = times(motion)
source
KomaMRIBase.amplsFunction
A = ampls(g::Grad)
A = ampls(r::RF)
A = ampls(d::ADC)

Get amplitude samples of MRI sequence event.

Arguments

  • gr: (::Grad) Gradient struct
  • rf: (::RF) RF struct
  • adc: (::ADC) ADC struct

Returns

  • A: (::Vector{Number}) vector with amplitude samples
source
KomaMRIBase.freqsFunction
f = freqs(r::RF)

Get frequency samples of MRI sequence event.

Arguments

  • rf: (::RF) RF struct

Returns

  • f: (::Vector{Number}) vector with frequency samples
source

Other functions

KomaMRIBase.trapzFunction
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
source
KomaMRIBase.cumtrapzFunction
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
source
KomaMRIBase.kfoldpermFunction
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
source

Sequence Building Blocks (SBB)

KomaMRIBase.PulseDesigner.RF_hardFunction
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)
source
KomaMRIBase.PulseDesigner.RF_sincFunction
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)
source
KomaMRIBase.PulseDesigner.EPIFunction
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)
source
KomaMRIBase.PulseDesigner.radial_baseFunction
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)
source
KomaMRIBase.PulseDesigner.spiral_baseFunction
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)
source
KomaMRIBase.PulseDesigner.EPI_exampleFunction
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)
source