API

In the following, you find the documentation of all exported functions of the MRFingerprintingRecon.jl package:

MRFingerprintingRecon.FFTNormalOpMethod
FFTNormalOp(img_shape, trj, U; cmaps)
FFTNormalOp(M, U; cmaps)
FFTNormalOp(Λ; cmaps)

Create normal operator of FFT operator. Differentiate between functions exploiting a pre-calculated kernel basis Λ and the functions which calculate Λ based on a passed trajectory trj or mask M.

Arguments

  • img_shape::Tuple{Int}: Image dimensions
  • traj::Vector{Matrix{Float32}}: Trajectory
  • U::Matrix{ComplexF32}=(1,): Basis coefficients of subspace
  • cmaps::Matrix{ComplexF32}: Coil sensitivities
  • M::Vector{Matrix{Float32}}: Mask
  • Λ::Array{Complex{T},3}: Toeplitz kernel basis
  • num_fft_threads::Int = round(Int, Threads.nthreads()/size(U, 2)) or `round(Int, Threads.nthreads()/size(Λ, 1)): Number of Threads for FFT
source
MRFingerprintingRecon.NFFTNormalOpMethod
NFFTNormalOp(img_shape, trj, U; cmaps, verbose, num_fft_threads)
NFFTNormalOp(img_shape, Λ, kmask_indcs; cmaps)

Create normal operator of NFFT operator. Differentiate between functions exploiting a pre-calculated Toeplitz kernel basis Λ and the function which calculates Λ based on a passed trajectory trj.

Arguments

  • img_shape::Tuple{Int}: Image dimensions
  • traj::AbstractVector{<:AbstractMatrix{T}}: Trajectory
  • U::AbstractMatrix{Tc}: Basis coefficients of subspace
  • cmaps::AbstractVector{Matrix{ComplexF32}}=[ones(T, img_shape)]: Coil sensitivities
  • Λ::Array{Complex{T},3}: Toeplitz kernel basis
  • kmask_indcs::Vector{Int}: Sampling indices of Toeplitz mask
  • verbose::Boolean=false: Verbose level
  • num_fft_threads::Int=round(Int, Threads.nthreads()/size(U, 2)) or round(Int, Threads.nthreads()/size(Λ, 1)): Number of threads for FFT
source
MRFingerprintingRecon.calcCoilMapsMethod
calcCoilMaps(data, trj, img_shape; U, density_compensation, kernel_size, calib_size, eigThresh_1, eigThresh_2, nmaps, verbose)

Estimate coil sensitivity maps using ESPIRiT [1].

Arguments

  • data::AbstractVector{<:AbstractMatrix{Complex{T}}}: Complex dataset either as AbstractVector of matrices or single matrix. The optional outer vector defines different time frames that are combined using the subspace defined in U
  • trj::AbstractVector{<:AbstractMatrix{T}}: Trajectory with samples corresponding to the dataset either as AbstractVector of matrices or single matrix.
  • img_shape::NTuple{N,Int}: Shape of image

Keyword Arguments

  • U::Matrix = N==3 ? ones(size(data,1)) : I(1): Basis coefficients of subspace (only defined if data and trj are vectors of matrices)
  • density_compensation=:radial_3D: Values of :radial_3D, :radial_2D, :none, or of type AbstractVector{<:AbstractVector}
  • kernel_size=ntuple(_ -> 6, N): Kernel size
  • calib_size=ntuple(_ -> 24, N): Size of calibration region
  • eigThresh_1=0.01: Threshold of first eigenvalue
  • eigThresh_2=0.9: Threshold of second eigenvalue
  • nmaps=1: Number of estimated maps
  • verbose::Boolean=false: Verbosity level

return

  • cmaps::::Vector{<:Array{Complex{T}}}: Coil sensitivities as Vector of arrays

References

[1] Uecker, M., Lai, P., Murphy, M.J., Virtue, P., Elad, M., Pauly, J.M., Vasanawala, S.S. and Lustig, M. (2014), ESPIRiT—an eigenvalue approach to autocalibrating parallel MRI: Where SENSE meets GRAPPA. Magn. Reson. Med., 71: 990-1001. https://doi.org/10.1002/mrm.24751

source
MRFingerprintingRecon.calculateBackProjectionMethod
calculateBackProjection(data, trj, img_shape; U, density_compensation, verbose)
calculateBackProjection(data, trj, cmaps::AbstractVector{<:AbstractArray{cT,N}}; U, density_compensation, verbose)
calculateBackProjection(data, trj, cmaps_img_shape; U, density_compensation, verbose)
calculateBackProjection(data, trj, cmaps; U)

Calculate (filtered) backprojection

Arguments

  • data <: Union{AbstractVector{<:AbstractMatrix{cT}},AbstractMatrix{cT}}: Complex dataset either as AbstractVector of matrices or single matrix. The optional outer matrix defines different time frames that are reconstructed in the subspace defined in U.
  • trj <: AbstractVector{<:AbstractMatrix{T}}: Trajectory with samples corresponding to the dataset. For a Cartesian reconstruction, use T <: Int and define trj[it][idim,ik] ∈ (1, img_shape[idim]). If T <: Float, the NFFT is used.

One of the following arguments needs to be supplied

  • img_shape::NTuple{N,Int}: Shape of image; in this case, the data is reconstructed coilwise.
  • cmaps::::AbstractVector{<:AbstractArray{T}}: Coil sensitivities; in this case, the coils are added up to a single backprojection.

Optional Keyword Arguments

  • U::Matrix = I(length(data)) or = I(1)`: Basis coefficients of subspace (only defined if data and trj have different timeframes)
  • density_compensation = :none: Values of :radial_3D, :radial_2D, :none, or of type AbstractVector{<:AbstractVector}
  • verbose::Boolean = false: Verbosity level
