sdcflows.transform module#

The \(B_0\) unwarping transform formalism.

This module implements a data structure to represent the displacements field associated to the deformations caused by susceptibility-derived distortions. This implementation attempts to provide a single representation of such distortions independently of the estimation strategy (see Methods and implementation).

That is achieved by implementing multi-level B-Spline cubic interpolators. For one given level, a function \(f(\mathbf{s})\) can be represented as a linear combination of tensor-product cubic B-Spline basis (\(\Psi^3(\mathbf{k}, \mathbf{s})\); see Eq. \(\eqref{eq:2}\)).

\[f(\mathbf{s}) = \sum_{k_1} \sum_{k_2} \sum_{k_3} c(\mathbf{k}) \Psi^3(\mathbf{k}, \mathbf{s}). \label{eq:1}\tag{1}\]
class sdcflows.transform.B0FieldTransform(coeffs=None)[source]#

Bases: object

Represents and applies the transform to correct for susceptibility distortions.

apply(moving: SpatialImage, pe_dir: str | Sequence[str], ro_time: float | Sequence[float], xfms: Sequence[ndarray] | None = None, xfm_data2fmap: ndarray | None = None, approx: bool = True, order: int = 3, mode: str = 'constant', cval: float = 0.0, prefilter: bool = True, output_dtype: str | dtype | None = None, num_threads: int = None, allow_negative: bool = False)[source]#

Apply a transformation to an image, resampling on the reference spatial object.

Handles parallelization to resample 4D images.

Parameters:
  • moving (SpatialImage) – The image object containing the data to be resampled in reference space

  • pe_dir (str or list of str) – A valid PhaseEncodingDirection metadata value.

  • ro_time (float or list of float) – The total readout time in seconds.

  • xfms (None or list) – A list of 4×4 matrices, each one formalizing the estimated head motion alignment to the scan’s reference. Therefore, each of these matrices express the transform of every voxel’s RAS (physical) coordinates in the image used as reference for realignment into the coordinates of each of the EPI series volume.

  • xfm_data2fmap (numpy.ndarray) – Transform that maps coordinates on the target_reference onto the fieldmap reference (that is, the linear transform through which the fieldmap can be resampled in register with the target_reference). In other words, xfm_data2fmap is the result of calling a registration tool such as ANTs configured for a linear transform with at most 12 degrees of freedom and called with the image carrying target_affine as reference and the fieldmap reference as moving. The result of such a registration framework is an affine (our xfm_data2fmap here) that maps coordinates in reference (target) RAS onto the fieldmap RAS.

  • approx (bool) – If True, do not reconstruct the B-Spline field directly on the target (which will likely not be aligned with the fieldmap’s grid), but rather use the fieldmap’s grid and then use just regular interpolation.

  • order (int, optional) – The order of the spline interpolation, default is 3. The order has to be in the range 0-5.

  • mode ({‘constant’, ‘reflect’, ‘nearest’, ‘mirror’, ‘wrap’}, optional) – Determines how the input image is extended when the resamplings overflows a border. Default is ‘constant’.

  • cval (float, optional) – Constant value for mode='constant'. Default is 0.0.

  • prefilter (bool, optional) – Determines if the image’s data array is prefiltered with a spline filter before interpolation. The default is True, which will create a temporary float64 array of filtered values if order > 1. If setting this to False, the output will be slightly blurred if order > 1, unless the input is prefiltered, i.e. it is the result of calling the spline filter on the original input.

  • output_dtype (str or dtype) – Override the output data type, instead of propagating it from the moving image.

  • num_threads (int) – The maximum number of parallel resamplings at any given time of execution. Use this parameter to set an upper bound to memory utilization.

  • allow_negative (bool) – Remove negative values introduced in interpolation (may happen for nonnegative data when order \(\gt\) 3). Set this value to True if your moving image does have negative values.

Returns:

resampled – The data imaged after resampling to reference space.

Return type:

SpatialImage

coeffs#

