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:
Attempt to automatically determine whether a deformation field is specified in absolute or relative coordinates. |
|
Convert a deformation field between storing absolute coordinates or relative displacements. |
|
Adjust the source and/or reference spaces of the given deformation field. |
|
Applies a |
|
Convert a |
-
class
fsl.transform.nonlinear.NonLinearTransform(*args, **kwargs)[source]¶ Bases:
fsl.data.image.ImageClass which represents a nonlinear transformation. This is just a base class for the
DeformationFieldandCoefficientFieldclasses.A nonlinear transformation is an
Imagewhich 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 theNifti.getAffine()method for an overview of the values thatsrcSpaceandrefSpacemay take.- Parameters
image – A string containing the name of an image file to load, or a
numpyarray, or anibabelimage object.src –
Niftirepresenting the source image.ref –
Niftirepresenting the reference image.srcSpace – Coordinate system in the source image that this
NonLinearTransformmaps to. Defaults to'fsl'.refSpace – Coordinate system in the reference image that this
NonLinearTransformmaps from. Defaults to'fsl'.
All other arguments are passed through to
Image.__init__().
-
property
srcSpace¶ Return the source image coordinate system this
NonLinearTransformmaps from - seeNifti.getAffine().
-
property
refSpace¶ Return the reference image coordinate system this
NonLinearTransformmaps to - seeNifti.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 thatfrom_andtomay take.- Parameters
coords – A sequence of XYZ coordinates, or
numpyarray of shape(n, 3)containingnsets of coordinates in the reference space.from – Reference image space that
coordsare defined into – Source image space to transform
coordsinto
- Returns
The corresponding coordinates in the source image space.
-
__module__= 'fsl.transform.nonlinear'¶
-
-
class
fsl.transform.nonlinear.DeformationField(*args, **kwargs)[source]¶ Bases:
fsl.transform.nonlinear.NonLinearTransformClass 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
DeformationFieldis aligned with the reference image in their world coordinate systems (i.e. theirsformaffines 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 thedetectDeformationType()function.
All other arguments are passed through to
NonLinearTransform.__init__().
-
property
deformationType¶ The type of this
DeformationField-'absolute'or'relative'.
-
property
absolute¶ Trueif thisDeformationFieldcontains absolute coordinates.
-
property
relative¶ Trueif thisDeformationFieldcontains 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
numpyarray of shape(n, 3)containingnsets of coordinates in the reference space.from – Reference image space that
coordsare defined into – Source image space to transform
coordsinto
- Returns
coords, transformed into the source image space
-
__module__= 'fsl.transform.nonlinear'¶
-
class
fsl.transform.nonlinear.CoefficientField(*args, **kwargs)[source]¶ Bases:
fsl.transform.nonlinear.NonLinearTransformClass 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
CoefficientFieldas thesrcToRefMatargument to__init__(), and is subsequently accessed via thesrcToRefMat()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
srcSpacecoordinates into reference imagerefSpacecoordinates (typically'fsl'coordinates, i.e. scaled voxels potentially with a left-right flip).
See the
NonLinearTransformclass 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
Noneif there isn’t one.
-
property
refToSrcMat¶ Return the inverse of the initial source-to-reference affine, or
Noneif there isn’t one.
-
asDeformationField(defType='relative', premat=True)[source]¶ Convert this
CoefficientFieldto aDeformationField.This method is a wrapper around
coefficientFieldToDeformationField()
-
transform(coords, from_=None, to=None, premat=True)[source]¶ Overrides
NonLinearTransform.transform(). Transforms the givencoordsfrom the reference image space into the source image space.- Parameters
coords – A sequence of XYZ coordinates, or
numpyarray of shape(n, 3)containingnsets of coordinates in the reference space.from – Reference image space that
coordsare defined into – Source image space to transform
coordsintopremat – If
True, the inversesrcToRefMat()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 forcoords
-
__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 likefieldcontains 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
DeformationFieldinstancedefType – Either
'absolute'or'relative'. If not provided, the opposite type tofield.deformationTypeis used.
- Returns
A
numpy.arraycontaining 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 thefrom_andtoarguments.- Parameters
field – A
DeformationFieldinstancefrom – New reference image coordinate system
to – New source image coordinate system
- Returns
A new
DeformationFieldwhich transforms between the referencefrom_coordinate system and the sourcetocoordinate system.
-
fsl.transform.nonlinear.applyDeformation(image, field, ref=None, order=1, mode=None, cval=None, premat=None)[source]¶ Applies a
DeformationFieldto anImage.The image is transformed into the space of the field’s reference image space. See the
scipy.ndimage.interpolation.map_coordinatesfunction for details on theorder,modeandcvaloptions.If an alternate reference image is provided via the
refargument, the deformation field is resampled into its space, and then applied to the input image. It is therefore assumed that an alternaterefis aligned in world coordinates with the field’s actual reference image.- Parameters
image –
Imageto be transformedfield –
DeformationFieldto useref – Alternate reference image - if not provided,
field.refis usedorder – Spline interpolation order, passed through to the
scipy.ndimage.affine_transformfunction -0corresponds to nearest neighbour interpolation,1(the default) to linear interpolation, and3to 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
imageis not in the same space asfield.src. Assumed to transform fromimagevoxel coordinates intofield.srcvoxel coordinates.
- Returns
numpy.arraycontaining the transformed image data.
-
fsl.transform.nonlinear.coefficientFieldToDeformationField(field, defType='relative', premat=True)[source]¶ Convert a
CoefficientFieldinto aDeformationField.- Parameters
field –
CoefficientFieldto convertdefType – The type of deformation field - either
'relative'(the default) or'absolute'.premat – If
True(the default), thesrcToRefMat()is encoded into the deformation field.
- Returns
DeformationFieldcalculated fromfield.