NPL
Neurological Programs and Libraries
Image registration algorithms

Classes

class  npl::RigidInfoComp
 The Rigid MI Computer is used to compute the mutual information and gradient of mutual information between two images. As the name implies, it is designed for 6 parameter rigid transforms. More...
 
class  npl::RigidCorrComp
 The Rigid Corr Computer is used to compute the correlation and gradient of correlation between two images. As the name implies, it is designed for 6 parameter rigid transforms. More...
 
class  npl::ProbDistCorrInfoComp
 The distortion correction MI Computer is used to compute the mutual information and gradient of mutual information between two images using nonrigid, unidirectional, B-spline transform. In this variant the moving image should be a probability map. More...
 
class  npl::DistCorrInfoComp
 The distortion correction MI Computer is used to compute the mutual information and gradient of mutual information between two images using nonrigid, unidirectional, B-spline transform. More...
 
struct  npl::Rigid3DTrans
 Struct for holding information about a rigid transform. Note that rotation R = Rx*Ry*Rz, where Rx, Ry, and Rz are the rotations about x, y and z aaxes, and the angles are stored (in radians) in the rotation member. More...
 

Functions

ptr< MRImage > npl::motionCorrect (ptr< const MRImage > input, size_t ref)
 Performs motion correction on a set of volumes. Each 3D volume is extracted and linearly registered with the ref volume. More...
 
Rigid3DTrans npl::corReg3D (ptr< const MRImage > fixed, ptr< const MRImage > moving, const std::vector< double > &sigmas)
 Performs correlation based registration between two 3D volumes. note that the two volumes should have identical sampling and identical orientation. If that is not the case, an exception will be thrown. More...
 
Rigid3DTrans npl::informationReg3D (ptr< const MRImage > fixed, ptr< const MRImage > moving, const std::vector< double > &sigmas, size_t nbins=128, size_t binradius=4, std::string metric="MI", double stopx=0.001)
 Performs information-based registration between two 3D volumes. note that the two volumes should have identical sampling and identical orientation. If that is not the case, an exception will be thrown. More...
 
ptr< MRImage > npl::infoDistCor (ptr< const MRImage > fixed, ptr< const MRImage > moving, bool otsu, int dir, double bspace, double jac, double tps, const std::vector< double > &sigmas, size_t nbins=128, size_t binradius=4, string metric="MI", size_t hist=8, double stopx=1e-5, double beta=0.5)
 Information based registration between two 3D volumes. note that the two volumes should have identical sampling and identical orientation. More...
 
int npl::cor3DDerivTest (double step, double tol, ptr< const MRImage > in1, ptr< const MRImage > in2)
 This function checks the validity of the derivative functions used to optimize between-image corrlation. More...
 
int npl::information3DDerivTest (double step, double tol, ptr< const MRImage > in1, ptr< const MRImage > in2)
 This function checks the validity of the derivative functions used to optimize between-image corrlation. More...
 
int npl::distcorDerivTest (double step, double tol, shared_ptr< const MRImage > in1, shared_ptr< const MRImage > in2, double regj=0, double regt=0)
 This function checks the validity of the derivative functions used to optimize between-image corrlation. More...
 
std::ostream & npl::operator<< (std::ostream &stream, const Rigid3DTrans &rigid)
 Prints a rigid transform. More...
 

Detailed Description

Registration functions are implemented with classes linked to optimization functions. All registration algorithms ultimately may be performed with a simple function call, but the Computer classes (of which there is currently just RigidCorrComp) are exposed in case you want to use them for your own registration algorithms.

A computer is needed for every pair of Metric and Transform type.

Note that for rigid computers the state_x variable will store the final result. The first 3 parameters are the rotation about the X, Y, Z axes, in degrees. The second 3 will be the shift in X, Y, Z in mm.

Function Documentation

int npl::cor3DDerivTest ( double  step,
double  tol,
ptr< const MRImage in1,
ptr< const MRImage in2 
)

This function checks the validity of the derivative functions used to optimize between-image corrlation.

Parameters
stepTest step size
tolTolerance in error between analytical and Numeric gratient
in1Image 1
in2Image 2
Returns
0 if success, -1 if failure
Rigid3DTrans npl::corReg3D ( ptr< const MRImage fixed,
ptr< const MRImage moving,
const std::vector< double > &  sigmas 
)