B-Spline coefficients (one value per control point).

fit(target_reference: SpatialImage, xfm_data2fmap: ndarray | None = None, approx: bool = True) bool[source]#

Generate the interpolation matrix (and the VSM with it).

Implements Eq. \(\eqref{eq:1}\), interpolating \(f(\mathbf{s})\) for all voxels in the target-image’s extent.

Parameters:
  • target_reference (SpatialImage) – The image object containing a reference grid (same as that of the data to be resampled). If a 4D dataset is provided, then the fourth dimension will be dropped.

  • xfm_data2fmap (numpy.ndarray) – Transform that maps coordinates on the target_reference onto the fieldmap reference (that is, the linear transform through which the fieldmap can be resampled in register with the target_reference). In other words, xfm_data2fmap is the result of calling a registration tool such as ANTs configured for a linear transform with at most 12 degrees of freedom and called with the image carrying target_affine as reference and the fieldmap reference as moving. The result of such a registration framework is an affine (our xfm_data2fmap here) that maps coordinates in reference (target) RAS onto the fieldmap RAS.

  • approx (bool) – If True, do not reconstruct the B-Spline field directly on the target (which will likely not be aligned with the fieldmap’s grid), but rather use the fieldmap’s grid and then use just regular interpolation.

Returns:

updatedTrue if the internal field representation was fit, False if cache was valid and will be reused.

Return type:

bool

mapped#

A cache of the interpolated field in Hz (i.e., the fieldmap mapped on to the target image we want to correct).

to_displacements(ro_time, pe_dir, itk_format=True)[source]#

Generate a NIfTI file containing a displacements field transform compatible with ITK/ANTs.

The displacements field can be calculated following Eq. (2) in the fieldmap fitting section.

Parameters:
  • ro_time (float) – The total readout time in seconds (only if vsm=False).

  • pe_dir (str) – The PhaseEncodingDirection metadata value (only if vsm=False).

Returns:

spatialimage – A NIfTI 1.0 object containing the distortion.

Return type:

nibabel.nifti.Nifti1Image

sdcflows.transform.disp_to_fmap(xyz_nii, epi_nii, ro_time, pe_dir, itk_format=True)[source]#

Convert a displacements field into a fieldmap in Hz.

This is the inverse operation to the previous function.

Parameters:
  • xyz_nii (nibabel.Nifti1Image) – Displacements field in NIfTI format.

  • epi_nii (nibabel.Nifti1Image) – Source EPI image.

  • ro_time (float) – The total readout time of the EPI image in seconds.

  • pe_dir (str) – The PhaseEncodingDirection metadata value of the EPI image.

Returns:

spatialimage – A NIfTI 1.0 object containing the field in Hz.

Return type:

nibabel.nifti.Nifti1Image

sdcflows.transform.fmap_to_disp(fmap_nii, ro_time, pe_dir, itk_format=True)[source]#

Convert a fieldmap in Hz into an ITK/ANTs-compatible displacements field.

The displacements field can be calculated following Eq. (2) in the fieldmap fitting section.

Parameters:
  • fmap_nii (os.pathlike) – Path to a voxel-shift-map (VSM) in NIfTI format

  • ro_time (float) – The total readout time in seconds

  • pe_dir (str) – The PhaseEncodingDirection metadata value

Returns:

spatialimage – A NIfTI 1.0 object containing the distortion.

Return type:

nibabel.nifti.Nifti1Image

sdcflows.transform.grid_bspline_weights(target_nii, ctrl_nii, dtype='float32')[source]#

Evaluate tensor-product B-Spline weights on a grid.

For each of the N input samples \((s_1, s_2, s_3)\) and K control points or knots \(\mathbf{k} =(k_1, k_2, k_3)\), the tensor-product cubic B-Spline kernel weights are calculated:

\[\Psi^3(\mathbf{k}, \mathbf{s}) = \beta^3(s_1 - k_1) \cdot \beta^3(s_2 - k_2) \cdot \beta^3(s_3 - k_3), \label{eq:2}\tag{2}\]

