fsl.transform.x5
This module contains functions for reading/writing linear/non-linear transformations from/to BIDS X5 files. The following functions are available:
Return the type of the given X5 file - either |
|
Read a linear X5 transformation file from |
|
Write a linear transformation to |
|
Read a nonlinear X5 transformation file from |
|
Write a nonlinear X5 transformation to |
Warning
This is a development release, and is subject to change.
An X5 file is a HDF5 container file which stores a linear or non-linear transformation from one NIfTI image to another.
Several terms may be used to refer to these images, such as source and reference, moving and fixed, from and to, etc. In an X5 file, the two images are simply referred to as A and B, where A refers to the starting point of the transformation, and B to the end point.
X5 files enable a transformation from the world coordinate system of image
A to the world coordinate system of image B. The world
coordinate system of an image is defined by its sform
or qform
(hereafter referred to as the sform
), which is contained in its NIfTI
header.
Custom HDF5 groups
HDF5 files are composed primarily of groups, attributes, and datasets:
Groups are simply containers for attributes, datasets, and other groups.
Datasets are strongly-typed, structured N-dimensional arrays.
Attributes are strongly-typed scalar values associated with a group or dataset.
To simplify the file format definitions below, we shall first define a few custom HDF5 groups. In the file format definitions, a HDF5 group which is listed as being of one of these custom types shall contain the attributes/datasets that are listed here.
affine
A HDF5 group which is listed as being of type affine contains an affine transformation, which can be used to transform coordinates from one space into another. Groups of type affine have the following fields:
Name |
Type |
Value/Description |
|
attribute |
|
|
dataset |
The affine transformation matrix - a |
|
dataset |
Optional pre-calculated inverse |
space
A HDF5 group which is listed as being of type space contains all of the information required to define the space of a NIfTI image, including its shape, dimensions, and voxel-to-world affine transformation.
Groups of type space have the following fields:
Name |
Type |
Value/Description |
|
attribute |
|
|
attribute |
|
|
attribute |
|
|
affine |
The image voxel-to-world transformation (its
|
deformation
A HDF5 group which is listed as being of type deformation contains a non-linear transformation, which can be used to transform coordinates from one space (space A) into another (space B).
The transformation is represented as a 3D deformation field which, at each voxel within the field, may contain:
relative displacements from space A to space B (i.e. for a given location in space A, you can add the displacement values to the coordinates of that location to obtain the coordinates of the corresponding location in space B).
absolute coordinates in space B.
The Mapping
affine can be used to calculate a correspondence between the
deformation field coordinate system and the coordinate system of space A -
it is assumed that space A and the deformation field share a common world
coordinate system.
Groups of type deformation have the following fields:
Name |
Type |
Value/Description |
|
attribute |
|
|
attribute |
|
|
dataset |
The deformation field - a |
|
affine |
The field voxel-to-world transformation (its
|
Linear X5 files
Linear X5 transformation files contain an affine transformation matrix of
shape (4, 4)
, which can be used to transform image A world
coordinates into image B world coordinates.
Linear X5 transformation files are assumed to adhere to the HDF5 structure defined in the table below. All fields are required unless otherwise noted.
Name |
Type |
Value/Description |
Metadata |
||
|
attribute |
|
|
attribute |
|
|
attribute |
JSON string containing unstructured metadata. |
Transformation |
||
|
attribute |
|
|
affine |
Affine transformation from image A world coordinates to image B world coordinates |
|
space |
Image A space |
|
space |
Image B space |
Storage of FSL FLIRT matrices in linear X5 files
FLIRT outputs the result of a linear registration from a source image to a
reference image as an affine matrix of shape (4, 4)
. This matrix encodes a
transformation from source image FSL coordinates to reference image FSL
coordinates [*].
In contrast, X5 matrices encode a transformation in world coordinates,
i.e. they can be used to transform coordinates from the source image to the
reference image, after both images have been transformed into a common
coordinate system via their respective sform
affines.
The fsl.transform.flirt
module contains functions for converting
between FLIRT-style matrices and X5 style matrices.
Non-linear X5 files
Non-linear X5 transformation files contain a non-linear transformation from image A world coordinates to image B world coordinates. The transformation is represented as a 3D deformation field which, at each voxel within the field, may contain:
relative displacements from image A to image B (i.e. for a given location in the image A world coordinate system, add the displacement values to the coordinates to obtain the corresponding location in the image B world coordinate system).
absolute coordinates in the image B world coordinate system.
File format specification
Non-linear X5 transformation files are assumed to adhere to the following HDF5 structure. All fields are required unless otherwise noted.
Name |
Type |
Value/Description |
||
Metadata |
||||
|
attribute |
|
||
|
attribute |
|
||
|
attribute |
JSON string containing unstructured metadata. |
||
Transformation |
||||
Name |
Type |
Value/Description |
||
|
attribute |
|
||
|
deformation |
The deformation field, encoding a nonlinear transformation from image A to image B |
||
|
deformation |
Optional pre-calculated inverse, encoding a nonlinear transformation from image B to image A |
||
|
space |
Image A space |
||
|
space |
Image B space |
Storage of FSL FNIRT warp fields in non-linear X5 files
FLIRT outputs the result of a non-linear registration between a source image and a reference image as either a warp field, or a coefficient field which can be used to generate a warp field. A warp field is defined in terms of the reference image - the warp field has the same shape and FOV as the reference image, and contains either:
relative displacements from the corresponding reference image location to the unwarped source image location
absolute unwarped source image coordinates
The reference image for a FNIRT warp field thus corresponds to image A in a X5 non-linear transform, and the FNIRT source image to image B.
FNIRT warp fields are defined in FSL coordinates - a relative warp contains displacements from reference image FSL coordinates to source image FSL coordinates, and an absolute warp contains source image FSL coordinates.
When a FNIRT warp field is stored in an X5 file, the displacements/coordinates must be adjusted so that they encode a transformation from reference image world coordinates to source image world coordinates.
Conversion of FNIRT coefficient fields
A FNIRT coefficient field can be used to generate a deformation field which contains relative displacements from reference image FSL coordinates to source image FSL coordinates. If an initial affine registration was used as the starting point for FNIRT, this generated displacement field contains relative displacements from the reference image to the aligned source image, i.e. after it has been transformed by the initial affine alignment.
When a FNIRT coefficient field is stored in an X5 file, it must first be converted to a displacement field. The displacements must then be adjusted so that they take into account the initial affine alignment (if relevant), and so that they encode displacements from reference image world coordinates to source image world coordinates.
The fsl.transform.fnirt
module contains functions which can be used to
perform all of the conversions and adjustments required to store FNIRT
transformations as X5 files.
- exception fsl.transform.x5.X5Error[source]
Bases:
Exception
Error raised if an invalid/incompatible file is detected.
- __module__ = 'fsl.transform.x5'
- __weakref__
list of weak references to the object (if defined)
- fsl.transform.x5.inferType(fname)[source]
Return the type of the given X5 file - either
'linear'
or'nonlinear'
.- Parameters:
fname – Name of a X5 file
- Returns:
'linear'
or'nonlinear'
- fsl.transform.x5.writeLinearX5(fname, xform, src, ref)[source]
Write a linear transformation to
fname
.
- fsl.transform.x5.readNonLinearX5(fname)[source]
Read a nonlinear X5 transformation file from
fname
.- Parameters:
fname – File name to read from
- Returns:
- fsl.transform.x5.writeNonLinearX5(fname, field)[source]
Write a nonlinear X5 transformation to
fname
.- Parameters:
fname – File name to write to
field – A
DeformationField
- fsl.transform.x5._readMetadata(group)[source]
Reads a metadata block from the given group, and raises a
X5Error
if it does not look valid.- Parameters:
group – A
h5py.Group
object.- Returns:
A
dict
containing the metadata.
- fsl.transform.x5._writeMetadata(group)[source]
Writes a metadata block into the given group.
- Parameters:
group – A
h5py.Group
object.
- fsl.transform.x5._readAffine(group)[source]
Reads an affine from the given group.
- Parameters:
group – A
h5py.Group
object.- Returns:
numpy
array containing a(4, 4)
affine transformation.
- fsl.transform.x5._writeAffine(group, xform)[source]
Writes the given affine transformation and its inverse to the given group.
- Parameters:
group – A
h5py.Group
object.xform –
numpy
array containing a(4, 4)
affine transformation.
- fsl.transform.x5._readSpace(group)[source]
Reads a space group, defining a source or reference space.
- Parameters:
group – A
h5py.Group
object.- Returns:
Nifti
object. defining the mapping
- fsl.transform.x5._writeSpace(group, img)[source]
Writes a space specified by
img
to the given group.- Parameters:
group – A
h5py.Group
object.img –
Nifti
object. defining the mapping
- fsl.transform.x5._readDeformation(group)[source]
Reads a deformation from the given group.
- Parameters:
group – A
h5py.Group
object- Returns:
A tuple containing
A
numpy.array
containing the deformation fieldA
numpy.array
of shape(4, 4)
containing the voxel to world affine for the deformation fieldThe deformation type - either
'absolute'
or'relative'
- fsl.transform.x5._writeDeformation(group, field)[source]
Write a deformation field to the given group.
- Parameters:
group – A
h5py.Group
objectfield – A
DeformationField
object