ants.registration package

Submodules

ants.registration.affine_initializer module

ants.registration.affine_initializer.affine_initializer(fixed_image, moving_image, search_factor=20, radian_fraction=0.1, use_principal_axis=False, local_search_iterations=10, mask=None, txfn=None)[source]

A multi-start optimizer for affine registration Searches over the sphere to find a good initialization for further registration refinement, if needed. This is a wrapper for the ANTs function antsAffineInitializer.

ANTsR function: affineInitializer

Parameters:
  • fixed_image (ANTsImage) – the fixed reference image

  • moving_image (ANTsImage) – the moving image to be mapped to the fixed space

  • search_factor (scalar) – degree of increments on the sphere to search

  • radian_fraction (scalar) – between zero and one, defines the arc to search over

  • use_principal_axis (boolean) – boolean to initialize by principal axis

  • local_search_iterations (scalar) – gradient descent iterations

  • mask (ANTsImage (optional)) – optional mask to restrict registration

  • txfn (string (optional)) – filename for the transformation

Returns:

transformation matrix

Return type:

ndarray

Example

>>> import ants
>>> fi = ants.image_read(ants.get_ants_data('r16'))
>>> mi = ants.image_read(ants.get_ants_data('r27'))
>>> txfile = ants.affine_initializer( fi, mi )
>>> tx = ants.read_transform(txfile, dimension=2)

ants.registration.apply_transforms module

ants.registration.apply_transforms.apply_transforms(fixed, moving, transformlist, interpolator='linear', imagetype=0, whichtoinvert=None, compose=None, defaultvalue=0, singleprecision=False, verbose=False, **kwargs)[source]

Apply a transform list to map an image from one domain to another. In image registration, one computes mappings between (usually) pairs of images. These transforms are often a sequence of increasingly complex maps, e.g. from translation, to rigid, to affine to deformation. The list of such transforms is passed to this function to interpolate one image domain into the next image domain, as below. The order matters strongly and the user is advised to familiarize with the standards established in examples.

ANTsR function: antsApplyTransforms

Parameters:
  • fixed (ANTsImage) – fixed image defining domain into which the moving image is transformed. The output will have the same pixel type as this image.

  • moving (AntsImage) – moving image to be mapped to fixed space.

  • transformlist (list of strings) – list of transforms generated by ants.registration where each transform is a filename.

  • interpolator (string) –

    Choice of interpolator. Supports partial matching.

    linear nearestNeighbor multiLabel for label images (deprecated, prefer genericLabel) gaussian bSpline cosineWindowedSinc welchWindowedSinc hammingWindowedSinc lanczosWindowedSinc genericLabel use this for label images

  • imagetype (python:integer) – choose 0/1/2/3 mapping to scalar/vector/tensor/time-series

  • whichtoinvert (list of booleans (optional)) – Must be same length as transformlist. whichtoinvert[i] is True if transformlist[i] is a matrix, and the matrix should be inverted. If transformlist[i] is a warp field, whichtoinvert[i] must be False. If the transform list is a matrix followed by a warp field, whichtoinvert defaults to (True,False). Otherwise it defaults to [False]*len(transformlist)).

  • compose (string (optional)) – if it is a string pointing to a valid file location, this will force the function to return a composite transformation filename.

  • defaultvalue (scalar) – Default voxel value for mappings outside the image domain.

  • singleprecision (boolean) – if True, use float32 for computations. This is useful for reducing memory usage for large datasets, at the cost of precision.

  • verbose (boolean) – print command and run verbose application of transform.

  • kwargs (keyword arguments) – extra parameters

Return type:

ANTsImage or string (transformation filename)

Example

>>> import ants
>>> fixed = ants.image_read( ants.get_ants_data('r16') )
>>> moving = ants.image_read( ants.get_ants_data('r64') )
>>> fixed = ants.resample_image(fixed, (64,64), 1, 0)
>>> moving = ants.resample_image(moving, (64,64), 1, 0)
>>> mytx = ants.registration(fixed=fixed , moving=moving ,
                             type_of_transform = 'SyN' )
