sdcflows.interfaces.bspline module#

Filtering of \(B_0\) field mappings with B-Splines.

class sdcflows.interfaces.bspline.ApplyCoeffsField(from_file=None, resource_monitor=None, **inputs)[source]#

Bases: SimpleInterface

Undistort a target, distorted image with a fieldmap, formalized by its B-Spline coefficients.

Preconditions:

  • We have a “target EPI” - a BOLD or DWI dataset (or even MPRAGE, same principle), without having gone through HMC or SDC.

  • We have also the list of HMC matrices that accounts for head-motion, so after resampling the dataset through this list of transforms the head does not move anymore.

  • We have estimated the fieldmap’s coefficients

  • We have the “fieldmap-to-data” affine transform that aligns the target dataset (e.g., EPI) and the fieldmap’s “magnitude” (phasediff et al.) or “reference” (pepolar, syn).

The algorithm is implemented in the B0FieldTransform object. First, we will call fit, which results in:

  1. The reference grid of the target dataset is projected onto the fieldmap space

  2. The B-Spline coefficients are applied to reconstruct the field on the grid resulting above.

After which, we can then call apply. This second step will:

  1. Find the location of every voxel on each timepoint (meaning, after the head moved) and progress (or recede) along the phase-encoding axis to find the actual (voxel) coordinates of each voxel. With those coordinates known, interpolation is trivial.

  2. Generate a spatial image with the new data.

Mandatory Inputs:
  • in_coeff (a list of items which are a pathlike object or string representing an existing file) – Input coefficients as calculated in the estimation stage.

  • in_data (a pathlike object or string representing a file) – Input EPI data to be corrected.

  • pe_dir (a list of items which are ‘i’ or ‘i-’ or ‘j’ or ‘j-’ or ‘k’ or ‘k-’) – The phase-encoding direction corresponding to in_data.

  • ro_time (a list of items which are a float) – EPI readout time (s).

Optional Inputs:
  • approx (a boolean) – Reconstruct the fieldmap on its original grid and then interpolate on the rotated grid, rather than reconstructing directly on the rotated grid. (Nipype default value: True)

  • data2fmap_xfm (a list of items which are a pathlike object or string representing an existing file) – The transform by which the fieldmap can be resampled on the target EPI’s grid. Mutually exclusive with inputs: f, m, a, p, 2, d, a, t, a, _, x, f, m.

  • fmap2data_xfm (a list of items which are a pathlike object or string representing an existing file) – The transform by which the target EPI can be resampled on the fieldmap’s grid. Mutually exclusive with inputs: d, a, t, a, 2, f, m, a, p, _, x, f, m.

  • in_xfms (a list of items which are a list of items which are a list of items which are a float) – List of head-motion correction matrices.

  • num_threads (an integer) – Number of threads.

Outputs:
  • out_corrected (a list of items which are a pathlike object or string representing an existing file)

  • out_field (a list of items which are a pathlike object or string representing an existing file)

class sdcflows.interfaces.bspline.BSplineApprox(from_file=None, resource_monitor=None, **inputs)[source]#

Bases: SimpleInterface

Approximate the \(B_0\) field using tensor-product B-Splines.

The approximation effectively smooths the data, removing spikes and other sources of noise, as well as enables the extrapolation of the \(B_0\) field beyond the brain mask, which alleviates boundary effects in correction.

This interface resolves the optimization problem of obtaining the B-Spline coefficients \(c(\mathbf{k})\) that best approximate the data samples within the brain mask \(f(\mathbf{s})\), following Eq. (17) – in that case for 2D – of [Unser1999]. Here, and for the case of 3D, the formulism is adapted in Eq. (1) of the transform module.

References

[Unser1999]

M. Unser, “Splines: A Perfect Fit for Signal and Image Processing,” IEEE Signal Processing Magazine 16(6):22-38, 1999.

See also

grid_bspline_weights() - for the evaluation of the tensor-product, cubic B-Splines (\(\Psi^3(\mathbf{k}, \mathbf{s})\)) formalized in Eq. (2) of the transform module.

