fsl.transform.nonlinear

This module contains data structures and functions for working with FNIRT-style nonlinear transformations.

The DeformationField and CoefficientField can be used to load and interact with FNIRT-style transformation images. The following utility functions are also available:

detectDeformationType

Attempt to automatically determine whether a deformation field is specified in absolute or relative coordinates.

convertDeformationType

Convert a deformation field between storing absolute coordinates or relative displacements.

convertDeformationSpace

Adjust the source and/or reference spaces of the given deformation field.

applyDeformation

Applies a DeformationField to an Image.

coefficientFieldToDeformationField

Convert a CoefficientField into a DeformationField.

class fsl.transform.nonlinear.NonLinearTransform(*args, **kwargs)[source]

Bases: Image

Class which represents a nonlinear transformation. This is just a base class for the DeformationField and CoefficientField classes.

A nonlinear transformation is an Image which contains some mapping between a source image coordinate system and a reference image coordinate system.

In FSL, non-linear transformations are defined in terms of the reference image coordinate system. At a given location in the reference image space, the non-linear mapping at that location can be used to calculate the corresponding location in the source image space. Therefore, these non-linear transformation effectively encode a transformation from the reference image to the (unwarped) source image.

__init__(image, src, ref, srcSpace=None, refSpace=None, **kwargs)[source]

Create a NonLinearTransform. See the Nifti.getAffine() method for an overview of the values that srcSpace and refSpace may take.

Parameters:
  • image – A string containing the name of an image file to load, or a numpy array, or a nibabel image object.

  • srcNifti representing the source image.

  • refNifti representing the reference image.

  • srcSpace – Coordinate system in the source image that this NonLinearTransform maps to. Defaults to 'fsl'.

  • refSpace – Coordinate system in the reference image that this NonLinearTransform maps from. Defaults to 'fsl'.

All other arguments are passed through to Image.__init__().

property src

Return a reference to the Nifti instance representing the source image.

property ref

Return a reference to the Nifti instance representing the reference image.

property srcSpace

Return the source image coordinate system this NonLinearTransform maps from - see Nifti.getAffine().

property refSpace

Return the reference image coordinate system this NonLinearTransform maps to - see Nifti.getAffine().

transform(coords, from_=None, to=None)[source]

Transform coordinates from the reference image space to the source image space. Implemented by sub-classes.

See the Nifti.getAffine() method for an overview of the values that from_ and to may take.

Parameters:
  • coords – A sequence of XYZ coordinates, or numpy array of shape (n, 3) containing n sets of coordinates in the reference space.

  • from – Reference image space that coords are defined in

  • to – Source image space to transform coords into

Returns:

The corresponding coordinates in the source image space.

__annotations__ = {}
__module__ = 'fsl.transform.nonlinear'
class fsl.transform.nonlinear.DeformationField(*args, **kwargs)[source]

Bases: NonLinearTransform

Class which represents a deformation (a.k.a. warp) field which, at each voxel, contains either:

  • a relative displacement from the reference image space to the source image space, or

  • absolute coordinates in the source space

It is assumed that the a DeformationField is aligned with the reference image in their world coordinate systems (i.e. their sform affines project the reference image and the deformation field into alignment).

__init__(image, src, ref=None, **kwargs)[source]

Create a DisplacementField.

Parameters:
  • ref – Optional. If not provided, it is assumed that the reference is defined in the same space as image.

  • defType – Either 'absolute' or 'relative', indicating the type of this displacement field. If not provided, will be inferred via the detectDeformationType() function.

All other arguments are passed through to NonLinearTransform.__init__().

property deformationType

The type of this DeformationField - 'absolute' or 'relative'.

property absolute

True if this DeformationField contains absolute coordinates.

property relative

True if this DeformationField contains relative displacements.

transform(coords, from_=None, to=None)[source]

Transform the given XYZ coordinates from the reference image space to the source image space.

Parameters:
  • coords – A sequence of XYZ coordinates, or numpy array of shape (n, 3) containing n sets of coordinates in the reference space.

  • from – Reference image space that coords are defined in

  • to – Source image space to transform coords into

Returns:

coords, transformed into the source image space

__annotations__ = {}
__module__ = 'fsl.transform.nonlinear'
class fsl.transform.nonlinear.CoefficientField(*args, **kwargs)[source]

Bases: NonLinearTransform

Class which represents a B-spline coefficient field generated by FNIRT.

The displacements() method can be used to calculate relative displacements for a set of reference space voxel coordinates.

A FNIRT nonlinear transformation often contains a premat, a global affine transformation from the source space to the reference space, which was calculated with FLIRT, and used as the starting point for the non-linear optimisation performed by FNIRT.

This affine may be provided when creating a CoefficientField as the srcToRefMat argument to __init__(), and is subsequently accessed via the srcToRefMat() attribute.

__init__(image, src, ref, srcSpace=None, refSpace=None, fieldType='cubic', knotSpacing=None, fieldToRefMat=None, srcToRefMat=None, **kwargs)[source]

