#include <itkMutualInformationImageToImageMetric.h>
Inheritance diagram for itk::MutualInformationImageToImageMetric:
Public Types | |
typedef MutualInformationImageToImageMetric | Self |
typedef ImageToImageMetric< TFixedImage, TMovingImage > | Superclass |
typedef SmartPointer< Self > | Pointer |
typedef SmartPointer< const Self > | ConstPointer |
typedef Superclass::TransformType | TransformType |
typedef Superclass::TransformPointer | TransformPointer |
typedef Superclass::TransformJacobianType | TransformJacobianType |
typedef Superclass::InterpolatorType | InterpolatorType |
typedef Superclass::MeasureType | MeasureType |
typedef Superclass::DerivativeType | DerivativeType |
typedef Superclass::ParametersType | ParametersType |
typedef Superclass::FixedImageType | FixedImageType |
typedef Superclass::MovingImageType | MovingImageType |
typedef Superclass::FixedImageConstPointer | FixedImageConstPointer |
typedef Superclass::MovingImageConstPointer | MovingImageCosntPointer |
typedef FixedImageType::IndexType | FixedImageIndexType |
typedef FixedImageIndexType::IndexValueType | FixedImageIndexValueType |
typedef MovingImageType::IndexType | MovingImageIndexType |
typedef TransformType::InputPointType | FixedImagePointType |
typedef TransformType::OutputPointType | MovingImagePointType |
Public Methods | |
virtual const char * | GetClassName () const |
itkStaticConstMacro (MovingImageDimension, unsigned int, MovingImageType::ImageDimension) | |
void | GetDerivative (const ParametersType ¶meters, DerivativeType &Derivative) const |
MeasureType | GetValue (const ParametersType ¶meters) const |
void | GetValueAndDerivative (const ParametersType ¶meters, MeasureType &Value, DerivativeType &Derivative) const |
void | SetNumberOfSpatialSamples (unsigned int num) |
virtual unsigned int | GetNumberOfSpatialSamples () const |
virtual void | SetMovingImageStandardDeviation (double _arg) |
virtual double | GetMovingImageStandardDeviation () const |
virtual void | SetFixedImageStandardDeviation (double _arg) |
virtual double | GetFixedImageStandardDeviation () |
virtual void | SetKernelFunction (KernelFunction *_arg) |
virtual KernelFunction * | GetKernelFunction () |
Static Public Methods | |
Pointer | New () |
Protected Methods | |
MutualInformationImageToImageMetric () | |
virtual | ~MutualInformationImageToImageMetric () |
void | PrintSelf (std::ostream &os, Indent indent) const |
MutualInformationImageToImageMetric computes the mutual information between a fixed and moving image to be registered.
This class is templated over the FixedImage type and the MovingImage type.
The fixed and moving images are set via methods SetFixedImage() and SetMovingImage(). This metric makes use of user specified Transform and Interpolator. The Transform is used to map points from the fixed image to the moving image domain. The Interpolator is used to evaluate the image intensity at user specified geometric points in the moving image. The Transform and Interpolator are set via methods SetTransform() and SetInterpolator().
The calculations are based on the method of Viola and Wells where the probability density distributions are estimated using Parzen windows.
By default a Gaussian kernel is used in the density estimation. Other option include Cauchy and spline-based. A user can specify the kernel passing in a pointer a KernelFunction using the SetKernelFunction() method.
Mutual information is estimated using two sample sets: one to calculate the singular and joint pdf's and one to calculate the entropy integral. By default 50 samples points are used in each set. Other values can be set via the SetNumberOfSpatialSamples() method.
Quality of the density estimate depends on the choice of the kernel's standard deviation. Optimal choice will depend on the images. It is can be shown that around the optimal variance, the mutual information estimate is relatively insensitive to small changes of the standard deviation. In our experiments, we have found that a standard deviation of 0.4 works well for images normalized to have a mean of zero and standard deviation of 1.0. The variance can be set via methods SetFixedImageStandardDeviation() and SetMovingImageStandardDeviation().
Implementaton of this class is based on: Viola, P. and Wells III, W. (1997). "Alignment by Maximization of Mutual Information" International Journal of Computer Vision, 24(2):137-154
Definition at line 90 of file itkMutualInformationImageToImageMetric.h.
|
Reimplemented from itk::ImageToImageMetric< TFixedImage, TMovingImage >. Definition at line 99 of file itkMutualInformationImageToImageMetric.h. |
|
Type of the derivative. Reimplemented from itk::ImageToImageMetric< TFixedImage, TMovingImage >. Definition at line 113 of file itkMutualInformationImageToImageMetric.h. |
|
Reimplemented from itk::ImageToImageMetric< TFixedImage, TMovingImage >. Definition at line 117 of file itkMutualInformationImageToImageMetric.h. |
|
Index and Point typedef support. Definition at line 121 of file itkMutualInformationImageToImageMetric.h. |
|
Definition at line 122 of file itkMutualInformationImageToImageMetric.h. |
|
Definition at line 124 of file itkMutualInformationImageToImageMetric.h. |
|
Type of the fixed Image. Reimplemented from itk::ImageToImageMetric< TFixedImage, TMovingImage >. Definition at line 115 of file itkMutualInformationImageToImageMetric.h. |
|
Type of the Interpolator Base class Reimplemented from itk::ImageToImageMetric< TFixedImage, TMovingImage >. Definition at line 111 of file itkMutualInformationImageToImageMetric.h. |
|
Type of the measure. Reimplemented from itk::ImageToImageMetric< TFixedImage, TMovingImage >. Definition at line 112 of file itkMutualInformationImageToImageMetric.h. |
|
Definition at line 118 of file itkMutualInformationImageToImageMetric.h. |
|
Definition at line 123 of file itkMutualInformationImageToImageMetric.h. |
|
Definition at line 125 of file itkMutualInformationImageToImageMetric.h. |
|
Type of the moving Image. Reimplemented from itk::ImageToImageMetric< TFixedImage, TMovingImage >. Definition at line 116 of file itkMutualInformationImageToImageMetric.h. |
|
Type of the parameters. Reimplemented from itk::ImageToImageMetric< TFixedImage, TMovingImage >. Definition at line 114 of file itkMutualInformationImageToImageMetric.h. |
|
Reimplemented from itk::ImageToImageMetric< TFixedImage, TMovingImage >. Definition at line 98 of file itkMutualInformationImageToImageMetric.h. |
|
Standard class typedefs. Reimplemented from itk::ImageToImageMetric< TFixedImage, TMovingImage >. Definition at line 96 of file itkMutualInformationImageToImageMetric.h. |
|
Reimplemented from itk::ImageToImageMetric< TFixedImage, TMovingImage >. Definition at line 97 of file itkMutualInformationImageToImageMetric.h. |
|
Reimplemented from itk::ImageToImageMetric< TFixedImage, TMovingImage >. Definition at line 110 of file itkMutualInformationImageToImageMetric.h. |
|
Reimplemented from itk::ImageToImageMetric< TFixedImage, TMovingImage >. Definition at line 109 of file itkMutualInformationImageToImageMetric.h. |
|
Types inherited from Superclass. Reimplemented from itk::ImageToImageMetric< TFixedImage, TMovingImage >. Definition at line 108 of file itkMutualInformationImageToImageMetric.h. |
|
|
|
Definition at line 180 of file itkMutualInformationImageToImageMetric.h. |
|
Run-time type information (and related methods). Reimplemented from itk::ImageToImageMetric< TFixedImage, TMovingImage >. |
|
Get the derivatives of the match measure. Implements itk::SingleValuedCostFunction. |
|
Set/Get the fixed image intensitiy standard deviation. This defines the kernel bandwidth used in the joint probability distribution calculation. Default value is 0.4 which works well for image intensities normalized to a mean of 0 and standard deviation of 1.0. Value is clamped to be always greater than zero. |
|
Set/Get the kernel function. This is used to calculate the joint probability distribution. Default is the GaussianKernelFunction. |
|
Set/Get the moving image intensitiy standard deviation. This defines the kernel bandwidth used in the joint probability distribution calculation. Default value is 0.4 which works well for image intensities normalized to a mean of 0 and standard deviation of 1.0. Value is clamped to be always greater than zero. |
|
Get the number of spatial samples. |
|
Get the value. Implements itk::SingleValuedCostFunction. |
|
Get the value and derivatives for single valued optimizers. Reimplemented from itk::SingleValuedCostFunction. |
|
Enum of the moving image dimension. |
|
Method for creation through the object factory. Reimplemented from itk::Object. |
|
Methods invoked by Print() to print information about the object including superclasses. Typically not called by the user (use Print() instead) but used in the hierarchical print process to combine the output of several classes. Reimplemented from itk::ImageToImageMetric< TFixedImage, TMovingImage >. |
|
Set/Get the fixed image intensitiy standard deviation. This defines the kernel bandwidth used in the joint probability distribution calculation. Default value is 0.4 which works well for image intensities normalized to a mean of 0 and standard deviation of 1.0. Value is clamped to be always greater than zero. |
|
Set/Get the kernel function. This is used to calculate the joint probability distribution. Default is the GaussianKernelFunction. |
|
Set/Get the moving image intensitiy standard deviation. This defines the kernel bandwidth used in the joint probability distribution calculation. Default value is 0.4 which works well for image intensities normalized to a mean of 0 and standard deviation of 1.0. Value is clamped to be always greater than zero. |
|
Set the number of spatial samples. This is the number of image samples used to calculate the joint probability distribution. The number of spatial samples is clamped to be a minimum of 1. Default value is 50. |