Fibers.ADCworkType
ADCwork{T}

Pre-allocated workspace for ADC fit computations

  • T::DataType : Data type for computations (default: Float32)
  • nvol::Int : Number of volumes in DWI series
  • ib0::Vector{Bool} : Indicator for b=0 volumes [nvol]
  • ipos::Vector{Vector{Bool}} : Indicator for volumes with positive DWI signal [nvol]
  • logs::Vector{Vector{T}} : Logarithm of DWI signal [nvol]
  • d::Vector{Vector{T}} : Linear system solution vector [2]
  • A::Matrix{T} : System matrix [nvol x 2]
  • pA::Matrix{T} : Pseudo-inverse of system matrix [2 x nvol]
Fibers.DSIType

Container for outputs of a DSI reconstruction

Fibers.DSIworkType
DSIwork{T}

Pre-allocated workspace for DSI reconstruction computations

  • T::DataType : Data type for computations (default: Float32)
  • nfft::Int : Size of FFT matrix
  • nvert::Int : Number of ODF vertices (on the half sphere)
  • dqr::T : Spacing of radial sampling for PDF integration
  • iq_ind::Vector{Int} : Indices of voxels in q-space sphere
  • H::Array{T, 3} : Hanning window [nfft^3]
  • p::Vector{Array{T, 3}} : Diffusion propagator [nfft^3]
  • X::Vector{Array{T, 3}} : DWI signals in q-space arrangement [nfft^3]
  • x::Vector{Array{Complex{T}, 3}} : FFT of q-space signal [nfft^3]
  • `xtmp::Vector{Array{Complex{T}, 3}} : FFT of q-space signal [nfft^3]
  • F::FFTW.cFFTWPlan : Linear operator for FFT
  • `qr2::Vector{T} : Q-space radii squared for PDF integration
  • `iqsubinterp::Array{T, 3} : Non-integer q-space indices for interpolation
  • isort::Vector{Vector{Int}} : Indices of ODF peak vertices (sorted) [nvert]
  • o::Vector{Vector{T}} : ODF amplitudes [nvert]
  • odf_peak::Vector{Vector{T}}: ODF amplitudes at local peaks [nvert]
  • faces::Matrix{Int} : ODF faces (on the half sphere) [nvert x 3]
Fibers.DTIType

Container for outputs of a DTI fit

Fibers.DTIworkType
DTIwork{T}

Pre-allocated workspace for DTI fit computations

  • T::DataType : Data type for computations (default: Float32)
  • nvol::Int : Number of volumes in DWI series
  • ib0::Vector{Bool} : Indicator for b=0 volumes [nvol]
  • ipos::Vector{Vector{Bool}} : Indicator for volumes with positive DWI signal [nvol]
  • logs::Vector{Vector{T}} : Logarithm of DWI signal [nvol]
  • d::Vector{Vector{T}} : Linear system solution vector [7]
  • A::Matrix{T} : System matrix [nvol x 7]
  • pA::Matrix{T} : Pseudo-inverse of system matrix [7 x nvol]
Fibers.GQIType

Container for outputs of a GQI fit

Fibers.GQIworkType
GQIwork{T}

Pre-allocated workspace for GQI reconstruction computations

  • T::DataType : Data type for computations (default: Float32)
  • nvol::Int : Number of volumes in DWI series
  • nvert::Int : Number of ODF vertices (on the half sphere)
  • isort::Vector{Vector{Int}} : Indices of ODF peak vertices (sorted) [nvert]
  • s::Vector{Vector{T}} : DWI signal [nvol]
  • o::Vector{Vector{T}} : ODF amplitudes [nvert]
  • odf_peak::Vector{Vector{T}}: ODF amplitudes at local peaks [nvert]
  • faces::Matrix{Int} : ODF faces (on the half sphere) [nvert x 3]
  • A::Matrix{T} : System matrix [nvert x nvol]
Fibers.LUTType

Container for segmentation and tract look-up tables

Fibers.MRIType

Container for header and image data of an MRI volume or volume series

Fibers.MRIMethod
MRI(vol::Array{T,N}) where T<:Number

Return an empty MRI structure

Fibers.MRIMethod
MRI(ref::MRI{A}, nframes::Integer=ref.nframes, datatype::DataType=eltype(A)) where A<:AbstractArray

Return an MRI structure whose header fields are populated based on a reference MRI structure ref, and whose image array is populated with zeros.

Optionally, the new MRI structure can be created with a different number of frames (nframes) than the reference MRI structure.

Fibers.ODFType

Vertices and faces for ODF computation

Fibers.RUMBAworkType
RUMBAwork{T}

Pre-allocated workspace for RUMBA-SD reconstruction computations

  • T::DataType : Data type for computations (default: Float32)
  • nmask::Int : Number of voxels in brain mask
  • ndir::Int : Number of unique diffusion-encoding directions
  • ncomp::Int : Number of components in fODF signal decomposition
  • nvox::NTuple{3,Int} : Dimensions of image volume (or nothing if TV not used)
  • isort::Vector{Vector{Int}} : Fiber ODF peak vertex indices (sorted)
  • fodf_peak::Vector{Vector{T}} : Fiber ODF amplitudes at local peaks
  • dodf_mat::Matrix{T} : Work matrices in diffusion ODF space
  • dodf_sig_mat::Matrix{T}
  • Iratio::Matrix{T}
  • fodf_mat::Matrix{T} : Work matrices in fiber ODF space
  • rl_mat::Matrix{T}
  • rl_mat2::Matrix{T}
  • tv_mat::Matrix{T}
  • tv_vol::Vector{Array{T,3}} : Work volumes in image space
  • Gx_vol::Vector{Array{T,3}}
  • Gy_vol::Vector{Array{T,3}}
  • Gz_vol::Vector{Array{T,3}}
  • Div_vol::Vector{Array{T,3}}
  • λ::Array{T,3}
  • σ2_vec::Matrix{T} : Noise variance and SNR by voxel
  • snr_vec::Matrix{T}
Fibers.StreamWorkType
StreamWork{T}

Pre-allocated workspace for streamline tractography

  • T::DataType : Data type (default: Float32)
  • len_min::Int : Minimum streamline length (default: 3)
  • len_max::Int : Maximum streamline length (default: max(nx,ny,nz))
  • cosang_thresh::T : Cosine of maximum bending angle (default: cosd(45))
  • step_size::T : Step length, in voxels (default: .5 voxel)
  • smooth_coeff::T : Vector smoothing coefficient, in 0-1
  • micro_search_cosang::T : Cosine of search angle (default: cosd(10))
  • micro_search_dist::Vector{Int} : Search distance, in voxels (default: 15)
  • strdims::Vector{Int} : In-plane dimensions for 2D LCM [2]
  • dxyz::Matrix{Int} : Coordinate increments for voxel jumps in LCM [3 4]
  • edgetype::Matrix{Int} : Voxel edges connected by i-th LCM element [2 10]
  • mask::BitArray{3} : Brain mask (voxel-wise) [nx ny nz]
  • ovecs::Array{T, 5} : Orientation vectors [3 nvec nx ny nz]
  • lcms::Array{T, 4} : LCM elements [nmat nx ny nz]
  • sublist::Vector{Vector{T}} : Subvoxel sampling offsets [nsub][3]
  • pos_now::Vector{Vector{T}} : Current (x, y, z) position [ncore][3]
  • vec_now::Vector{Vector{T}} : Current orientation vector [ncore][3]
  • pos_next::Vector{Vector{T}} : Next (x, y, z) position [ncore][3]
  • vec_next::Vector{Vector{T}} : Next orientation vector [ncore][3]
  • ivec_next::Vector{Int} : Index of next orientation vector [ncore]
  • cosang::Vector{Vector{T}} : Cosine of angle b/w vectors [ncore][nvec]
  • cosangabs::Vector{Vector{T}} : Absolute value of above [ncore][nvec]
  • dvox::Vector{Vector{Int}} : Change in voxel position [ncore][3]
  • lcm::Vector{Vector{T}} : LCM vector at single voxel [ncore][10]
  • isdiff::Vector{Bool} : Indicator of method difference [ncore]
  • str::Vector{Vector{Matrix{T}}} : Streamlines [ncore][nstr][3 npts]
  • flag::Vector{Vector{Vector{Bool}}} : Flag on streamlines [ncore][nstr][npts]
Fibers.TractType

Container for header and streamline data stored in .trk format

Fibers.TractMethod
Tract{T}(ref::MRI) where T<:Number

Return a Tract structure with data type T and header fields populated based on the reference MRI structure ref

Fibers.TractMethod
Tract{T}()

Return an empty Tract structure with data type T

Fibers.XformMethod
Xform{T}()

Return an empty Xform structure with data type T

Base.invMethod
Base.inv(xfm::Xform{T})

Invert a transform and return a new Xform structure

Base.showMethod
Base.show(mri::MRI; plane::Char='a', z::Union{Int64, Nothing}=nothing, t::Union{Int64, Nothing}=nothing, title::Union{String, Nothing}=nothing)

Show the z-th slice from the t-th frame of an MRI structure, along the specified plane ('a': axial; 's': sagittal; 'c': coronal).

Fibers.adc_fitMethod
adc_fit(dwi::MRI, mask::MRI)

Fit the apparent diffusion coefficient (ADC) to DWIs and return it as an MRI structure.

Fibers.adc_fitMethod
adc_fit(dwi::Vector{T}, W::ADCwork{T}) where T<:AbstractFloat

Fit the ADC for a single voxel

Fibers.ang2rotMethod
ang2rot(φ, θ)

Convert polar `φ` and azimuthal `θ` angles to a rotation matrix.

`φ` and `θ` must be in radians.
Fibers.besseli_ratioMethod
besseli_ratio(nu::Integer, z::Float32)

Evaluate the ratio of the modified Bessel functions of the first kind, of orders nu and nu-1, at the points in z.

Instead of direct evaluation, i.e., besseli(nu,z) / besseli(nu-1,z), compute with Perron's continued fraction equation, which is an order of magnitude faster.

Perron's continued fraction is substantially superior than Gauss' continued fraction when z >> nu, and only moderately inferior otherwise (Walter Gautschi and Josef Slavik, Math. Comp., 32(143):865-875, 1978).

Fibers.cart2polMethod
cart2pol(x, y)

Transform Cartesian coordinates (`x`, `y`) to polar coordinates (`φ`, `ρ`),
where `φ` is in radians.
Fibers.cart2sphMethod
cart2sph(x, y, z)

Transform Cartesian coordinates (`x`, `y`, `z`) to spherical coordinates
(`φ`, `θ`, `ρ`), where `φ` and `θ` are in radians.
Fibers.dispFunction
disp(mri::MRI, mrimod::Union{MRI, Nothing}=nothing)

Quick display of an MRI structure (an image slice and a summary of header info) in the terminal window

Arguments:

  • mri::MRI: the main image to display
  • mrimod:MRI: an optional image to modulate the main image by (e.g., an FA map to modulate a vector map)
Fibers.dsi_recFunction
dsi_rec(dwi::MRI, mask::MRI, odf_dirs::ODF=sphere_642, σ::Float32=Float32(1.25))

Perform diffusion spectrum imaging (DSI) reconstruction of DWIs, and return a DSI structure.

If you use this method, please cite: Wedeen, V. J., et al. (2005). Mapping complex tissue architecture with diffusion spectrum magnetic resonance imaging. Magnetic resonance in medicine, 54(6), 1377–1386. https://doi.org/10.1002/mrm.20642

Arguments

  • dwi::MRI: A series of DWIs, stored in an MRI structure with valid .bvec and .bval fields
  • mask::MRI: A brain mask volume, stored in an MRI structure
  • odf_dirs::ODF=sphere_642: The vertices and faces of the ODF tessellation, stored in an ODF structure
  • hann_width::Int=32: Width of Hanning window (in q-space voxels)

Output

In the DSI structure:

  • .pdf: PDF amplitudes on the half sphere
  • .odf: ODF amplitudes on the half sphere
  • .peak: Orientation vectors of the 3 peak ODF amplitudes
  • .qa: Quantitative anisotropy for each of the 3 peak orientations
Fibers.dsi_writeMethod
dsi_write(dsi::DSI, basename::String)

Write the volumes from a DSI structure that was created by dsi_rec() to files whose names start with the specified base name.

Fibers.dti_fitMethod
dti_fit(dwi::MRI, mask::MRI)

Fit tensors to DWIs and return a DTI structure.

Fibers.dti_fit_lsMethod
dti_fit_ls(dwi::MRI, mask::MRI)

Perform least-squares fitting of tensors from DWIs and return a DTI structure.

If you use this method, please cite: Peter Basser et al. (1994). Estimation of the effective self-diffusion tensor from the NMR spin echo. Journal of Magnetic Resonance Series B, 103(3), 247–254. https://doi.org/10.1006/jmrb.1994.1037

Fibers.dti_fit_lsMethod
dti_fit_ls(dwi::Vector{T}, W::DTIwork{T}) where T<:AbstractFloat

Perform least-squares fitting of tensor for a single voxel

Fibers.dti_mapsMethod
dti_maps(eigval1::T, eigval2::T, eigval3::T) where T<:AbstractFloat

Return the radial diffusivity (RD), mean diffusivity (MD), and fractional anisotropy (FA) given the 3 eigenvalues the diffusion tensor.

Fibers.dti_writeMethod
dti_write(dti::DTI, basename::String)

Write the volumes from a DTI structure that was created by dti_fit() to files whose names start with the specified base name.

Fibers.find_peaks!Method
find_peaks!(W::Union{GQIwork, DSIwork})

Find the vertices whose amplitudes are local peaks and return them sorted (assume that ODF amplitudes have been computed)

Fibers.get_tmp_pathFunction
get_tmp_path(tmpdir::String="")

Return path to a directory where temporary files can be stored.

Search for candidate directories in the following order:

  1. $TMPDIR: Check if environment variable is defined and directory exists
  2. $TEMPDIR: Check if environment variable is defined and directory exists
  3. /scratch: Check if directory exists
  4. /tmp: Check if directory exists
  5. tmpdir: Check if tmpdir argument was passed and directory exists

If none of the above exist, use current directory (./) and print warning.

Fibers.gqi_recFunction
gqi_rec(dwi::MRI, mask::MRI, odf_dirs::ODF=sphere_642, σ::Float32=Float32(1.25))

Perform generalized q-space imaging (GQI) reconstruction of DWIs, and return a GQI structure.

If you use this method, please cite: Fang-Cheng Yeh, et al. (2010). Generalized q-sampling imaging. IEEE Transactions on Medical Imaging, 29(9), 1626–1635. https://doi.org/10.1109/TMI.2010.2045126

Arguments

  • dwi::MRI: A series of DWIs, stored in an MRI structure with valid .bvec and .bval fields
  • mask::MRI: A brain mask volume, stored in an MRI structure
  • odf_dirs::ODF=sphere_642: The vertices and faces of the ODF tessellation, stored in an ODF structure
  • σ::Float32=Float32(1.25): Diffusion sampling length factor

Output

In the GQI structure:

  • .odf: ODF amplitudes on the half sphere
  • .peak: Orientation vectors of the 3 peak ODF amplitudes
  • .qa: Quantitative anisotropy for each of the 3 peak orientations
Fibers.gqi_writeMethod
gqi_write(gqi::GQI, basename::String)

Write the volumes from a GQI structure that was created by gqi_rec() to files whose names start with the specified base name.

Fibers.infoMethod
info(mri::MRI)

Show basic info from the header of an MRI structure

Fibers.isinmaskMethod
isinmask(point::Vector{T}, mask::BitArray)

Check if a point is in a mask array

Fibers.load_brukerMethod
load_bruker(indir::String; headeronly::Bool=false, reco::Integer=1)

Read Bruker image data from disk and return an MRI structure similar to the FreeSurfer MRI struct defined in mri.h.

Arguments

  • indir::String: Path to a Bruker scan directory (files called method, acqp,

pdata/1/reco, and pdata/1/2dseq are expected to be found in it)

Optional arguments

  • headeronly::Bool=false: If true, the pixel data are not read in.

  • reco::Integer=1: Number of image reconstruction to load, where 1 denotes the online reconstruction and higher numbers denote any additional offline reconstructions performed after acquisition.

Fibers.load_mghMethod
load_mgh(fname::String; slices::Union{Vector{Unsigned}, Nothing}=nothing, frames::Union{Vector{Unsigned}, Nothing}=nothing, headeronly::Bool=false)

Load a .mgh or .mgz file from disk.

Return:

  • The image data as an array

  • The 4x4 vox2ras matrix that transforms 0-based voxel indices to coordinates such that: [x y z]' = M * [i j k 1]'

  • MRI parameters as a vector: [tr, flip_angle, te, ti]

  • The image volume dimensions as a vector (useful when headeronly is true, in which case the image array will be empty)

Arguments

  • fname::String: Path to the .mgh/.mgz file

  • slices::Vector: 1-based slice numbers to load (default: read all slices)

  • frames::Vector: 1-based frame numbers to load (default: read all frames)

  • headeronly::Bool=false: If true, the pixel data are not read in.

Fibers.load_niftiMethod
load_nifti(fname::Stringr; headeronly::Bool=false) -> NIfTIheader, Array

Load a NIfTI (.nii or .nii.gz) volume from disk and return an array containing the image data and a NIfTIheader structure and the image data as an array.

Handle compressed NIfTI (nii.gz) by issuing an external Unix call to uncompress the file to a temporary file, which is then deleted.

The output NIfTIheader structure contains:

  • the units for each dimension of the volume [mm or msec], in the .pixdim field
  • the sform and qform matrices, in the .sform and .qform fields
  • the vox2ras matrix, which is the sform (if valid), otherwise the qform, in the .vox2ras field
Fibers.load_nifti_hdrMethod
load_nifti_hdr(fname::String) -> NIfTIheader

Load the header of a .nii volume from disk.

Return a NIfTIheader structure, where units have been converted to mm and msec, the sform and qform matrices have been computed and stored in the .sform and .qform fields, and the .vox2ras field has been set to the sform (if valid), then the qform (if valid).

Assume that the input file is uncompressed (compression is handled in the wrapper load_nifti()).

Handle data with more than 32k cols by looking at the .dim field of the header. If dim[2] = -1, then the .glmin field contains the number of columns. This is FreeSurfer specific, for handling surfaces. When the total number of spatial voxels equals 163842, then reshape the volume to 163842x1x1xnframes. This is for handling the 7th order icosahedron used by FS group analysis.

Fibers.mri_filenameFunction
mri_filename(fstring::String, checkdisk::Bool=true)

Return a valid file name, file stem, and file extension, given a string fstring that can be either a file name or a file stem.

Valid extensions are: mgh, mgz, nii, nii.gz. Thus a file name is expected to have the form stem.{mgh, mgz, nii, nii.gz}.

If fstring is a file name, then the stem and extension are determined from fstring.

If fstring is a file stem and checkdisk is true (default), then the file name and extension are determined by searching for a file on disk named fstring.{mgh, mgz, nii, nii.gz}, where the possible extensions are searched for in this order. If no such file is found, then empty strings are returned.

Fibers.mri_readMethod
mri_read(inbase::String, type::DataType; headeronly::Bool=false, permutedata::Bool=false)

Read a set of image files with base name inbase, which were generated by an analysis, into a struct of a custom type (e.g., DTI)

Fibers.mri_readMethod
mri_read(infile::String; headeronly::Bool=false, permutedata::Bool=false, reco::Integer=1)

Read an image volume from disk and return an MRI structure similar to the FreeSurfer MRI struct defined in mri.h.

Arguments

  • infile::String: Path to the input, which can be:
    1. An MGH file, e.g., f.mgh or f.mgz
    2. A NIfTI file, e.g., f.nii or f.nii.gz
    3. A file stem, in which case the format and full file name are determined by finding a file on disk named infile.{mgh, mgz, nii, nii.gz}
    4. A Bruker scan directory, which is expected to contain the following files: method, acqp, pdata/1/reco, pdata/1/2dseq

Optional arguments

  • headeronly::Bool=false: If true, the pixel data are not read in.

  • permutedata::Bool==false: If true, the first two dimensions of the image volume are permuted in the .vol, .volsize, and .volres fields of the output structure (this is the default behavior of the MATLAB version). The permutedata will not affect the vox2ras matrices, which always map indices in the (column, row, slice) convention.

  • reco::Integer=1: (Only relevant to reading Bruker scan directories) Number of image reconstruction to load, where 1 denotes the online reconstruction and higher numbers denote any additional offline reconstructions performed after acquisition.

Output

In the MRI structure:

  • Times are in ms and angles are in radians.

  • vox2ras0: This field contains the vox2ras matrix that converts 0-based (column, row, slice) voxel indices to (x, y, z) coordinates. This is the also the matrix that mri_write() uses to derive all geometry information (e.g., direction cosines, voxel resolution, P0). Thus if any other geometry fields of the structure are modified, the change will not be reflected in the output volume.

  • vox2ras1: This field contains the vox2ras matrix that converts 1-based (column, row, slice) voxel indices to (x, y, z) coordinates.

  • niftihdr: If the input file is NIfTI, this field contains a NIfTIheader structure.

Fibers.mri_read_bfiles!Method
mri_read_bfiles!(dwi::MRI, infile1::String, infile2::String)

Set the .bval and .bvec fields of the MRI structure dwi, by reading a DWI b-value table and gradient table from the text files infile1 and infile2. The two input files can be specified in any order. The gradient table file must contain 3 times as many entries as the b-value table file.

Return the b-value table as a vector of size n and the gradient table as a matrix of size (n, 3).

Fibers.mri_read_bfilesMethod
mri_read_bfiles(infile1::String, infile2::String)

Read a DWI b-value table and gradient table from text files infile1 and infile2. The two input files can be specified in any order. The gradient table file must contain 3 times as many entries as the b-value table file.

Return the b-value table as a vector of size n and the gradient table as a matrix of size (n, 3).

Fibers.mri_set_geometry!Method
mri_set_geometry!(mri::MRI)

Set header fields related to volume geometry in an MRI structure.

These are redundant fields that can be derived from the vox2ras0, width, height, depth fields. They are in the MRI struct defined in mri.h, as well as the MATLAB version of this structure, so they have been added here for completeness.

Note: mri_write() derives all geometry information (i.e., direction cosines, voxel resolution, and P0 from vox2ras0. If any of the redundant geometry elements below are modified, the change will not be reflected in the output volume.

Fibers.mri_writeFunction
mri_write(mri::MRI, outfile::String, datatype::DataType=Float32)

Write an MRI volume to disk. Return true is an error occurred (i.e., the number of bytes written were not as expected based on the size of the volume).

Arguments

  • mri::MRI: A structure like that returned by mri_read(). The geometry (i.e., direction cosines, voxel resolution, and P0) are all recomputed from mri.vox2ras0. So, if a method has changed one of the other fields, e.g., mri.x_r, this change will not be reflected in the output volume.

  • outfile::String: Path to the output file, which can be:

    1. An MGH file, e.g., f.mgh or f.mgz (uncompressed or compressed)
    2. A NIfTI file, e.g., f.nii or f.nii.gz (uncompressed or compressed).
  • datatype::DataType=eltype(mri.vol): Only applies to NIfTI and can be UInt8, Int16, Int32, Float32, Float64, Int8, UInt16, UInt32. By default, the native data type of the volume array is preserved when writing to disk.

Fibers.pol2cartMethod
pol2cart(φ, ρ)

Transform polar coordinates (`φ`, `ρ`) to Cartesian coordinates (`x`, `y`),
where `φ` is in radians.
Fibers.rumba_peaks!Method

Find fODF peaks, given an fODF amplitude vector and an amplitude threshold.

Return the number of peaks found.

Fibers.rumba_recMethod
rumba_rec(dwi::MRI{T}, mask::MRI, odf_dirs::ODF=sphere_724, niter::Integer=600, λ_para::T=T(1.7*10^-3), λ_perp::T=T(0.2*10^-3), λ_csf::T=T(3.0*10^-3), λ_gm::T=T(0.8*10^-4), ncoils::Integer=1, coil_combine::String="SMF-SENSE", ipat_factor::Integer=1, use_tv::Bool=true) where T <: AbstractFloat

Perform robust and unbiased model-based spherical deconvolution (RUMBA-SD) reconstruction of DWIs, and return a RUMBASD structure.

If you use this method, please cite: Erick J. Canales-Rodríguez, et al. (2015). Spherical deconvolution of multichannel diffusion MRI data with non-Gaussian noise models and spatial regularization. PLoS ONE, 10(10), e0138910. https://doi.org/10.1371/journal.pone.0138910

Arguments

  • dwi::MRI{T}: A series of DWIs, stored in an MRI structure with valid .bvec and .bval fields
  • mask::MRI: A brain mask volume, stored in an MRI structure
  • odf_dirs::ODF=sphere_724: The vertices and faces of the ODF tessellation, stored in an ODF structure
  • λ_para::T=T(1.7*10^-3): Axial diffusivity in white-matter voxels with a single fiber population
  • λ_perp::T=T(0.2*10^-3): Radial diffusivity in white-matter voxels with a single fiber population
  • λ_csf::T=T(3.0*10^-3): Mean diffusivity in CSF voxels
  • λ_gm::T=T(0.8*10^-4): Mean diffusivity in gray-matter voxels
  • ncoils::Integer=1: Number of receive coil elements, if the DWIs were collected with parallel imaging, or 1 otherwise
  • coil_combine::String="SMF-SENSE": Method that was used to combine signals from receive coil elements (possible values: "SMF-SENSE" or "SoS-GRAPPA"; has no effect if ncoils equals 1)
  • ipat_factor::Integer=1: Acceleration factor, if the DWIs were collected with parallel imaging, or 1 otherwise
  • use_tv::Bool=true: If true, include a total-variation regularization term in the FOD reconstruction

Output

In the RUMBASD structure:

  • .fodf: Volume fractions of anisotropic compartments (one per vertex)
  • .fgm: Volume fraction of GM isotropic compartment
  • .fcsf: Volume fraction of CSF isotropic compartment
  • .peak: Orientation vectors of the 5 peak ODF amplitudes
  • .gfa: Generalized fractional anisotropy
  • .var: Estimated noise variance
  • .snr_mean: Estimated mean SNR
  • .snr_std: Estimated SNR standard deviation
Fibers.rumba_sd_iterate!Method
fodf:     Volume fractions of anisotropic compartments
fgm:      Volume fraction of GM isotropic component
fcsf:     Volume fraction of CSF isotropic component
var:      Noise variance
snr_mean: Estimated mean snr
snr_std:  Estimated SNR standard deviation
Fibers.rumba_tv!Method
Total variation (TV) regularization term

Read the fODF amplitudes for a single vertex from a 3D volume in the workspace and compute the TV regularization term in place, overwriting the volume

Fibers.rumba_writeMethod
rumba_write(rumba::RUMBASD, basename::String)

Write the volumes from a RUMBASD structure that was created by rumba_rec() to files whose names start with the specified base name.

Fibers.save_mghFunction
save_mgh(vol::Array, fname::String, M::Matrix=Diagonal(ones(4)), mr_parms::Vector=zeros(4))

Write an MRI volume to a .mgh or .mgz file. Return true is an error occurred (i.e., the number of bytes written were not as expected based on the size of the volume).

Arguments

  • vol::Array: the image data

  • fname::String: path to the output file

  • M::Matrix: the 4x4 vox2ras transform such that [x y z]' = M * [i j k 1]', where the voxel indices (i, j, k) are 0-based.

  • mr_parms::Vector: a vector of MRI parameters, [tr, flip_angle, te, ti]

Fibers.save_niftiMethod
save_nifti(hdr::NIfTIheader, vol::Array{T}, fname::String) where T<:Number

Write an MRI volume to a .nii or .nii.gz file. Return true is an error occurred (i.e., the number of bytes written were not as expected based on the size of the volume).

Arguments

  • hdr::NIfTIheader: a NIfTI header structure

  • vol::Array{T}: an array that contains the image data

  • fname::String: path to the output file

Handle data structures with more than 32k cols by setting hdr.dim[2] = -1 and hdr.glmin = ncols. This is FreeSurfer specific, for handling surfaces. The exception to this is when the total number of spatial voxels equals 163842, then the volume is reshaped to 27307x1x6xnframes. This is for handling the 7th order icosahedron used by FS group analysis.

Fibers.sph2cartMethod
sph2cart(φ, θ, ρ)

Transform spherical coordinates (`φ`, `θ`, `ρ`) to Cartesian coordinates
(`x`, `y`, `z`).

`φ` and `θ` must be in radians.
Fibers.str_add!Method
str_add!(tr::Tract{T}, xyz::Vector{Matrix{Txyz}}, scalars::Union{Vector{Matrix{Ts}}, Vector{Vector{Ts}}, Nothing}=nothing, properties::Union{Matrix{Tp}, Vector{Tp}, Nothing}=nothing) where {T<:Number, Txyz<:Number, Ts<:Number, Tp<:Number}

Append new streamlines to a Tract structure of data type T

Required inputs

  • tr::Tract{T}: Tract structure that the streamlines will be added to
  • xyz::Vector{Matrix}: Voxel coordinates of the points on the new streamlines [nstr][3 x npts]

Optional inputs (required only if Tract structure contains them)

  • scalars::Union{Vector{Matrix}, Vector{Vector}, Nothing}: Scalars associated with each point on the new streamlines [nstr][nscalar x npts] or (if nscalar == 1) [nstr][npts]
  • properties::Union{Matrix, Vector, Nothing}: Properties associated with each of the new streamlines [nprop x nstr] or (if nprop == 1) [nstr]
Fibers.str_mergeMethod
str_merge(tr1::Tract{T}, tr2::Tract{T}...) where T<:Number

Merge streamlines from two or more Tract structures and return a new Tract structure

Fibers.str_xformMethod
str_xform(xfm::Xform{T}, tr::Tract{T}) where T<:Number

Apply a transform to streamline coordinates and return a new Tract structure

Fibers.streamMethod
stream(ovec::Union{MRI,Vector{MRI}}; odf::Union{MRI,Nothing}=nothing}, f::Union{MRI,Vector{MRI},Nothing}=nothing, f_thresh::Real=.03, fa::Union{MRI,Nothing}=nothing, fa_thresh::Real=.1, mask::Union{MRI,Nothing}=nothing, lcms::Union{MRI,Nothing}=nothing, lcm_thresh::Real=.099, seed::Union{MRI,Nothing}=nothing, nsub::Union{Integer,Nothing}=3, len_min::Integer=3, len_max::Integer=(isa(ovec,MRI) ? maximum(ovec.volsize) : maximum(ovec[1].volsize)), ang_thresh::Union{Real,Nothing}=45, step_size::Union{Real,Nothing}=.5, smooth_coeff::Union{Real,Nothing}=.2, verbose::Bool=false, search_dist::Integer=15, search_ang::Number=10)

Streamline tractography

Arguments

  • ovec::Union{MRI,Vector{MRI}} : Orientation vectors [nvec][nx ny nz 3]

Optional arguments

  • odf::MRI : ODFs nx ny nz 3 nvert
  • f::Union{MRI,Vector{MRI}} : Vector amplitudes (or volume fractions), for vector-wise masking [nvec][nx ny nz]
  • f_thresh::Real : Minimum vector amplitude (or volume fraction), for vector-wise masking (default: .03)
  • fa::MRI : FA (or other microstructural measure), for voxel-wise masking (default: none)
  • fa_thresh::Real : Minimum FA (or other microstructural measure), for voxel-wise masking (default: .1)
  • mask::MRI : Brain mask [nx ny nz]
  • seed::Union{MRI,Nothing} : Seed mask nx ny nz
  • nsub::Integer : Number of subvoxel samples per seed voxel (default: 3)
  • len_min::Integer : Minimum streamline length (default: 3)
  • len_max::Integer : Maximum streamline length (default: max(nx,ny,nz))
  • ang_thresh::Real : Maximum bending angle, in degrees (default: 45)
  • step_size::Real : Step length, in voxels (default: .5 voxel)
  • smooth_coeff::Real : Vector smoothing coefficient, in 0-1
  • search_dist::Integer : Micro search distance, in voxels (default: 15)
  • search_ang::Integer : Micro search angle, in degrees (default: 10)
  • lcms::Union{MRI,Nothing} : LCMs (default: none) [ny nx nz nmat]
  • lcm_thresh::Real : Minimum LCM coefficient (default: .099)
  • verbose::Bool : Count differences b/w propagation methods (default: no)

Outputs

In the Tract structure

  • .str: Voxel coordinates of points along streamlines
  • .scalar: Indicator function of a method difference (only if verbose==true)
Fibers.tensor_modelMethod
tensor_model(φ::Number, θ::Number, λ::Vector, b::Vector, g::Matrix, s0::Number)

Compute the expected diffusion-weighted signal in a voxel, assuming that
diffusion can be modeled by a tensor with orientation angles `φ`, `θ` and
eigenvalues `λ` and that the signal was acquired with b-values `b` and 
gradient vectors `g`, and that the non-diffusion-weighted signal is `s0`.
Fibers.trk_readMethod
trk_read(infile::String)

Read tractography streamlines from infile and return a Tract structure

Input file must be in .trk format, see: http://trackvis.org/docs/?subsect=fileformat

Fibers.trk_writeMethod
 trk_write(tr::Tract, outfile::String)

Write a Tract structure to a file in the .trk format

Return true if an error occurred (i.e., the number of bytes written was not the expected based on the size of the Tract structure)

Fibers.view_axesMethod
view_axes(vox2ras::Matrix, plane::Char)

Given the vox2ras matrix of an image volume, return the axes along which the volume has to be displayed to view the specified plane ('a': axial; 's': sagittal; 'c': coronal).

Fibers.vol_to_rgbFunction
vol_to_rgb(vol::Array, maxint::Union{Number, Nothing}=nothing)

Convert an image array to an RGB/Gray array for display.

Determine how the image should be displayed:

  • If all the values are IDs in the FreeSurfer color look-up table, the image is assumed to be is a segmentation map and is converted to RGB based on the FreeSurfer color look-up table.
  • If the image has size 3 in any dimension, and the sum of squares of the values in that dimension is approx. 1 everywhere, the image is assumed to be a vector map and is converted to RGB based on vector orientation.
  • Otherwise, the image is assumed to be a generic intensity map and is converted to grayscale, optionally clamping intensities above maxint.

Return an array of RGB/Gray values the same size as the input vol.

Fibers.vox2ras_0to1Method
vox2ras_0to1(M0::Matrix)

Convert a 0-based vox2ras matrix M0 to a 1-based vox2ras matrix such that:

Pxyz = M0 * [c r s 1]' = M1 * [c+1 r+1 s+1 1]'

Fibers.vox2ras_tkregMethod
vox2ras_tkreg(voldim::Vector, voxres::Vector)

Return a 0-based vox2ras transform of a volume that is compatible with the registration matrix produced by tkregister. May work with MINC xfm.

Arguments

  • voldim = [ncols; nrows; nslices ...]
  • volres = [colres; rowres; sliceres ...]
Fibers.vox2ras_to_orientMethod
vox2ras_to_orient(vox2ras::Matrix)

Convert a vox2ras matrix to a 3-character array indicating image orientation, e.g., RAS, LIA, etc.

Fibers.vox2ras_to_qformMethod
vox2ras_to_qform(vox2ras::Matrix)

Convert a vox2ras matrix to NIfTI qform parameters. The vox2ras should be 6 DOF.

Return the following NIfTI header fields:

  • hdr.quatern_b
  • hdr.quatern_c
  • hdr.quatern_d
  • hdr.qoffset_x
  • hdr.qoffset_y
  • hdr.qoffset_z
  • hdr.pixdim[1]

From DG's vox2rasToQform.m: This code mostly just follows CH's mriToNiftiQform() in mriio.c

Fibers.xfm_apply!Method
xfm_apply!(outpoint::Array{To}, xfm::Xform{Tx}, inpoint::Array{Ti})

Apply a transform specified in an Xform structure to N points specified in a voxel coordinate array of length 3*N, in place

Fibers.xfm_applyMethod
xfm_apply(xfm::Xform{Tx}, point::Array{Ti})

Apply a transform specified in an Xform structure to N points specified in a voxel coordinate array of length 3*N, and return the transformed points as a voxel coordinate array of the same shape as the input array

Fibers.xfm_composeMethod
xfm_compose(xfm1::Xform{T}, xfm2::Xform{T}...)

Compose two transforms and return a new Xform structure

Note that the last argument of this function is the "innermost" transform, i.e., the one that would be applied first to the input data, and the first argument is the "outermost" transform, i.e., the one that would be applied last: output = xfm1 * xfm2 * ... * input

Fibers.xfm_readFunction
xfm_read(ltafile::String, T::DataType=Float32)

Read a transform from a .lta file and return an Xform structure

Fibers.xfm_readFunction
xfm_read(matfile::String, inref::MRI, outref::MRI, T::DataType=Float32)

Read a transform from a .mat file and return an Xform structure

Reference volumes in the input and output space of the transform must also be specified, as MRI structures

Fibers.xfm_rotate!Method
xfm_rotate!(outpoint::Vector{To}, xfm::Xform{Tx}, inpoint::Vector{Ti})

Apply the rotation component of a transform specified in an Xform structure to a point specified in a voxel coordinate vector of length 3, in place

Fibers.xfm_rotateMethod
xfm_rotate(xfm::Xform{Tx}, point::Vector{Ti})

Apply the rotation component of a transform specified in an Xform structure to a point specified in a voxel coordinate vector of length 3, and return the transformed point as a voxel coordinate vector of length 3