Performs correlation based registration between two 3D volumes. note that the two volumes should have identical sampling and identical orientation. If that is not the case, an exception will be thrown.

Todo:
make it v = Ru + s, then u = INV(R)*(v - s)
Parameters
fixedImage which will be the target of registration.
movingImage which will be rotated then shifted to match fixed.
sigmasStandard deviation of smoothing at each level
Returns
Output rigid transform
int npl::distcorDerivTest ( double  step,
double  tol,
shared_ptr< const MRImage in1,
shared_ptr< const MRImage in2,
double  regj = 0,
double  regt = 0 
)

This function checks the validity of the derivative functions used to optimize between-image corrlation.

Parameters
stepTest step size
tolTolerance in error between analytical and Numeric gratient
in1Image 1
in2Image 2
regjJackobian weight
regtThin-plate-spline weight
Returns
0 if success, -1 if failure
ptr<MRImage> npl::infoDistCor ( ptr< const MRImage fixed,
ptr< const MRImage moving,
bool  otsu,
int  dir,
double  bspace,
double  jac,
double  tps,
const std::vector< double > &  sigmas,
size_t  nbins = 128,
size_t  binradius = 4,
string  metric = "MI",
size_t  hist = 8,
double  stopx = 1e-5,
double  beta = 0.5 
)

Information based registration between two 3D volumes. note that the two volumes should have identical sampling and identical orientation.

Parameters
fixedImage which will be the target of registration.
movingImage which will be rotated then shifted to match fixed.
otsuperform OTSU thresholding on images after downsampling
dirdirection/dimension of distortion
bspaceSpacing of B-Spline knots (in mm)
jacjacobian regularizer weight
tpsthin-plate-spline regularizer weight
sigmasStandard deviation of smoothing at each level
nbinsNumber of bins in marginal PDF
binradiusradius of parzen window, to smooth pdf
metricType of information based metric to use
histlength of history to keep
stopxMinimum step size before stopping
betaFraction of previous step size to consider for next step size. Essentially this decides how quickly to reduce step sizes (note 1 would mean NO reduction and would continue forever, 0 would do 1 step then stop)
Returns
parameters of bspline
int npl::information3DDerivTest ( double  step,
double  tol,
ptr< const MRImage in1,
ptr< const MRImage in2 
)

This function checks the validity of the derivative functions used to optimize between-image corrlation.

Parameters
stepTest step size
tolTolerance in error between analytical and Numeric gratient
in1Image 1
in2Image 2
Returns
0 if success, -1 if failure
Rigid3DTrans npl::informationReg3D ( ptr< const MRImage fixed,
ptr< const MRImage moving,
const std::vector< double > &  sigmas,
size_t  nbins = 128,
size_t  binradius = 4,
std::string  metric = "MI",
double  stopx = 0.001 
)

Performs information-based registration between two 3D volumes. note that the two volumes should have identical sampling and identical orientation. If that is not the case, an exception will be thrown.

Parameters
fixedImage which will be the target of registration.
movingImage which will be rotated then shifted to match fixed.
sigmasStandard deviation of smoothing kernel at each level
nbinsNumber of bins during marginal density estimation (joint with have nbins*nbins)
binradiusDuring parzen window, the radius of the smoothing kernel
metricmetric to use (default is MI)
stopxSmallest step to take, stop after steps reach this size
Returns
Rigid transform.
ptr<MRImage> npl::motionCorrect ( ptr< const MRImage input,
size_t  ref 
)

Performs motion correction on a set of volumes. Each 3D volume is extracted and linearly registered with the ref volume.

Parameters
input3+D volume (set of 3D volumes, all higher dimensions are treated equally as separate volumes).
refReference t, all images will be registered to the specified timepoint
Returns
Motion corrected volume.
std::ostream& npl::operator<< ( std::ostream &  stream,
const Rigid3DTrans rigid 
)
inline

Prints a rigid transform.

Parameters
streamOutput stream
rigidRigid transform
Returns
after this is inserted, stream

Definition at line 1088 of file registration.h.