ITK
5.0.0
Insight Segmentation and Registration Toolkit
|
#include <itkMutualInformationImageToImageMetric.h>
Computes the mutual information between two images to be registered.
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 method GetValue() computes of the mutual information while method GetValueAndDerivative() computes both the mutual information and its derivatives with respect to the transform parameters.
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 KernelFunctionBase 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 94 of file itkMutualInformationImageToImageMetric.h.
Classes | |
class | SpatialSample |
Static Public Member Functions | |
static Pointer | New () |
Static Public Member Functions inherited from itk::Object | |
static bool | GetGlobalWarningDisplay () |
static void | GlobalWarningDisplayOff () |
static void | GlobalWarningDisplayOn () |
static Pointer | New () |
static void | SetGlobalWarningDisplay (bool flag) |
Static Public Member Functions inherited from itk::LightObject | |
static void | BreakOnError () |
static Pointer | New () |
Static Public Attributes | |
static constexpr unsigned int | MovingImageDimension = MovingImageType::ImageDimension |
Static Public Attributes inherited from itk::ImageToImageMetric< TFixedImage, TMovingImage > | |
static constexpr unsigned int | FixedImageDimension = TFixedImage::ImageDimension |
static constexpr unsigned int | MovingImageDimension = TMovingImage::ImageDimension |
Private Types | |
using | CoordinateRepresentationType = typename Superclass::CoordinateRepresentationType |
using | DerivativeFunctionType = CentralDifferenceImageFunction< MovingImageType, CoordinateRepresentationType > |
using | SpatialSampleContainer = std::vector< SpatialSample > |
Private Member Functions | |
void | CalculateDerivatives (const FixedImagePointType &, DerivativeType &, TransformJacobianType &) const |
virtual void | SampleFixedImageDomain (SpatialSampleContainer &samples) const |
using itk::MutualInformationImageToImageMetric< TFixedImage, TMovingImage >::ConstPointer = SmartPointer< const Self > |
Definition at line 104 of file itkMutualInformationImageToImageMetric.h.
|
private |
Definition at line 238 of file itkMutualInformationImageToImageMetric.h.
|
private |
Definition at line 240 of file itkMutualInformationImageToImageMetric.h.
using itk::MutualInformationImageToImageMetric< TFixedImage, TMovingImage >::DerivativeType = typename Superclass::DerivativeType |
Definition at line 118 of file itkMutualInformationImageToImageMetric.h.
using itk::MutualInformationImageToImageMetric< TFixedImage, TMovingImage >::FixedImageConstPointer = typename Superclass::FixedImageConstPointer |
Definition at line 122 of file itkMutualInformationImageToImageMetric.h.
using itk::MutualInformationImageToImageMetric< TFixedImage, TMovingImage >::FixedImageIndexType = typename FixedImageType::IndexType |
Index and Point type alias support
Definition at line 126 of file itkMutualInformationImageToImageMetric.h.
using itk::MutualInformationImageToImageMetric< TFixedImage, TMovingImage >::FixedImageIndexValueType = typename FixedImageIndexType::IndexValueType |
Definition at line 127 of file itkMutualInformationImageToImageMetric.h.
using itk::MutualInformationImageToImageMetric< TFixedImage, TMovingImage >::FixedImagePointType = typename TransformType::InputPointType |
Definition at line 129 of file itkMutualInformationImageToImageMetric.h.
using itk::MutualInformationImageToImageMetric< TFixedImage, TMovingImage >::FixedImageType = typename Superclass::FixedImageType |
Definition at line 120 of file itkMutualInformationImageToImageMetric.h.
using itk::MutualInformationImageToImageMetric< TFixedImage, TMovingImage >::InterpolatorType = typename Superclass::InterpolatorType |
Definition at line 116 of file itkMutualInformationImageToImageMetric.h.
using itk::MutualInformationImageToImageMetric< TFixedImage, TMovingImage >::KernelFunctionType = KernelFunctionBase<double> |
Definition at line 132 of file itkMutualInformationImageToImageMetric.h.
using itk::MutualInformationImageToImageMetric< TFixedImage, TMovingImage >::MeasureType = typename Superclass::MeasureType |
Definition at line 117 of file itkMutualInformationImageToImageMetric.h.
using itk::MutualInformationImageToImageMetric< TFixedImage, TMovingImage >::MovingImageCosntPointer = typename Superclass::MovingImageConstPointer |
Definition at line 123 of file itkMutualInformationImageToImageMetric.h.
using itk::MutualInformationImageToImageMetric< TFixedImage, TMovingImage >::MovingImageIndexType = typename MovingImageType::IndexType |
Definition at line 128 of file itkMutualInformationImageToImageMetric.h.
using itk::MutualInformationImageToImageMetric< TFixedImage, TMovingImage >::MovingImagePointType = typename TransformType::OutputPointType |
Definition at line 130 of file itkMutualInformationImageToImageMetric.h.
using itk::MutualInformationImageToImageMetric< TFixedImage, TMovingImage >::MovingImageType = typename Superclass::MovingImageType |
Definition at line 121 of file itkMutualInformationImageToImageMetric.h.
using itk::MutualInformationImageToImageMetric< TFixedImage, TMovingImage >::ParametersType = typename Superclass::ParametersType |
Definition at line 119 of file itkMutualInformationImageToImageMetric.h.
using itk::MutualInformationImageToImageMetric< TFixedImage, TMovingImage >::Pointer = SmartPointer< Self > |
Definition at line 103 of file itkMutualInformationImageToImageMetric.h.
using itk::MutualInformationImageToImageMetric< TFixedImage, TMovingImage >::Self = MutualInformationImageToImageMetric |
Standard class type aliases.
Definition at line 101 of file itkMutualInformationImageToImageMetric.h.
|
private |
SpatialSampleContainer type alias support
Definition at line 209 of file itkMutualInformationImageToImageMetric.h.
using itk::MutualInformationImageToImageMetric< TFixedImage, TMovingImage >::Superclass = ImageToImageMetric< TFixedImage, TMovingImage > |
Definition at line 102 of file itkMutualInformationImageToImageMetric.h.
using itk::MutualInformationImageToImageMetric< TFixedImage, TMovingImage >::TransformJacobianType = typename Superclass::TransformJacobianType |
Definition at line 115 of file itkMutualInformationImageToImageMetric.h.
using itk::MutualInformationImageToImageMetric< TFixedImage, TMovingImage >::TransformPointer = typename Superclass::TransformPointer |
Definition at line 114 of file itkMutualInformationImageToImageMetric.h.
using itk::MutualInformationImageToImageMetric< TFixedImage, TMovingImage >::TransformType = typename Superclass::TransformType |
Types inherited from Superclass.
Definition at line 113 of file itkMutualInformationImageToImageMetric.h.
|
protected |
|
overrideprotecteddefault |
|
private |
Calculate the intensity derivatives at a point
|
virtual |
Create an object from an instance, potentially deferring to a factory. This method allows you to create an instance of an object that is exactly the same type as the referring object. This is useful in cases where an object has been cast back to a base class.
Reimplemented from itk::Object.
|
overridevirtual |
Get the derivatives of the match measure.
Implements itk::SingleValuedCostFunction.
|
virtual |
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.
|
virtual |
Set/Get the kernel function. This is used to calculate the joint probability distribution. Default is the GaussianKernelFunction.
|
virtual |
Set/Get the kernel function. This is used to calculate the joint probability distribution. Default is the GaussianKernelFunction.
|
virtual |
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.
|
virtual |
Run-time type information (and related methods).
Reimplemented from itk::ImageToImageMetric< TFixedImage, TMovingImage >.
|
virtual |
Get the number of spatial samples.
|
overridevirtual |
Get the value.
Implements itk::SingleValuedCostFunction.
|
overridevirtual |
Get the value and derivatives for single valued optimizers.
Reimplemented from itk::SingleValuedCostFunction.
|
static |
Method for creation through the object factory.
|
overrideprotectedvirtual |
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 >.
|
privatevirtual |
Uniformly select samples from the fixed image buffer.
|
virtual |
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.
|
virtual |
Set/Get the kernel function. This is used to calculate the joint probability distribution. Default is the GaussianKernelFunction.
|
virtual |
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.
void itk::MutualInformationImageToImageMetric< TFixedImage, TMovingImage >::SetNumberOfSpatialSamples | ( | unsigned int | num | ) |
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.
|
private |
Definition at line 242 of file itkMutualInformationImageToImageMetric.h.
|
private |
Definition at line 221 of file itkMutualInformationImageToImageMetric.h.
|
private |
Definition at line 224 of file itkMutualInformationImageToImageMetric.h.
|
private |
Definition at line 222 of file itkMutualInformationImageToImageMetric.h.
|
private |
Definition at line 220 of file itkMutualInformationImageToImageMetric.h.
|
private |
Definition at line 219 of file itkMutualInformationImageToImageMetric.h.
|
mutableprivate |
Container to store sample set A - used to approximate the probability density function (pdf).
Definition at line 213 of file itkMutualInformationImageToImageMetric.h.
|
mutableprivate |
Container to store sample set B - used to approximate the mutual information value.
Definition at line 217 of file itkMutualInformationImageToImageMetric.h.
|
static |
Enum of the moving image dimension.
Definition at line 135 of file itkMutualInformationImageToImageMetric.h.