>>> mywarpedimage = ants.apply_transforms( fixed=fixed, moving=moving,
                                           transformlist=mytx['fwdtransforms'] )
ants.registration.apply_transforms.apply_transforms_to_points(dim, points, transformlist, whichtoinvert=None, verbose=False)[source]

Apply a transform list to map a pointset from one domain to another. In registration, one computes mappings between pairs of domains. These transforms are often a sequence of increasingly complex maps, e.g. from translation, to rigid, to affine to deformation. The list of such transforms is passed to this function to interpolate one image domain into the next image domain, as below. The order matters strongly and the user is advised to familiarize with the standards established in examples. Importantly, point mapping goes the opposite direction of image mapping, for both reasons of convention and engineering.

ANTsR function: antsApplyTransformsToPoints

Parameters:
  • dim (python:integer) – dimensionality of the transformation.

  • points (data frame) – moving point set with n-points in rows of at least dim columns - we maintain extra information in additional columns. this should be a data frame with columns names x, y, z, t.

  • transformlist (list of strings) – list of transforms generated by ants.registration where each transform is a filename.

  • whichtoinvert (list of booleans (optional)) – Must be same length as transformlist. whichtoinvert[i] is True if transformlist[i] is a matrix, and the matrix should be inverted. If transformlist[i] is a warp field, whichtoinvert[i] must be False. If the transform list is a matrix followed by a warp field, whichtoinvert defaults to (True,False). Otherwise it defaults to [False]*len(transformlist)).

  • verbose (boolean) –

Return type:

data frame of transformed points

Example

>>> import ants
>>> fixed = ants.image_read( ants.get_ants_data('r16') )
>>> moving = ants.image_read( ants.get_ants_data('r27') )
>>> reg = ants.registration( fixed, moving, 'Affine' )
>>> d = {'x': [128, 127], 'y': [101, 111]}
>>> pts = pd.DataFrame(data=d)
>>> ptsw = ants.apply_transforms_to_points( 2, pts, reg['fwdtransforms'])

ants.registration.build_template module

ants.registration.build_template.build_template(initial_template=None, image_list=None, iterations=3, gradient_step=0.2, blending_weight=0.75, weights=None, useNoRigid=True, output_dir=None, **kwargs)[source]

Estimate an optimal template from an input image_list

ANTsR function: N/A

Parameters:
  • initial_template (ANTsImage) – initialization for the template building

  • image_list (ANTsImages) – images from which to estimate template

  • iterations (python:integer) – number of template building iterations

  • gradient_step (scalar) – for shape update gradient

  • blending_weight (scalar) – weight for image blending

  • weights (vector) – weight for each input image

  • useNoRigid (boolean) – equivalent of -y in the script. Template update step will not use the rigid component if this is True.

  • output_dir (path) – directory name where intermediate transforms are written

  • kwargs (keyword args) – extra arguments passed to ants registration

Return type:

ANTsImage

Example

>>> import ants
>>> image = ants.image_read( ants.get_ants_data('r16') )
>>> image2 = ants.image_read( ants.get_ants_data('r27') )
>>> image3 = ants.image_read( ants.get_ants_data('r85') )
>>> timage = ants.build_template( image_list = ( image, image2, image3 ) ).resample_image( (45,45))
>>> timagew = ants.build_template( image_list = ( image, image2, image3 ), weights = (5,1,1) )

ants.registration.create_jacobian_determinant_image module

ants.registration.create_jacobian_determinant_image.create_jacobian_determinant_image(domain_image, tx, do_log=False, geom=False)[source]

Compute the jacobian determinant from a transformation file

ANTsR function: createJacobianDeterminantImage

Parameters:
  • domain_image (ANTsImage) – image that defines transformation domain

  • tx (string) – deformation transformation file name

  • do_log (boolean) – return the log jacobian

  • geom (boolean) – use the geometric jacobian calculation (boolean)

