Title: | Image Registration Using the 'NiftyReg' Library |
---|---|
Description: | Provides an 'R' interface to the 'NiftyReg' image registration tools <https://github.com/KCL-BMEIS/niftyreg>. Linear and nonlinear registration are supported, in two and three dimensions. |
Authors: | Jon Clayden [cre, aut] , Marc Modat [aut], Benoit Presles [aut], Thanasis Anthopoulos [aut], Pankaj Daga [aut] |
Maintainer: | Jon Clayden <[email protected]> |
License: | GPL-2 |
Version: | 2.8.4 |
Built: | 2024-11-18 03:23:04 UTC |
Source: | https://github.com/jonclayden/rniftyreg |
This function allows a precomputed transformation to be applied to a new image or set of points.
applyTransform(transform, x, interpolation = 3L, nearest = FALSE, internal = FALSE)
applyTransform(transform, x, interpolation = 3L, nearest = FALSE, internal = FALSE)
transform |
|
x |
A numeric vector, representing a pixel/voxel location in source space, or a matrix with rows representing such points, or an image with the same dimensions as the original source image. |
interpolation |
A single integer specifying the type of interpolation to be applied to the final resampled image. May be 0 (nearest neighbour), 1 (trilinear) or 3 (cubic spline). No other values are valid. |
nearest |
Logical value: if |
internal |
If |
Points may be transformed from source to target space exactly under an affine transformation, but nonlinear transformation is inexact. Its accuracy will depend to some extent on the density of the control point grid and the geometry of the deformation in the vicinity of the points of interest. Nevertheless, it should be quite sufficient for most purposes.
The method is to first convert the control points to a deformation field
(cf. deformationField
), which encodes the location of each
target space voxel in the source space. The target voxel closest to the
requested location is found by searching through this deformation field, and
returned if nearest
is TRUE
or it coincides exactly with the
requested location. Otherwise, a block of four voxels in each dimension
around the point of interest is extracted from the deformation field, and
the final location is estimated by local cubic spline regression.
A resampled image or matrix of transformed points.
Jon Clayden <[email protected]>
niftyreg.linear
, niftyreg.nonlinear
This function does the opposite to decomposeAffine
, building
up an affine matrix from its components. It can be useful for testing, or
for rescaling images.
buildAffine(translation = c(0, 0, 0), scales = c(1, 1, 1), skews = c(0, 0, 0), angles = c(0, 0, 0), source = NULL, target = NULL, anchor = c("none", "origin", "centre", "center"))
buildAffine(translation = c(0, 0, 0), scales = c(1, 1, 1), skews = c(0, 0, 0), angles = c(0, 0, 0), source = NULL, target = NULL, anchor = c("none", "origin", "centre", "center"))
translation |
Translations along each axis, in |
scales |
Scale factors along each axis. |
skews |
Skews in the XY, XZ and YZ planes. |
angles |
Roll, pitch and yaw rotation angles, in radians. If
|
source |
The source image for the transformation (required). |
target |
The target image for the transformation. If |
anchor |
The fixed point for the transformation. Setting this parameter
to a value other than |
A 4x4 affine matrix representing the composite transformation. Note that NiftyReg affines logically transform backwards, from target to source space, so the matrix may be the inverse of what is expected.
Jon Clayden <[email protected]>
Compute the composition of two or more transforms, the single transform that combines their effects in order.
composeTransforms(...)
composeTransforms(...)
... |
Affine or nonlinear transforms, possibly obtained from
|
The composed transform. If all arguments are affines then the result will also be an affine; otherwise it will be a deformation field.
The source image for the composed transform is generally the source
image from the first transform, and the target is the target image from
the second transform. However, the target image attached to half
transforms (as calculated by halfTransform
) generally has a
modified xform, compared to the original target. Therefore, composing a
half transform with itself may not be exactly equivalent to the original.
Jon Clayden <[email protected]>
niftyreg.linear
, niftyreg.nonlinear
,
deformationField
An affine matrix is composed of translation, scale, skew and rotation transformations. This function extracts these components, after first inverting the matrix so that it transforms from source to target space.
decomposeAffine(affine)
decomposeAffine(affine)
affine |
A 4x4 matrix representing an affine transformation matrix. |
A list with components:
A 3x3 matrix representing only the scale operation embodied in the full affine transformation.
A 3x3 matrix representing only the skew operation embodied in the full affine transformation.
A 3x3 matrix representing only the rotation operation embodied in the full affine transformation.
A length-3 named numeric vector representing the
translations (in pixunits
units) in each of the X, Y and
Z directions.
A length-3 named numeric vector representing the scale factors in each of the X, Y and Z directions. Scale factors of 1 represent no effect.
A length-3 named numeric vector representing the skews in each of the XY, XZ and YZ planes.
A length-3 named numeric vector representing the rotation angles (in radians) about each of the X, Y and Z directions, i.e., roll, pitch and yaw.
The decomposition is not perfect, and there is one particular
degenerate case when the pitch angle is very close to pi/2
radians,
known as “Gimbal lock”. In this case the yaw angle is arbitrarily set to
zero.
Affine matrices embodying rigid-body transformations include only 6 degrees of freedom, rather than the full 12, so skews will always be zero and scales will always be unity (to within rounding error). Likewise, affine matrices derived from 2D registration will not include components relating to the Z direction.
Jon Clayden <[email protected]>
This function is used to calculate the deformation field corresponding to a specified linear or nonlinear transformation. The deformation field gives the location in source image space corresponding to the centre of each voxel in target space. It is used as a common form for linear and nonlinear transformations, and allows them to be visualised.
deformationField(transform, jacobian = TRUE)
deformationField(transform, jacobian = TRUE)
transform |
|
jacobian |
A logical value: if |
An "internalImage"
representing the deformation field. If
requested, the Jacobian map is stored in an attribute, which can be
extracted using the jacobian
accessor function.
Jon Clayden <[email protected]>
niftyreg.linear
, niftyreg.nonlinear
These functions extract forward and reverse transformations in a form
compatible with applyTransform
and other functions. They are
(S3) generic, but only methods for "niftyreg"
objects currently
exist.
forward(object, ...) ## S3 method for class 'niftyreg' forward(object, i = 1L, ...) reverse(object, ...) ## S3 method for class 'niftyreg' reverse(object, i = 1L, ...)
forward(object, ...) ## S3 method for class 'niftyreg' forward(object, i = 1L, ...) reverse(object, ...) ## S3 method for class 'niftyreg' reverse(object, i = 1L, ...)
object |
An R object. |
... |
Additional arguments. Not currently used. |
i |
The transformation number to extract. There will only be more than one in the case of multiple registration. |
A transformation object, an image or affine matrix, with suitable
attributes giving pointers to source and target images. If there is no
transformation information in the object then NULL
is returned.
Jon Clayden <[email protected]>
This function calculates the half-way transformation corresponding to its argument. Applying this transformation results in points or images in a space halfway between the original source and target images, which can be a useful common space in some applications.
halfTransform(transform)
halfTransform(transform)
transform |
The half-way transform, in a similar format to transform
.
Jon Clayden <[email protected]>
niftyreg.linear
, niftyreg.nonlinear
This function is used to invert an affine matrix. It is a wrapper around
solve
, which additionally sets appropriate attributes.
invertAffine(affine)
invertAffine(affine)
affine |
An existing 4x4 affine matrix. |
The inverted affine matrix.
Jon Clayden <[email protected]>
affine <- readAffine(system.file("extdata","affine.txt",package="RNiftyReg")) print(affine) print(invertAffine(affine))
affine <- readAffine(system.file("extdata","affine.txt",package="RNiftyReg")) print(affine) print(invertAffine(affine))
isAffine
returns a logical value indicating whether its argument is,
or resembles, a 4x4 affine matrix. asAffine
converts other objects to
the affine class, attaching or updating the source and target image
attributes. Affine transformations are a class of linear transformations
which preserve points, straight lines and planes, and may consist of a
combination of rotation, translation, scale and skew operations.
isAffine(object, strict = FALSE) asAffine(object, source = NULL, target = NULL, ...) ## S3 method for class 'niftyreg' asAffine(object, source = NULL, target = NULL, i = 1L, ...) ## S3 method for class 'niftyregRDS' asAffine(object, source = NULL, target = NULL, ...) ## S3 method for class 'affine' asAffine(object, source = NULL, target = NULL, ...) ## S3 method for class 'niftiImage' asAffine(object, source = attr(object, "source"), target = attr(object, "target"), ...) ## Default S3 method: asAffine(object, source = NULL, target = NULL, ...) ## S3 method for class 'affine' print(x, ...)
isAffine(object, strict = FALSE) asAffine(object, source = NULL, target = NULL, ...) ## S3 method for class 'niftyreg' asAffine(object, source = NULL, target = NULL, i = 1L, ...) ## S3 method for class 'niftyregRDS' asAffine(object, source = NULL, target = NULL, ...) ## S3 method for class 'affine' asAffine(object, source = NULL, target = NULL, ...) ## S3 method for class 'niftiImage' asAffine(object, source = attr(object, "source"), target = attr(object, "target"), ...) ## Default S3 method: asAffine(object, source = NULL, target = NULL, ...) ## S3 method for class 'affine' print(x, ...)
object |
An R object. |
strict |
If |
source , target
|
New source and target images for the transformation. |
... |
Additional parameters to methods. |
i |
The transformation number, for |
x |
An |
NiftyReg's convention is for affine matrices to transform world coordinates
(in the sense of voxelToWorld
) from TARGET to SOURCE space, although
transforms are logically applied the other way.
For isAffine
, a logical value, which is TRUE
if
object
appears to be an affine matrix. For asAffine
, a
classed affine object with source and target attributes set appropriately.
2D affines are a subset of 3D affines, and are stored in a 4x4 matrix for internal consistency, even though a 3x3 matrix would suffice.
Jon Clayden <[email protected]>
This function tried to determine whether an object is an image that the
package knows how to handle. If its class is "nifti"
,
"niftiImage"
, "internalImage"
or "MriImage"
, then the
result is always TRUE
. Likewise if it has an internal image pointer.
If it has no dim
attribute, or looks like an affine matrix, then the
result is FALSE
. Otherwise the value of the unsure
argument
is returned.
isImage(object, unsure = NA)
isImage(object, unsure = NA)
object |
An R object. |
unsure |
The value to return if the function can't tell whether or not
the |
Jon Clayden <[email protected]>
This function extracts the Jacobian determinant map associated with a deformation field.
jacobian(x)
jacobian(x)
x |
An R object, probably a deformation field. |
Jon Clayden <[email protected]>
The niftyreg
function performs linear or nonlinear registration for
two and three dimensional images. 4D images may also be registered
volumewise to a 3D image, or 3D images slicewise to a 2D image. This
function is a common wrapper for niftyreg.linear
and
niftyreg.nonlinear
.
niftyreg(source, target, scope = c("affine", "rigid", "nonlinear"), init = NULL, sourceMask = NULL, targetMask = NULL, symmetric = TRUE, interpolation = 3L, estimateOnly = FALSE, sequentialInit = FALSE, internal = NA, precision = c("double", "single"), threads = getOption("RNiftyReg.threads"), ...) ## S3 method for class 'niftyreg' asNifti(x, ...) ## S3 method for class 'niftyreg' as.array(x, ...)
niftyreg(source, target, scope = c("affine", "rigid", "nonlinear"), init = NULL, sourceMask = NULL, targetMask = NULL, symmetric = TRUE, interpolation = 3L, estimateOnly = FALSE, sequentialInit = FALSE, internal = NA, precision = c("double", "single"), threads = getOption("RNiftyReg.threads"), ...) ## S3 method for class 'niftyreg' asNifti(x, ...) ## S3 method for class 'niftyreg' as.array(x, ...)
source |
The source image, an object of class |
target |
The target image, an object of class |
scope |
A string describing the scope, or number of degrees of freedom
(DOF), of the registration. The currently supported values are
|
init |
Transformation(s) to be used for initialisation, which may be
|
sourceMask |
An optional mask image in source space, whose nonzero
region will be taken as the region of interest for the registration.
Ignored when |
targetMask |
An optional mask image in target space, whose nonzero region will be taken as the region of interest for the registration. |
symmetric |
Logical value. Should forward and reverse transformations be estimated simultaneously? |
interpolation |
A single integer specifying the type of interpolation to be applied to the final resampled image. May be 0 (nearest neighbour), 1 (trilinear) or 3 (cubic spline). No other values are valid. |
estimateOnly |
Logical value: if |
sequentialInit |
If |
internal |
If |
precision |
Working precision for the registration. Using single- precision may be desirable to save memory when coregistering large images. |
threads |
For OpenMP-capable builds of the package, the maximum number of threads to use. |
... |
Further arguments to |
x |
A |
A list of class "niftyreg"
with components:
An array or internal image representing the registered and
resampled source
image in the space of the target
image.
This element is NULL
if the estimateOnly
parameter is
TRUE
.
A list of (linear or nonlinear) transformations from source to target space.
A list of (linear or nonlinear) transformations from target to source space.
A list of integer vectors, giving the number of iterations completed at each “level” of the algorithm. Note that for the first level of the linear algorithm specifically, twice the specified number of iterations is allowed.
An internal representation of the source image for each registration.
An internal representation of the target image.
The as.array
method for this class returns the image
element.
If substantial parts of the target image are zero-valued, for example
because the target image has been brain-extracted, it can be useful to
pass it as a target mask as well as the target image, viz.
niftyreg(source, target, targetMask=target)
.
Jon Clayden <[email protected]>
Please see niftyreg.linear
or
niftyreg.nonlinear
for references relating to each type of
registration.
niftyreg.linear
and niftyreg.nonlinear
,
which do most of the work. Also, forward
and
reverse
to extract transformations, and
applyTransform
to apply them to new images or points.
## Not run: source <- readNifti(system.file("extdata", "epi_t2.nii.gz", package="RNiftyReg")) target <- readNifti(system.file("extdata", "flash_t1.nii.gz", package="RNiftyReg")) result <- niftyreg(source, target, scope="affine") ## End(Not run)
## Not run: source <- readNifti(system.file("extdata", "epi_t2.nii.gz", package="RNiftyReg")) target <- readNifti(system.file("extdata", "flash_t1.nii.gz", package="RNiftyReg")) result <- niftyreg(source, target, scope="affine") ## End(Not run)
The niftyreg.linear
function performs linear registration for two and
three dimensional images. 4D images may also be registered volumewise to a
3D image, or 3D images slicewise to a 2D image. Rigid-body (6 degrees of
freedom) and affine (12 degrees of freedom) registration can currently be
performed.
niftyreg.linear(source, target, scope = c("affine", "rigid"), init = NULL, sourceMask = NULL, targetMask = NULL, symmetric = TRUE, nLevels = 3L, maxIterations = 5L, useBlockPercentage = 50L, interpolation = 3L, verbose = FALSE, estimateOnly = FALSE, sequentialInit = FALSE, internal = NA, precision = c("double", "single"), threads = getOption("RNiftyReg.threads"))
niftyreg.linear(source, target, scope = c("affine", "rigid"), init = NULL, sourceMask = NULL, targetMask = NULL, symmetric = TRUE, nLevels = 3L, maxIterations = 5L, useBlockPercentage = 50L, interpolation = 3L, verbose = FALSE, estimateOnly = FALSE, sequentialInit = FALSE, internal = NA, precision = c("double", "single"), threads = getOption("RNiftyReg.threads"))
source |
The source image, an object of class |
target |
The target image, an object of class |
scope |
A string describing the scope, or number of degrees of freedom
(DOF), of the registration. The currently supported values are
|
init |
Transformation(s) to be used for initialisation, which may be
|
sourceMask |
An optional mask image in source space, whose nonzero
region will be taken as the region of interest for the registration.
Ignored when |
targetMask |
An optional mask image in target space, whose nonzero region will be taken as the region of interest for the registration. |
symmetric |
Logical value. Should forward and reverse transformations be estimated simultaneously? |
nLevels |
A single integer specifying the number of levels of the algorithm that should be applied. If zero, no optimisation will be performed, and the final affine matrix will be the same as its initialisation value. |
maxIterations |
A single integer specifying the maximum number of iterations to be used within each level. Fewer iterations may be used if a convergence test deems the process to have completed. |
useBlockPercentage |
A single integer giving the percentage of blocks to use for calculating correspondence at each step of the algorithm. The blocks with the highest intensity variance will be chosen. |
interpolation |
A single integer specifying the type of interpolation to be applied to the final resampled image. May be 0 (nearest neighbour), 1 (trilinear) or 3 (cubic spline). No other values are valid. |
verbose |
A single logical value: if |
estimateOnly |
Logical value: if |
sequentialInit |
If |
internal |
If |
precision |
Working precision for the registration. Using single- precision may be desirable to save memory when coregistering large images. |
threads |
For OpenMP-capable builds of the package, the maximum number of threads to use. |
This function performs the dual operations of finding a transformation to optimise image alignment, and resampling the source image into the space of the target image.
The algorithm is based on a block-matching approach and Least Trimmed Squares (LTS) fitting. Firstly, the block matching provides a set of corresponding points between a target and a source image. Secondly, using this set of corresponding points, the best rigid or affine transformation is evaluated. This two-step loop is repeated until convergence to the best transformation is achieved.
In the NiftyReg implementation, normalised cross-correlation between the target and source blocks is used to evaluate correspondence. The block width is constant and has been set to 4 voxels. A coarse-to-fine approach is used, where the registration is first performed on down-sampled images (using a Gaussian filter to resample images), and finally performed on full resolution images.
The source image may have 2, 3 or 4 dimensions, and the target 2 or 3. The dimensionality of the target image determines whether 2D or 3D registration is applied, and source images with one more dimension than the target (i.e. 4D to 3D, or 3D to 2D) will be registered volumewise or slicewise, as appropriate. In the latter case the last dimension of the resulting image is taken from the source image, while all other dimensions come from the target. One affine matrix is returned for each registration performed.
See niftyreg
.
Jon Clayden <[email protected]>
The algorithm used by this function is described in the following publication.
M. Modat, D.M. Cash, P. Daga, G.P. Winston, J.S. Duncan & S. Ourselin (2014). Global image registration using a symmetric block-matching approach. Journal of Medical Imaging 1(2):024003.
niftyreg
, which can be used as an interface to this
function, and niftyreg.nonlinear
for nonlinear registration.
Also, forward
and reverse
to extract
transformations, and applyTransform
to apply them to new
images or points.
The niftyreg.nonlinear
function performs nonlinear registration for
two and three dimensional images. 4D images may also be registered
volumewise to a 3D image, or 3D images slicewise to a 2D image. The warping
is based on free-form deformations, parameterised using an image of control
points.
niftyreg.nonlinear(source, target, init = NULL, sourceMask = NULL, targetMask = NULL, symmetric = TRUE, nLevels = 3L, maxIterations = 150L, nBins = 64L, bendingEnergyWeight = 0.001, linearEnergyWeight = 0.01, jacobianWeight = 0, finalSpacing = c(5, 5, 5), spacingUnit = c("voxel", "world"), interpolation = 3L, verbose = FALSE, estimateOnly = FALSE, sequentialInit = FALSE, internal = NA, precision = c("double", "single"), threads = getOption("RNiftyReg.threads"))
niftyreg.nonlinear(source, target, init = NULL, sourceMask = NULL, targetMask = NULL, symmetric = TRUE, nLevels = 3L, maxIterations = 150L, nBins = 64L, bendingEnergyWeight = 0.001, linearEnergyWeight = 0.01, jacobianWeight = 0, finalSpacing = c(5, 5, 5), spacingUnit = c("voxel", "world"), interpolation = 3L, verbose = FALSE, estimateOnly = FALSE, sequentialInit = FALSE, internal = NA, precision = c("double", "single"), threads = getOption("RNiftyReg.threads"))
source |
The source image, an object of class |
target |
The target image, an object of class |
init |
Transformation(s) to be used for initialisation, which may be
|
sourceMask |
An optional mask image in source space, whose nonzero
region will be taken as the region of interest for the registration.
Ignored when |
targetMask |
An optional mask image in target space, whose nonzero region will be taken as the region of interest for the registration. |
symmetric |
Logical value. Should forward and reverse transformations be estimated simultaneously? |
nLevels |
A single integer specifying the number of levels of the algorithm that should be applied. If zero, no optimisation will be performed, and the final control-point image will be the same as its initialisation value. |
maxIterations |
A single integer specifying the maximum number of iterations to be used within each level. Fewer iterations may be used if a convergence test deems the process to have completed. |
nBins |
A single integer giving the number of bins to use for the joint histogram created by the algorithm. |
bendingEnergyWeight |
A numeric value giving the weight of the bending energy term in the cost function. |
linearEnergyWeight |
A numeric value giving the weight of the linear energy term in the cost function. |
jacobianWeight |
A numeric value giving the weight of the Jacobian determinant term in the cost function. |
finalSpacing |
A numeric vector giving the spacing of control points in the final grid, along the X, Y and Z directions respectively. This is set from the initial control point image, if one is supplied. |
spacingUnit |
A character string giving the units in which the
|
interpolation |
A single integer specifying the type of interpolation to be applied to the final resampled image. May be 0 (nearest neighbour), 1 (trilinear) or 3 (cubic spline). No other values are valid. |
verbose |
A single logical value: if |
estimateOnly |
Logical value: if |
sequentialInit |
If |
internal |
If |
precision |
Working precision for the registration. Using single- precision may be desirable to save memory when coregistering large images. |
threads |
For OpenMP-capable builds of the package, the maximum number of threads to use. |
This function performs the dual operations of finding a transformation to
optimise image alignment, and resampling the source image into the space of
the target image (and vice-versa, if symmetric
is TRUE
).
Unlike niftyreg.linear
, this transformation is nonlinear, and
the degree of deformation may vary across the image.
The nonlinear warping is based on free-form deformations. A lattice of equally-spaced control points is defined over the target image, each of which can be moved to locally modify the mapping to the source image. In order to assess the quality of the warping between the two images, an objective function based on the normalised mutual information is used, with penalty terms based on the bending energy or the squared log of the Jacobian determinant. The objective function value is optimised using a conjugate gradient scheme.
The source image may have 2, 3 or 4 dimensions, and the target 2 or 3. The dimensionality of the target image determines whether 2D or 3D registration is applied, and source images with one more dimension than the target (i.e. 4D to 3D, or 3D to 2D) will be registered volumewise or slicewise, as appropriate. In the latter case the last dimension of the resulting image is taken from the source image, while all other dimensions come from the target. One image of control points is returned for each registration performed.
See niftyreg
.
Performing a linear registration first, and then initialising the
nonlinear transformation with the result (via the init
parameter),
is highly recommended in most circumstances.
Jon Clayden <[email protected]>
The algorithm used by this function is described in the following publication.
M. Modat, G.R. Ridgway, Z.A. Taylor, M. Lehmann, J. Barnes, D.J. Hawkes, N.C. Fox & S. Ourselin (2010). Fast free-form deformation using graphics processing units. Computer Methods and Programs in Biomedicine 98(3):278-284.
niftyreg
, which can be used as an interface to this
function, and niftyreg.linear
for linear registration. Also,
forward
and reverse
to extract
transformations, and applyTransform
to apply them to new
images or points.
This function is used to read a 4x4 numeric matrix representing an affine
transformation from a file. It is a wrapper around read.table
which
additionally ensures that required attributes are set. The type of the
matrix must be specified, as there are differing conventions across
software packages.
readAffine(fileName, source = NULL, target = NULL, type = NULL)
readAffine(fileName, source = NULL, target = NULL, type = NULL)
fileName |
A string giving the file name to read the affine matrix from. |
source |
The source image for the transformation. If |
target |
The target image for the transformation. If |
type |
The type of the affine matrix, which describes what convention
is it is stored with. Currently valid values are |
An matrix with class "affine"
, converted to the NiftyReg
convention and with source
and target
attributes set
appropriately.
Jon Clayden <[email protected]>
print(readAffine(system.file("extdata","affine.txt",package="RNiftyReg")))
print(readAffine(system.file("extdata","affine.txt",package="RNiftyReg")))
This function reads files in NIfTI-1 or NIfTI-2 format into R, using the
standard NIfTI C library. It extends the equivalent function from the
RNifti
package with source and target image parameters.
readNifti(file, source = NULL, target = NULL, internal = FALSE)
readNifti(file, source = NULL, target = NULL, internal = FALSE)
file |
A character vector of file names. |
source , target
|
If the specified |
internal |
Logical value. If |
An array or internal image, with class "niftiImage"
, and
possibly also "internalImage"
.
Jon Clayden <[email protected]>
These objects save a full transformation object, including source and target image metadata, to a self-contained RDS file, or load it back from such a file.
saveTransform(transform, fileName = NULL) loadTransform(x)
saveTransform(transform, fileName = NULL) loadTransform(x)
transform |
|
fileName |
The file name to save to. If |
x |
A file name to read from, or a serialised transform object. |
saveTransform
returns a serialised transform object, if no
filename is given; otherwise it is called for its side-effect of writing
to file. loadTransform
returns a deserialised transform object.
Jon Clayden <[email protected]>
This function calculates a similarity measure between two images, after resampling one into the space of the other. The only supported measure is currently normalised mutual information, which is also used as a cost function by the registration algorithms.
similarity(source, target, targetMask = NULL, interpolation = 3L, threads = getOption("RNiftyReg.threads"))
similarity(source, target, targetMask = NULL, interpolation = 3L, threads = getOption("RNiftyReg.threads"))
source |
The source image, in any acceptable form. |
target |
The target image. Must have the same dimensionality as the source image. |
targetMask |
An optional mask image in target space, whose nonzero region will be the area over which the measure is calculated. |
interpolation |
A single integer specifying the type of interpolation to be applied to the source image when resampling it into the space of the target image. May be 0 (nearest neighbour), 1 (trilinear) or 3 (cubic spline). No other values are valid. |
threads |
For OpenMP-capable builds of the package, the maximum number of threads to use. |
A single numeric value representing the similarity between the images.
Jon Clayden <[email protected]>
These functions allow simple transformations to be applied quickly, or in a
chosen order. They represent simplified interfaces to the
buildAffine
and applyTransform
functions, and
are compatible with the chaining operator from the popular magrittr
package (although performing one single transformation may be preferable).
translate(source, translation, ...) rescale(source, scales, anchor = c("none", "origin", "centre", "center"), ...) skew(source, skews, anchor = c("none", "origin", "centre", "center"), ...) rotate(source, angles, anchor = c("none", "origin", "centre", "center"), ...)
translate(source, translation, ...) rescale(source, scales, anchor = c("none", "origin", "centre", "center"), ...) skew(source, skews, anchor = c("none", "origin", "centre", "center"), ...) rotate(source, angles, anchor = c("none", "origin", "centre", "center"), ...)
source |
A 2D or 3D image, in the sense of |
translation |
Translations along each axis, in |
... |
Additional arguments to |
scales |
Scale factors along each axis. |
anchor |
The fixed point for the transformation. Setting this parameter
to a value other than |
skews |
Skews in the XY, XZ and YZ planes. |
angles |
Roll, pitch and yaw rotation angles, in radians. If
|
The transformed image.
Jon Clayden <[email protected]>
This function is used to write a 4x4 numeric matrix representing an affine
transformation to a file. A comment is also (optionally) written, which
specifies the matrix as using the NiftyReg convention, for the benefit of
readAffine
.
writeAffine(affine, fileName, comments = TRUE)
writeAffine(affine, fileName, comments = TRUE)
affine |
A 4x4 affine matrix. |
fileName |
A string giving the file name to write the matrix to. |
comments |
Logical value: if |
Jon Clayden <[email protected]>