OCTOPUS package

OCTOPUS.ORC module

Methods for ORC (Off-Resonance Correction) and off-resonance simulation

Author: Marina Manso Jimeno

Last updated: 01/20/2021

OCTOPUS.ORC.CPR(dataIn, dataInType, kt, df, nonCart=None, params=None)

Off-resonance Correction by Conjugate Phase Reconstruction Maeda, A., Sano, K. and Yokoyama, T. (1988), Reconstruction by weighted correlation for MRI with time-varying gradients. IEEE Transactions on Medical Imaging, 7(1): 26-31. doi: 10.1109/42.3926

Parameters:
  • dataIn (numpy.ndarray) – k-space raw data or image data
  • dataInType (str) – Can be either ‘raw’ or ‘im’
  • kt (numpy.ndarray) – k-space trajectory.
  • df (numpy.ndarray) – Field map
  • nonCart (int) – Non-cartesian trajectory option. Default is None (Cartesian).
  • params (dict) – Sequence parameters. Default is None (Cartesian).
Returns:

M_hat – Corrected image data.

Return type:

numpy.ndarray

OCTOPUS.ORC.MFI(dataIn, dataInType, kt, df, Lx, nonCart=None, params=None)

Off-resonance Correction by Multi-Frequency Interpolation Man, L., Pauly, J. M. and Macovski, A. (1997), Multifrequency interpolation for fast off‐resonance correction. Magn. Reson. Med., 37: 785-792. doi:10.1002/mrm.1910370523

Parameters:
  • dataIn (numpy.ndarray) – k-space raw data or image data
  • dataInType (str) – Can be either ‘raw’ or ‘im’
  • kt (numpy.ndarray) – k-space trajectory
  • df (numpy.ndarray) – Field map
  • Lx (int) – L (frequency bins) factor
  • nonCart (int) – Non-cartesian trajectory option. Default is None (Cartesian).
  • params (dict) – Sequence parameters. Default is None.
Returns:

M_hat – Corrected image data.

Return type:

numpy.ndarray

OCTOPUS.ORC.add_or(M, kt, df, nonCart=None, params=None)

Forward model for off-resonance simulation

Parameters:
  • M (numpy.ndarray) – Image data
  • kt (numpy.ndarray) – k-space trajectory
  • df (numpy.ndarray) – Field map
  • nonCart (int , optional) – Cartesian/Non-Cartesian trajectory option. Default is None.
  • params (dict , optional) – Sequence parameters. Default is None.
Returns:

M_or – Off-resonance corrupted image data

Return type:

numpy.ndarray

OCTOPUS.ORC.add_or_CPR(M, kt, df, nonCart=None, params=None)

Forward model for off-resonance simulation. The number of fourier transforms = number of unique values in the field map.

Parameters:
  • M (numpy.ndarray) – Image data
  • kt (numpy.ndarray) – k-space trajectory
  • df (numpy.ndarray) – Field map
  • nonCart (int , optional) – Cartesian/Non-Cartesian trajectory option. Default is None.
  • params (dict , optional) – Sequence parameters. Default is None.
Returns:

M_or – Off-resonance corrupted image data

Return type:

numpy.ndarray

OCTOPUS.ORC.check_inputs_cartesian(dataInShape, dataInType, ktShape, dfShape)

Check dimensions of the inputs

Parameters:
  • dataInShape (tuple) – Shape of input data. Must be [NxN] for image data and [NlinesxNcolumns] for raw data.
  • dataInType (str) – Type of data: im or raw
  • ktShape (tuple) – Shape of k-space trajectory. Must be [NlinesxNcolumns].
  • dfShape (tuple) – Shape of field map. Must be [NxN] and match the image matrix size.
OCTOPUS.ORC.check_inputs_noncartesian(dataInShape, dataInType, ktShape, dfShape, params)

Check dimensions of the inputs

Parameters:
  • dataInShape (tuple) – Shape of input data. Must be [NxN] for image data and [NpointsxNshots] for raw data.
  • dataInType (str) – Type of data: im or raw
  • ktShape (tuple) – Shape of k-space trajectory. Must be [NpointsxNshots].
  • dfShape (tuple) – Shape of field map. Must be [NxN] and match the image matrix size.
  • params (dict) – Sequence parameters. Required elements are Npoints, Nshots, and N.
OCTOPUS.ORC.coeffs_MFI_lsq(kt, f_L, df_range, t_vector)

Finds the coefficients for Multi-frequency interpolation method by least squares approximation.

Parameters:
  • kt (numpy.ndarray) – K-space trajectory
  • f_L (numpy.ndarray) – Frequency segments array.
  • df_range (tuple) – Frequency range of the field map (minimum and maximum).
  • params (dict) – Sequence parameters.
Returns:

cL – Coefficients look-up-table.

Return type:

dict

OCTOPUS.ORC.correct_from_kdat(method, kdat, ktraj, fmap, params, cart_opt)

Corrects numerical simulation data directly from k-space data to avoid artifacts from ksp-im and im-ksp conversion.

Parameters:
  • method (str) – Correction method to apply: CPR, fsCPR or MFI
  • kdat (numpy.ndarray) – Corrupted k-space data
  • ktraj (numpy.ndarray) – k-space trajectory
  • fmap (numpy.ndarray) – Field map in Hertz
  • params (dict) – Sequence acquisition parameters
  • cart_opt – Cartesian option: 1 (spiral), ‘EPI’
Returns:

M_corr – Corrected image

Return type:

numpy.ndarray

OCTOPUS.ORC.find_nearest(array, value)

Finds the index of the value’s closest array element

Parameters:
  • array (numpy.ndarray) – Array of values
  • value (float) – Value for which the closest element has to be found
Returns:

idx – Index of the closest element of the array

Return type:

int

OCTOPUS.ORC.fs_CPR(dataIn, dataInType, kt, df, Lx, nonCart=None, params=None)

Off-resonance Correction by frequency-segmented Conjugate Phase Reconstruction Noll, D. C., Pauly, J. M., Meyer, C. H., Nishimura, D. G. and Macovskj, A. (1992), Deblurring for non‐2D fourier transform magnetic resonance imaging. Magn. Reson. Med., 25: 319-333. doi:10.1002/mrm.1910250210

Parameters:
  • dataIn (numpy.ndarray) – k-space raw data or image data
  • dataInType (str) – Can be either ‘raw’ or ‘im’
  • kt (numpy.ndarray) – k-space trajectory
  • df (numpy.ndarray) – Field map
  • Lx (int) – L (frequency bins) factor
  • nonCart (int) – Non-cartesian trajectory option. Default is None (Cartesian).
  • params (dict) – Sequence parameters. Default is None (Cartesian).
Returns:

M_hat – Corrected image data.

Return type:

numpy.ndarray

OCTOPUS.ORC.orc(M, kt, df)

Off-resonance correction for Cartesian trajectories

Parameters:
  • M (numpy.ndarray) – Cartesian k-space data
  • kt (numpy.ndarray) – Cartesian k-space trajectory
  • df (numpy.ndarray) – Field map
Returns:

M_hat – Off-resonance corrected image data

Return type:

numpy.ndarray