GIRFReco.adjust_header!Method
adjust_header!(raw::RawAcquisitionData, recon_size, num_samples, interleave_number, single_slice)

Adjusts the header data for each interleave and slice of spiral diffusion RawAcquisitionData

Arguments

  • raw::RawAcquisitionData - RawAcquisitionData object
  • recon_size::Vector - Reconstruction matrix size
  • num_samples::Int - Number of samples per interleave
  • interleave_number::Int - Index of interleave for multi-shot acquisitionNumbers
  • single_slice::Bool - flag for single-slice reconstruction/acquisition
GIRFReco.apply_girf!Method
apply_girf!(a::AcquisitionData{T}, g::GirfApplier)

Applies the GIRF to the trajectories inside of a::AcquisitionData

Arguments

  • a::AcquisitionData{T} - AcquisitionData object
  • g::GirfApplier - GirfApplier object containing GIRF definition
GIRFReco.apply_k0!Method
apply_k0!(a::AcquisitionData{T}, g::GirfApplier)

Applies the K0 modulation due to imaging gradients to the data inside of a::AcquisitionData

Arguments

  • a::AcquisitionData{T} - AcquisitionData object
  • g::GirfApplier - GirfApplier containing GIRF definition
GIRFReco.calculate_b0_mapsMethod
calculate_b0_maps(me_data, slices, echotime_1, echotime_2)

Calculate B₀ map from the two images with different echo times via their phase difference (phase of imgTE2.*conj(imgTE1))

Arguments

  • me_data - [nX nY nZ 2 num_coils] 5D image array, 4th dim is echo time
  • slices::NTuple{num_slices,Int} - slice index vector (tuple?) for which map is computed
  • echotime_1::AbstractFloat - TE1 [ms]
  • echotime_2::AbstractFloat - TE2 [ms]
GIRFReco.check_acquisition_nodes!Method
check_acquisition_nodes!(a::AcquisitionData)

Validates processed AcquisitionData object to make sure that |kᵢ| < 0.5 ∀ i ∈ [1, Nₛ]

Arguments

  • a::AcquisitionData - AcquisitionData object
GIRFReco.check_profilesMethod
check_profiles(raw_data::RawAcquisitionData)

Sanity check of RawAcqData object by ploting all profiles to confirm its consistency with ISMRMRD file

Arguments

  • raw_data::RawAcquisitionData - RawAcquisitionData object
  • num_profiles_display - The number of profiles to be displayed
GIRFReco.do_k0_correction!Method
do_k0_correction!(raw_data, k0_phase_modulation, interleave)

Applies phase modulation due to 0th-order field fluctuations during the acquisition

Arguments

  • raw_data::RawAcquisitionData - RawAcquisitionData object
  • k0_phase_modulation::Matrix{Complex{T}} - Vector containing phase modulation measurements
  • interleave::Int - index of interleave
GIRFReco.get_slice_orderMethod
get_slice_order(r::RawAcquisitionData, sliceNum::Int, startProfile::Int, incProfile:Int)

Return a array of slice order index with ascending order of Z position.

e.g. For an interleaved pattern of slice position in RawAcquisitionData given below (in unit of mm): [-7, -3, 1, 5, 9, -9, -5, -1, 3, 7], the output will be [6, 1, 7, 2, 8, 3, 9, 4, 10, 5]

Arguments

  • r::RawAcquisitionData - A RawAcquisitionData that directly reads from original MRD file
  • sliceNum::Int - Total slice number that included in the RawAcquisitionData
  • startProfile::Int - Starting index of the profile in the RawAcqData for the first valid slice to be processed
  • incProfile::Int - Increment of profile index for the next valid slices

Output

  • orderedIndex - array with slice index of RawAcquisitionData with ascending order of position in Z.
GIRFReco.load_mapMethod
load_map(filename; do_split_phase::Bool = false)

Load calibration maps (sensitivity or B₀) from 4D NIfTI file(s)

For complex-valued data, magnitude and phase parts are stored in two files with suffix "magn" and "phase".

Arguments

  • filename::String - string filename with extension .nii, example "sensemap.nii"
  • do_split_phase::Bool=false - if true, data is saved in two nifti files with suffix "magn" and "phase", respectively to enable display in typical NIfTI viewers

Output

  • calib_map - [nX nY nZ {nChannels}] 4-D sensitivity or 3D B₀ map array
GIRFReco.merge_raw_interleavesMethod
merge_raw_interleaves(params, output_raw)

Merges multiple interleave data together from individually acquired interleave scans

Arguments

  • params - Dictionary
  • output_raw - Bool
GIRFReco.plot_reconstructionMethod
plot_reconstruction(images, slices_index, b0; is_slice_interleaved = false, rotation = 0)

Plots the magnitude and phase of the reconstructed images for a given slice or slices, along with a B₀ map if applicable

Arguments

  • images - Complex-valued images reconstructed using MRIReco.jl
  • slices_index::Vector{Int} - slices to plot
  • b0 - off-resonance map to plot along with images
  • is_slice_interleaved::Bool - for 2D scanning, indicate this value as true to make sure the slice order on the displayed results is correct
  • rotation::Int - Counterclock-wise rotation angle for each slice, should be a value from 0, 90, 180, 270 degrees
GIRFReco.plot_sense_mapsMethod
plot_sense_maps(sensitivity, num_channels; slice_index = 1)

Plots coil sensitivity maps from the channels, for a given number of num_channels plots on a given slice index.

Arguments

  • sensitivity - sensitivity maps, a 4D array: [nX, nY, nZ, nCoil]
  • num_channels - number of coils to be displayed.
  • slice_index - The index of the slice to be displayed (if multislice)