Mandatory Inputs:

in_data (a pathlike object or string representing an existing file) – Path to a fieldmap.

Optional Inputs:
  • bs_spacing (a list of items which are a tuple of the form: (a float, a float, a float)) – Spacing between B-Spline control points. (Nipype default value: [(40.0, 40.0, 20.0)])

  • debug (a boolean) – Generate extra assets for debugging. (Nipype default value: False)

  • extrapolate (a boolean) – Generate a field, extrapolated outside the brain mask. (Nipype default value: True)

  • in_mask (a pathlike object or string representing an existing file) – Path to a brain mask.

  • recenter (‘mode’ or ‘median’ or ‘mean’ or False) – Strategy to recenter the distribution of the input fieldmap. (Nipype default value: mode)

  • ridge_alpha (a float) – Controls the regularization. (Nipype default value: 0.01)

  • zooms_min (a float or a legal value) – Limit minimum image zooms, set 0.0 to use the original image. (Nipype default value: 4.0)

Outputs:
  • out_coeff (a list of items which are a pathlike object or string representing an existing file)

  • out_error (a pathlike object or string representing an existing file)

  • out_extrapolated (a pathlike object or string representing a file)

  • out_field (a pathlike object or string representing an existing file)

class sdcflows.interfaces.bspline.TOPUPCoeffReorient(from_file=None, resource_monitor=None, **inputs)[source]#

Bases: SimpleInterface

Revise the orientation of TOPUP-generated B-Spline coefficients.

TOPUP-generated “fieldcoeff” files are just B-Spline fields, where the shape of the field is fixated to be a decimated grid of the original image by an integer factor and added 3 pixels on each dimension. This is one root reason why TOPUP errors (FSL 6) or segfaults (FSL 5), when the input image has odd number of voxels along one or more directions.

These “fieldcoeff” are fixated to be zero-centered, and have “plumb” orientation (as in, aligned with cardinal/imaging axes). The q-form of these NIfTI files is always diagonal, with the decimation factors set on the diagonal (and hence, the voxel zooms). The origin of the q-form is set to the reference image’s shape. All the director cosines of the output coefficients will be positive. In other words, the output orientation is either RAS, ARS, ASR, SAR, or SRA.

This interface modifies these coefficient files to be fully-fledged NIfTI images aligned with the reference image. Therefore, the s-form header of the coefficients file is updated to match that of the reference file. The s-form header is used because the imaging axes may be oblique.

The q-form retains the original header and is marked with code 0.

Mandatory Inputs:
  • fmap_ref (a pathlike object or string representing an existing file) – The fieldmap reference.

  • in_coeff (a list of items which are a pathlike object or string representing a file) – Input coefficients file(s) from TOPUP.

  • pe_dir (‘i’ or ‘i-’ or ‘j’ or ‘j-’ or ‘k’ or ‘k-’ or ‘x’ or ‘x-’ or ‘y’ or ‘y-’ or ‘z’ or ‘z-’) – Phase encoding direction.

Outputs:

out_coeff (a list of items which are a pathlike object or string representing an existing file) – Patched coefficients.

class sdcflows.interfaces.bspline.TransformCoefficients(from_file=None, resource_monitor=None, **inputs)[source]#

Bases: SimpleInterface

Project coefficients files to another space through a rigid-body transform.

Mandatory Inputs:
  • fmap_ref (a pathlike object or string representing an existing file) – The fieldmap reference.

  • in_coeff (a list of items which are a pathlike object or string representing a file) – Input coefficients file(s).

  • transform (a pathlike object or string representing an existing file) – Rigid-body transform file.

Optional Inputs:

fmap_target (a pathlike object or string representing an existing file) – The distorted EPI target (feed to set debug mode on).

Outputs:

out_coeff (a list of items which are a pathlike object or string representing an existing file) – Moved coefficients.

sdcflows.interfaces.bspline.bspline_grid(img, control_zooms_mm=(40.0, 40.0, 20.0))[source]#

Create a Nifti1Image embedding the location of control points.