source
MRFingerprintingRecon.calculateGoldenMeansMethod
calculateGoldenMeans()

Function to calculate 3D golden means [1].

References

[1] Chan, R.W., Ramsay, E.A., Cunningham, C.H. and Plewes, D.B. (2009), Temporal stability of adaptive 3D radial MRI using multidimensional golden means. Magn. Reson. Med., 61: 354-363. https://doi.org/10.1002/mrm.21837

source
MRFingerprintingRecon.grog_calibMethod
grog_calib(data, trj, Nr)

Perform GROG kernel calibration based on whole radial trajectory and passed data. Calibration follows the work on self-calibrating radial GROG (https://doi.org/10.1002/mrm.21565).

Arguments

  • data::AbstractVector{<:AbstractMatrix{cT}}: Complex dataset passed as AbstractVector of matrices
  • trj::Vector{Matrix{Float32}}: Trajectory with samples corresponding to the dataset passed as AbstractVector of matrices with Float32 entries
  • Nr::Int: Number of samples per read out
source
MRFingerprintingRecon.grog_gridding!Method
grog_gridding!(data, trj, lnG, Nr, img_shape)

Perform gridding of data based on pre-calculated GROG kernel.

Arguments

  • data::AbstractVector{<:AbstractMatrix{cT}}: Complex dataset passed as AbstractVector of matrices
  • trj::Vector{Matrix{Float32}}: Trajectory with samples corresponding to the dataset passed as AbstractVector of matrices with Float32 entries
  • lnG::Vector{Matrix{Float32}}: Natural logarithm of GROG kernel in all dimensions
  • Nr::Int: Number of samples per read out
  • img_shape::Tuple{Int}: Image dimensions

Output

  • trj::Vector{Matrix{Int32}}: Cartesian trajectory with the elements trj[it][idim,ik] ∈ (1, img_shape[idim])

Dimensions:

  • data: [timesteps][samples, spokes, coils, repetitions of sampling pattern]
  • trj: [timesteps][dims, samples, repetitions]
  • lnG: [dims][Ncoils, Ncoils]
source
MRFingerprintingRecon.kooshballMethod
kooshball(Nr, theta, phi; thetaRot, phiRot, delay)

Function to calculate kooshball trajectory.

Arguments

  • Nr::Int: Number of read out samples
  • theta::Array{Float,2}: Array with dimensions: Ncyc, Nt defining the angles theta for each cycle and timestep.
  • phi::Array{Float,2}: Array with dimensions: Ncyc, Nt defining the angles phi for each cycle and timestep.
  • thetaRot::Float = 0: Fixed rotation angle along theta
  • phiRot::Float = 0: Fixed rotation angle along phi
  • delay::Tuple{Float, Float, Float} = (0, 0, 0): Gradient delays in (HF, AP, LR)
source
MRFingerprintingRecon.kooshballGAMethod
kooshballGA(Nr, Ncyc, Nt; thetaRot, phiRot, delay)

Function to calculate golden means [1] based kooshball trajectory.

Arguments

  • Nr::Int: Number of read out samples
  • Ncyc::Int: Number of cycles
  • Nt::Int: Number of time steps in the trajectory
  • thetaRot::Float = 0: Fixed rotation angle along theta
  • phiRot::Float = 0: Fixed rotation angle along phi
  • delay::Tuple{Float, Float, Float}= (0, 0, 0): Gradient delays in (HF, AP, LR)

References

[1] Chan, R.W., Ramsay, E.A., Cunningham, C.H. and Plewes, D.B. (2009), Temporal stability of adaptive 3D radial MRI using multidimensional golden means. Magn. Reson. Med., 61: 354-363. https://doi.org/10.1002/mrm.21837

source
MRFingerprintingRecon.radial_grog!Method
radial_grog!(data, trj, Nr, img_shape)

Perform GROG kernel calibration and gridding [1] of data in-place. The trajectory is returned with integer values.

Arguments

  • data::AbstractVector{<:AbstractMatrix{cT}}: Complex dataset passed as AbstractVector of matrices
  • trj::Vector{Matrix{Float32}}: Trajectory with samples corresponding to the dataset passed as AbstractVector of matrices with Float32 entries
  • Nr::Int: Number of samples per read out
  • img_shape::Tuple{Int}: Image dimensions

Output

  • trj::Vector{Matrix{Int32}}: Cartesian trajectory with the elements trj[it][idim,ik] ∈ (-img_shape[idim]/2, img_shape[idim]/2-1)

References

[1] Seiberlich, N., Breuer, F., Blaimer, M., Jakob, P. and Griswold, M. (2008), Self-calibrating GRAPPA operator gridding for radial and spiral trajectories. Magn. Reson. Med., 59: 930-935. https://doi.org/10.1002/mrm.21565

source
MRFingerprintingRecon.traj_2d_radial_goldenratioMethod
traj_2d_radial_goldenratio(Nr, Ncyc, Nt; thetaRot, phiRot, delay, N)

Function to calculate 2D golden ratio based trajectory [1]. By modifying N also tiny golden angles [2] are supported.

Arguments

  • Nr::Int: Number of read out samples
  • Ncyc::Int: Number of cycles
  • Nt::Int: Number of time steps in the trajectory
  • thetaRot::Float = 0: Fixed rotation angle along theta
  • phiRot::Float = 0: Fixed rotation angle along phi
  • delay::Tuple{Float, Float, Float} = (0, 0, 0): Gradient delays in (HF, AP, LR)
  • N::Int = 1: Number of tiny golden angle

References

[1] Winkelmann S, Schaeffter T, Koehler T, Eggers H, Doessel O. An optimal radial profile order based on the Golden Ratio for time-resolved MRI. IEEE TMI 26:68–76 (2007) [2] Wundrak S, Paul J, Ulrici J, Hell E, Geibel MA, Bernhardt P, Rottbauer W, Rasche V. Golden ratio sparse MRI using tiny golden angles. Magn Reson Med 75:2372-2378 (2016)

source
MRFingerprintingRecon.traj_cartesianMethod
traj_cartesian(T, Nx, Ny, Nz, Nt)

Generate a 3D Cartesian trajectory. As input to functions using NFFTs, the trajectory can be defined as floats ∈ [-0.5, 0.5). For use with Cartesian FFTs, the trajectory consists of integers ∈ [1, N].

Arguments

  • Nx::Int: Number of readout samples
  • Ny::Int: Number of phase encoding lines
  • Nz::Int: Number of phase encoding lines (third dimension)
  • Nt::Int: Number of time steps in the trajectory

Optional Keyword Argument

  • T::Type=Int: Type of output trajectory. If T <: Float, trajectory is defined ∈ (-0.5, 0.5). If T <: Int, trajectory consists of values ∈ (1, N) instead.
source