GIRFReco.preprocess_cartesian_dataMethod
preprocess_cartesian_data(r::RawAcquisitionData, do_save; filename = "data/processed_cartesian_file.h5")

Prepares Cartesian for reconstruction

Arguments

  • r::RawAcquisitionData{T} - RawAcquisitionData object
  • do_save::Boolean - Save the processed Cartesian data as a HDF5 file
  • filename - filename to save the preprocessed data
GIRFReco.remove_oversampling!Method
remove_oversampling!(raw::RawAcquisitionData)

Removes 2x readout oversampling in raw data along read-out dimension.

Arguments

  • raw::RawAcquisitionData{T} - RawAcquisitionData object
GIRFReco.save_mapMethod
save_map(filename, calib_map, resolution_mm; offset_mm = [0.0, 0.0, 0.0], do_split_phase::Bool = false, do_normalize::Bool = true)

Saves calibration maps (sensitivity or B₀) as 4D NIfTI file(s)

For complex-valued data, magnitude and phase can be split into separate files

Arguments

  • filename::String - string filename with extension .nii, example "sensemap.nii"
  • calib_map - [nX nY nZ {nChannels}] 4-D sensitivity or 3D B₀ map array
  • resolution_mm - resolution in mm, 3 element vector, e.g., [1.0, 1.0, 2.0]
  • offset_mm - isocenter offset in mm, default: [0.0, 0.0, 0.0]
  • do_split_phase::Bool = false - if true, data is saved in two nifti files with suffix "magn" and "phase", respectively to enable display in typical NIfTI viewers
  • do_normalize::Bool = true - if true, normalize the image by its magnitude maxima
GIRFReco.shift_kspace!Method
shift_kspace!(acqdata, shift)

This function applys additional phase ramps to k-space data to achive a given shift of center of image FOV in X and Y directions.

Perhaps this should be called shiftfov; however, since this function is modifying kspace data, it is named shiftkspace for now.

Arguments

  • acqdata::AcquisitionData{T} - AcquisitionData object to be modified
  • shift::AbstractVector - Vector containing shift with size [shiftX, shiftY]
GIRFReco.sync_traj_and_data!Method
sync_traj_and_data!(a::AcquisitionData)

Synchronizes k-space trajectory and sampled data as they do not usually have a common sampling rate

Arguments

  • raw_data::RawAcquisitionData - RawAcquisitionData object
  • traj::Trajectory - Trajectory object to be synchronized with data contained in raw_data
  • idx_crop::Int - Trajectory and Data may contain samples we don't want in the recon, usually at the end of acquisition. Ignore samples after idx_crop
  • interleave::Int - index of interleave
GIRFReco.validate_acq_data!Method
validate_acq_data!(a::AcquisitionData)

Validates processed AcquisitionData object after manipulation, etc...

Arguments

  • a::AcquisitionData - AcquisitionData object
GIRFReco.validate_siemens_mrd!Method
validate_siemens_mrd!(r::RawAcquisitionData)

Validates RawAcquisitionData object created from ISMRMRD format object

Arguments

  • r::RawAcquisitionData - RawAcquisitionData object
GIRFReco.read_gradient_text_fileMethod
read_gradient_text_file(filename, reconsize, delay)

Reads in text file containing gradient waveform information

Arguments

  • filename - filename (with full path) of text file with gradient waveform information
  • reconsize::Tuple{Int64,Int64,Int64} - size of reconstructed image (trailing dimension 1 for 2D acquisitions)
  • delay - delay in seconds from the nominal first sampling point to the actual first sampling point
GIRFReco.RMethod

R(x::Matrix{T}) Regularization function which penalizes roughness

Arguments

  • x::Matrix{T} - fieldmap estimate (in radians)
GIRFReco.estimate_b0_mapsFunction

estimateb0maps(im_data,slices, TE1,TE2,β,isrotated) Processes 3D volume data as output from MRIReco.reconstruction to estimate fieldmaps using the method presented by Funai and Fessler

Required Arguments

  • im_data - 5-D array with complex image data -> first
  • slices - vector of slices to process (must be within range of 3rd dimension of im_data)
  • TE1 - Echo time 1 [ms]
  • TE2 - Echo time 2 [ms]

Optional Arguments

  • isrotated - Boolean controlling whether to rotate the B0 maps to match the images or not (legacy feature)

Keyword Arguments

  • β - Regularization parameter controlling roughness penalty (larger = smoother)
  • reltol - early stopping criteria (exit if subsequent cost function change < reltol)
GIRFReco.ml_costMethod

ml_cost(x::Matrix{T},y::Matrix{Complex{T}},z::Matrix{Complex{T}}, β) Calculates the ML estimator cost between an estimated phase map and the underlying multi-echo scan data

Arguments

  • x::Matrix{T} - fieldmap estimate (in radians)
  • y::Matrix{Complex{T}} - Complex first-echo image data
  • z::Matrix{Complex{T}} - Complex second-echo image data
  • m::Matrix{T} - normalized weighting data (m ∈ [0,1]:= abs.(conj.(y).z)./maximum(abs.(conj.(y).z))), precomputed for speed
  • β - Regularization parameter controlling roughness penalty
GIRFReco.pcg_ml_est_fieldmapMethod

pcgmlest_fieldmap(y::AbstractMatrix{Complex{T}},z::AbstractMatrix{Complex{T}},β) Estimates the fieldmap using the method presented in https://doi.org/10.1109/tmi.2008.923956

Required Arguments

  • y::AbstractMatrix{Complex{T}} - Complex first-echo image data
  • z::AbstractMatrix{Complex{T}} - Complex second-echo image data

Optional Arguments

  • β - Regularization parameter controlling roughness penalty
  • reltol - early stopping criteria (exit if subsequent cost function change < reltol)