Create a CoefficientField.

Parameters:
  • fieldType – Must currently be 'cubic'

  • knotSpacing – A tuple containing the spline knot spacings along each axis.

  • fieldToRefMat – Affine transformation which can transform reference image voxel coordinates into coefficient field voxel coordinates.

  • srcToRefMat – Optional initial global affine transformation from the source image to the reference image. This is assumed to be a FLIRT-style matrix, i.e. it transforms from source image srcSpace coordinates into reference image refSpace coordinates (typically 'fsl' coordinates, i.e. scaled voxels potentially with a left-right flip).

See the NonLinearTransform class for details on the other arguments.

property fieldType

Return the type of the coefficient field, currently always 'cubic'.

property knotSpacing

Return a tuple containing spline knot spacings along the x, y, and z axes.

property fieldToRefMat

Return an affine transformation which can transform coefficient field voxel coordinates into reference image voxel coordinates.

property refToFieldMat

Return an affine transformation which can transform reference image voxel coordinates into coefficient field voxel coordinates.

property srcToRefMat

Return the initial source-to-reference affine, or None if there isn’t one.

property refToSrcMat

Return the inverse of the initial source-to-reference affine, or None if there isn’t one.

asDeformationField(defType='relative', premat=True)[source]

Convert this CoefficientField to a DeformationField.

This method is a wrapper around coefficientFieldToDeformationField()

transform(coords, from_=None, to=None, premat=True)[source]

Overrides NonLinearTransform.transform(). Transforms the given coords from the reference image space into the source image space.

Parameters:
  • coords – A sequence of XYZ coordinates, or numpy array of shape (n, 3) containing n sets of coordinates in the reference space.

  • from – Reference image space that coords are defined in

  • to – Source image space to transform coords into

  • premat – If True, the inverse srcToRefMat() is applied to the coordinates after the displacements have been added.

Returns:

coords, transformed into the source image space

displacements(coords)[source]

Calculate the relative displacements for the given coordinates.

Parameters:

coords(N, 3) array of reference image voxel coordinates.

Returns:

A (N, 3) array of relative displacements to the source image for coords

__annotations__ = {}
__module__ = 'fsl.transform.nonlinear'
fsl.transform.nonlinear.detectDeformationType(field)[source]

Attempt to automatically determine whether a deformation field is specified in absolute or relative coordinates.

Parameters:

field – A DeformationField

Returns:

'absolute' if it looks like field contains absolute coordinates, 'relative' otherwise.

fsl.transform.nonlinear.convertDeformationType(field, defType=None)[source]

Convert a deformation field between storing absolute coordinates or relative displacements.

Parameters:
  • field – A DeformationField instance

  • defType – Either 'absolute' or 'relative'. If not provided, the opposite type to field.deformationType is used.

Returns:

A numpy.array containing the adjusted deformation field.

fsl.transform.nonlinear.convertDeformationSpace(field, from_, to)[source]

Adjust the source and/or reference spaces of the given deformation field. See the Nifti.getAffine() method for the valid values for the from_ and to arguments.

Parameters:
  • field – A DeformationField instance

  • from – New reference image coordinate system

  • to – New source image coordinate system

Returns:

A new DeformationField which transforms between the reference from_ coordinate system and the source to coordinate system.

fsl.transform.nonlinear.applyDeformation(image, field, ref=None, order=1, mode=None, cval=None, premat=None)[source]

Applies a DeformationField to an Image.

The image is transformed into the space of the field’s reference image space. See the scipy.ndimage.interpolation.map_coordinates function for details on the order, mode and cval options.

If an alternate reference image is provided via the ref argument, the deformation field is resampled into its space, and then applied to the input image. It is therefore assumed that an alternate ref is aligned in world coordinates with the field’s actual reference image.

Parameters:
  • imageImage to be transformed

  • fieldDeformationField to use

  • ref – Alternate reference image - if not provided, field.ref is used

  • order – Spline interpolation order, passed through to the scipy.ndimage.affine_transform function - 0 corresponds to nearest neighbour interpolation, 1 (the default) to linear interpolation, and 3 to cubic interpolation.

  • mode – How to handle regions which are outside of the image FOV. Defaults to ‘’nearest’`.

  • cval – Constant value to use when mode='constant'.

  • premat – Optional affine transform which can be used if image is not in the same space as field.src. Assumed to transform from image voxel coordinates into field.src voxel coordinates.

Returns:

numpy.array containing the transformed image data.

fsl.transform.nonlinear.coefficientFieldToDeformationField(field, defType='relative', premat=True)[source]

Convert a CoefficientField into a DeformationField.

Parameters:
  • fieldCoefficientField to convert

  • defType – The type of deformation field - either 'relative' (the default) or 'absolute'.

  • premat – If True (the default), the srcToRefMat() is encoded into the deformation field.

Returns:

DeformationField calculated from field.