Return type:

ANTsImage

Example

>>> import ants
>>> fi = ants.image_read( ants.get_ants_data('r16'))
>>> mi = ants.image_read( ants.get_ants_data('r64'))
>>> fi = ants.resample_image(fi,(128,128),1,0)
>>> mi = ants.resample_image(mi,(128,128),1,0)
>>> mytx = ants.registration(fixed=fi , moving=mi, type_of_transform = ('SyN') )
>>> jac = ants.create_jacobian_determinant_image(fi,mytx['fwdtransforms'][0],1)
ants.registration.create_jacobian_determinant_image.deformation_gradient(warp_image, to_rotation=False, to_inverse_rotation=False, py_based=False)[source]

Compute the deformation gradient from an image containing a warp (deformation).

This function now includes a highly optimized pure Python/NumPy implementation.

ANTsR function: NA

Parameters:
  • warp_image (ANTsImage (or filename if not py_based)) – image that defines the deformation field (vector pixels)

  • to_rotation (boolean) – maps deformation gradient to a rotation matrix using polar decomposition.

  • to_inverse_rotation (boolean) – map the deformation gradient to a rotation matrix, and return its inverse. This is useful for reorienting tensors and vectors after resampling.

  • py_based (boolean) – If True, uses the optimized pure Python/NumPy implementation. If False, uses the classic C++ backend.

Returns:

where U is the x-component of deformation and xyz are spatial.

Return type:

ANTsImage with dimension*dimension components indexed in order U_xyz, V_xyz, W_xyz

Example

>>> import ants
>>> fi = ants.image_read( ants.get_ants_data('r16'))
>>> mi = ants.image_read( ants.get_ants_data('r64'))
>>> fi = ants.resample_image(fi,(128,128),1,0)
>>> mi = ants.resample_image(mi,(128,128),1,0)
>>> mytx = ants.registration(fixed=fi , moving=mi, type_of_transform = ('SyN') )
>>> warp = ants.image_read( mytx['fwdtransforms'][0] )
>>> # Use the fast, optimized Python implementation
>>> dg_py = ants.deformation_gradient( warp, py_based=True )
>>> dg_rot_py = ants.deformation_gradient( warp, to_rotation=True, py_based=True )

ants.registration.create_warped_grid module

ants.registration.create_warped_grid.create_warped_grid(image, grid_step=10, grid_width=2, grid_directions=(True, True), fixed_reference_image=None, transform=None, foreground=1, background=0)[source]

Deforming a grid is a helpful way to visualize a deformation field. This function enables a user to define the grid parameters and apply a deformable map to that grid.

ANTsR function: createWarpedGrid

Parameters:
  • image (ANTsImage) – input image

  • grid_step (scalar) – width of grid blocks

  • grid_width (scalar) – width of grid lines

  • grid_directions (tuple of booleans) – directions in which to draw grid lines, boolean vector

  • fixed_reference_image (ANTsImage (optional)) – reference image space

  • transform (list/tuple of strings (optional)) – vector of transforms

  • foreground (scalar) – intensity value for grid blocks

  • background (scalar) – intensity value for grid lines

Return type:

ANTsImage

Example

>>> import ants
>>> fi = ants.image_read( ants.get_ants_data( 'r16' ) )
>>> mi = ants.image_read( ants.get_ants_data( 'r64' ) )
>>> mygr = ants.create_warped_grid( mi )
>>> mytx = ants.registration(fixed=fi, moving=mi, type_of_transform = ('SyN') )
>>> mywarpedgrid = ants.create_warped_grid( mygr, grid_directions=(False,True),
                        transform=mytx['fwdtransforms'], fixed_reference_image=fi )

ants.registration.fsl2antstransform module

ants.registration.interface module

ants.registration.make_points_image module

ants.registration.metrics module

ants.registration.reflect_image module

ants.registration.reorient_image module

ants.registration.resample_image module

ants.registration.symmetrize_image module

Module contents