Linear Transforms

Linear transforms.

class nitransforms.linear.Affine(matrix=None, reference=None)[source]

Represents linear transforms on image data.

classmethod from_filename(filename, fmt=None, reference=None, moving=None)[source]

Create an affine from a transform file.

classmethod from_matvec(mat=None, vec=None, reference=None)[source]

Create an affine from a matrix and translation pair.

Example

>>> Affine.from_matvec(vec=(4, 0, 0))  
array([[1., 0., 0., 4.],
       [0., 1., 0., 0.],
       [0., 0., 1., 0.],
       [0., 0., 0., 1.]])
map(x, inverse=False)[source]

Apply \(y = f(x)\).

Parameters
  • x (N x D numpy.ndarray) – Input RAS+ coordinates (i.e., physical coordinates).

  • inverse (bool) – If True, apply the inverse transform \(x = f^{-1}(y)\).

Returns

y – Transformed (mapped) RAS+ coordinates (i.e., physical coordinates).

Return type

N x D numpy.ndarray

Examples

>>> xfm = Affine([[1, 0, 0, 1], [0, 1, 0, 2], [0, 0, 1, 3], [0, 0, 0, 1]])
>>> xfm.map((0,0,0))
array([[1., 2., 3.]])
>>> xfm.map((0,0,0), inverse=True)
array([[-1., -2., -3.]])
property matrix

Access the internal representation of this affine.

property ndim

Access the internal representation of this affine.

to_filename(filename, fmt='X5', moving=None)[source]

Store the transform in the requested output format.

class nitransforms.linear.LinearTransformsMapping(transforms, reference=None)[source]

Represents a series of linear transforms.

apply(spatialimage, reference=None, order=3, mode='constant', cval=0.0, prefilter=True, output_dtype=None)[source]

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

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

  • reference (spatial object, optional) – The image, surface, or combination thereof containing the coordinates of samples that will be sampled.

  • 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.

Returns

resampled – The data imaged after resampling to reference space.

Return type

spatialimage or ndarray

map(x, inverse=False)[source]

Apply \(y = f(x)\).

Parameters
  • x (N x D numpy.ndarray) – Input RAS+ coordinates (i.e., physical coordinates).

  • inverse (bool) – If True, apply the inverse transform \(x = f^{-1}(y)\).

Returns

y – Transformed (mapped) RAS+ coordinates (i.e., physical coordinates).

Return type

N x D numpy.ndarray

Examples

>>> xfm = LinearTransformsMapping([
...     [[1., 0, 0, 1.], [0, 1., 0, 2.], [0, 0, 1., 3.], [0, 0, 0, 1.]],
...     [[1., 0, 0, -1.], [0, 1., 0, -2.], [0, 0, 1., -3.], [0, 0, 0, 1.]],
... ])
>>> xfm.matrix
array([[[ 1.,  0.,  0.,  1.],
        [ 0.,  1.,  0.,  2.],
        [ 0.,  0.,  1.,  3.],
        [ 0.,  0.,  0.,  1.]],

       [[ 1.,  0.,  0., -1.],
        [ 0.,  1.,  0., -2.],
        [ 0.,  0.,  1., -3.],
        [ 0.,  0.,  0.,  1.]]])
>>> y = xfm.map([(0, 0, 0), (-1, -1, -1), (1, 1, 1)])
>>> y[0, :, :3]
array([[1., 2., 3.],
       [0., 1., 2.],
       [2., 3., 4.]])
>>> y = xfm.map([(0, 0, 0), (-1, -1, -1), (1, 1, 1)], inverse=True)
>>> y[0, :, :3]
array([[-1., -2., -3.],
       [-2., -3., -4.],
       [ 0., -1., -2.]])
to_filename(filename, fmt='X5', moving=None)[source]

Store the transform in the requested output format.

nitransforms.linear.load(filename, fmt=None, reference=None, moving=None)[source]

Load a linear transform file.

Examples

>>> xfm = load(regress_dir / "affine-LAS.itk.tfm")
>>> isinstance(xfm, Affine)
True
>>> xfm = load(regress_dir / "itktflist.tfm")
>>> isinstance(xfm, LinearTransformsMapping)
True