where each \(\beta^3\) represents the cubic B-Spline for one dimension. The 1D B-Spline kernel implementation uses numpy.piecewise, and is based on the closed-form given by Eq. (6) of [Unser1999].

By iterating over dimensions, the data samples that fall outside of the compact support of the tensor-product kernel associated to each control point can be filtered out and dismissed to lighten computation.

Finally, the resulting weights matrix \(\Psi^3(\mathbf{k}, \mathbf{s})\) can easily be identified in Eq. (1), and used as the design matrix for approximation of data.

Parameters:
  • target_nii (nibabel.spatialimages) – An spatial image object (typically, a Nifti1Image) embedding the target EPI image to be corrected. Provides the location of the N samples (total number of voxels) in the space.

  • ctrl_nii (nibabel.spatialimages) – An spatial image object (typically, a Nifti1Image) embedding the location of the control points of the B-Spline grid. The data array should contain a total of \(K\) knots (control points).

Returns:

weights – A sparse matrix of interpolating weights \(\Psi^3(\mathbf{k}, \mathbf{s})\) for the N voxels of the target EPI, for each of the total K knots. This sparse matrix can be directly used as design matrix for the fitting step of approximation/extrapolation.

Return type:

numpy.ndarray (\(N \times K\))

async sdcflows.transform.unwarp_parallel(fulldataset: ndarray, coordinates: ndarray, fmap_hz: ndarray, pe_info: Sequence[Tuple[int, float]], xfms: Sequence[ndarray], order: int = 3, mode: str = 'constant', cval: float = 0.0, prefilter: bool = True, output_dtype: str | dtype | None = None, max_concurrent: int = 4) ndarray[source]#

Unwarp an EPI dataset parallelizing across volumes.

Parameters:
  • fulldataset (ndarray) – An array of shape (I, J, K, T), where I, J, K are the dimensions of spatial axes and T is the number of volumes. The full data array of the EPI that are wanted after correction.

  • coordinates (ndarray) – An array of shape (3, I, J, K) array providing the voxel (index) coordinates of the reference image (i.e., interpolated points) before SDC/HMC.

  • fmap_hz (ndarray) – An array of shape (I, J, K) containing the displacement of each voxel in voxel units.

  • pe_info (tuple of (int, float)) – A tuple containing the index of the phase-encoding axis in the data array and the readout time (including sign, if displacements must be reversed)

  • xfms (list of obj:~numpy.ndarray) – A list of 4×4 matrices, each one formalizing the estimated head motion alignment to the scan’s reference. Therefore, each of these matrices express the transform of every voxel’s RAS (physical) coordinates in the image used as reference for realignment into the coordinates of each of the EPI series volume.

  • order (int, optional) – The order of the spline interpolation, default is 3. The order has to be in the range 0-5.

  • mode ({‘constant’, ‘reflect’, ‘nearest’, ‘mirror’, ‘wrap’}, optional) – Determines how the input image is extended when the resamplings overflows a border. Default is ‘constant’.

  • cval (float, optional) – Constant value for mode='constant'. Default is 0.0.

  • prefilter (bool, optional) – Determines if the image’s data array is prefiltered with a spline filter before interpolation. The default is True, which will create a temporary float64 array of filtered values if order > 1. If setting this to False, the output will be slightly blurred if order > 1, unless the input is prefiltered, i.e. it is the result of calling the spline filter on the original input.

  • output_dtype (str or dtype) – Override the output data type, instead of propagating it from the moving image.

  • max_concurrent (int) – The maximum number of parallel resamplings at any given time of execution. Use this parameter to set an upper bound to memory utilization.

async sdcflows.transform.worker(data: ndarray, coordinates: ndarray, pe_info: Tuple[int, float], hmc_xfm: ndarray, func: Callable, semaphore: Semaphore) ndarray[source]#

Create one worker and attach it to the execution loop.