ITK  4.0.0
Insight Segmentation and Registration Toolkit
Public Types | Public Member Functions | Static Public Member Functions | Static Public Attributes | Protected Member Functions | Private Types | Private Member Functions | Private Attributes
itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage > Class Template Reference

Computes the mutual information between two images to be registered using the method of Mattes et al. More...

#include <itkMattesMutualInformationImageToImageMetric.h>

Inheritance diagram for itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >:
Collaboration diagram for itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >:

List of all members.

Public Types

typedef
Superclass::BSplineTransformIndexArrayType 
BSplineTransformIndexArrayType
typedef
Superclass::BSplineTransformWeightsType 
BSplineTransformWeightsType
typedef SmartPointer< const SelfConstPointer
typedef
Superclass::CoordinateRepresentationType 
CoordinateRepresentationType
typedef Superclass::DerivativeType DerivativeType
typedef
Superclass::FixedImageConstPointer 
FixedImageConstPointer
typedef
Superclass::FixedImageSampleContainer 
FixedImageSampleContainer
typedef Superclass::FixedImageType FixedImageType
typedef
Superclass::ImageDerivativesType 
ImageDerivativesType
typedef Superclass::IndexValueType IndexValueType
typedef
Superclass::InterpolatorType 
InterpolatorType
typedef Image< PDFValueType, 3 > JointPDFDerivativesType
typedef Image< PDFValueType, 2 > JointPDFType
typedef Superclass::MeasureType MeasureType
typedef
Superclass::MovingImageConstPointer 
MovingImageConstPointer
typedef
Superclass::MovingImagePointType 
MovingImagePointType
typedef Superclass::MovingImageType MovingImageType
typedef
FixedImageType::OffsetValueType 
OffsetValueType
typedef Superclass::ParametersType ParametersType
typedef float PDFValueType
typedef SmartPointer< SelfPointer
typedef
MattesMutualInformationImageToImageMetric 
Self
typedef ImageToImageMetric
< TFixedImage, TMovingImage > 
Superclass
typedef
Superclass::TransformJacobianType 
TransformJacobianType
typedef
Superclass::TransformPointer 
TransformPointer
typedef Superclass::TransformType TransformType
typedef
Superclass::WeightsValueType 
WeightsValueType

Public Member Functions

virtual ::itk::LightObject::Pointer CreateAnother (void) const
void GetDerivative (const ParametersType &parameters, DerivativeType &Derivative) const
virtual const char * GetNameOfClass () const
MeasureType GetValue (const ParametersType &parameters) const
void GetValueAndDerivative (const ParametersType &parameters, MeasureType &Value, DerivativeType &Derivative) const
virtual void Initialize (void) throw ( ExceptionObject )
virtual void SetNumberOfHistogramBins (SizeValueType _arg)
virtual const SizeValueTypeGetNumberOfHistogramBins ()
virtual void SetUseExplicitPDFDerivatives (bool _arg)
virtual const bool & GetUseExplicitPDFDerivatives ()
virtual void UseExplicitPDFDerivativesOn ()
virtual void UseExplicitPDFDerivativesOff ()
const JointPDFType::Pointer GetJointPDF () const
const
JointPDFDerivativesType::Pointer 
GetJointPDFDerivatives () const

Static Public Member Functions

static Pointer New ()

Static Public Attributes

static const unsigned int MovingImageDimension = MovingImageType::ImageDimension

Protected Member Functions

 MattesMutualInformationImageToImageMetric ()
void PrintSelf (std::ostream &os, Indent indent) const
virtual ~MattesMutualInformationImageToImageMetric ()

Private Types

typedef
BSplineDerivativeKernelFunction< 3 > 
CubicBSplineDerivativeFunctionType
typedef BSplineKernelFunction< 3 > CubicBSplineFunctionType
typedef
JointPDFDerivativesType::IndexType 
JointPDFDerivativesIndexType
typedef
JointPDFDerivativesType::RegionType 
JointPDFDerivativesRegionType
typedef
JointPDFDerivativesType::SizeType 
JointPDFDerivativesSizeType
typedef
JointPDFDerivativesType::PixelType 
JointPDFDerivativesValueType
typedef JointPDFType::IndexType JointPDFIndexType
typedef JointPDFType::RegionType JointPDFRegionType
typedef JointPDFType::SizeType JointPDFSizeType
typedef JointPDFType::PixelType JointPDFValueType
typedef Array2D< PRatioTypePRatioArrayType
typedef double PRatioType

Private Member Functions

void ComputeFixedImageParzenWindowIndices (FixedImageSampleContainer &samples)
void ComputePDFDerivatives (ThreadIdType threadID, unsigned int sampleNumber, int movingImageParzenWindowIndex, const ImageDerivativesType &movingImageGradientValue, double cubicBSplineDerivativeValue) const
void GetValueAndDerivativeThreadPostProcess (ThreadIdType threadID, bool withinSampleThread) const
void GetValueAndDerivativeThreadPreProcess (ThreadIdType threadID, bool withinSampleThread) const
bool GetValueAndDerivativeThreadProcessSample (ThreadIdType threadID, SizeValueType fixedImageSample, const MovingImagePointType &mappedPoint, double movingImageValue, const ImageDerivativesType &movingImageGradientValue) const
void GetValueThreadPostProcess (ThreadIdType threadID, bool withinSampleThread) const
void GetValueThreadPreProcess (ThreadIdType threadID, bool withinSampleThread) const
bool GetValueThreadProcessSample (ThreadIdType threadID, SizeValueType fixedImageSample, const MovingImagePointType &mappedPoint, double movingImageValue) const
 MattesMutualInformationImageToImageMetric (const Self &)
void operator= (const Self &)

Private Attributes

CubicBSplineDerivativeFunctionType::Pointer m_CubicBSplineDerivativeKernel
CubicBSplineFunctionType::Pointer m_CubicBSplineKernel
double m_FixedImageBinSize
double m_FixedImageNormalizedMin
double m_FixedImageTrueMax
double m_FixedImageTrueMin
bool m_ImplicitDerivativesSecondPass
double m_MovingImageBinSize
std::vector< PDFValueTypem_MovingImageMarginalPDF
double m_MovingImageNormalizedMin
double m_MovingImageTrueMax
double m_MovingImageTrueMin
SizeValueType m_NumberOfHistogramBins
PRatioArrayType m_PRatioArray
std::vector< std::vector
< PDFValueType > > 
m_ThreaderFixedImageMarginalPDF
std::vector
< JointPDFType::Pointer
m_ThreaderJointPDF
std::vector
< JointPDFDerivativesType::Pointer
m_ThreaderJointPDFDerivatives
std::vector< int > m_ThreaderJointPDFEndBin
std::vector< int > m_ThreaderJointPDFStartBin
std::vector< double > m_ThreaderJointPDFSum
std::vector< DerivativeTypem_ThreaderMetricDerivative
bool m_UseExplicitPDFDerivatives

Detailed Description

template<class TFixedImage, class TMovingImage>
class itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >

Computes the mutual information between two images to be registered using the method of Mattes et al.

MattesMutualInformationImageToImageMetric 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().

If a BSplineInterpolationFunction is used, this class obtain image derivatives from the BSpline interpolator. Otherwise, image derivatives are computed using central differencing.

Warning:
This metric assumes that the moving image has already been connected to the interpolator outside of this class.

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 Mattes et al [1,2] where the probability density distribution are estimated using Parzen histograms. Since the fixed image PDF does not contribute to the derivatives, it does not need to be smooth. Hence, a zero order (box car) BSpline kernel is used for the fixed image intensity PDF. On the other hand, to ensure smoothness a third order BSpline kernel is used for the moving image intensity PDF.

On Initialize(), the FixedImage is uniformly sampled within the FixedImageRegion. The number of samples used can be set via SetNumberOfSpatialSamples(). Typically, the number of spatial samples used should increase with the image size.

The option UseAllPixelOn() disables the random sampling and uses all the pixels of the FixedImageRegion in order to estimate the joint intensity PDF.

During each call of GetValue(), GetDerivatives(), GetValueAndDerivatives(), marginal and joint intensity PDF's values are estimated at discrete position or bins. The number of bins used can be set via SetNumberOfHistogramBins(). To handle data with arbitray magnitude and dynamic range, the image intensity is scale such that any contribution to the histogram will fall into a valid bin.

One the PDF's have been contructed, the mutual information is obtained by doubling summing over the discrete PDF values.

Notes: 1. This class returns the negative mutual information value. 2. This class in not thread safe due the private data structures used to the store the sampled points and the marginal and joint pdfs.

References: [1] "Nonrigid multimodality image registration" D. Mattes, D. R. Haynor, H. Vesselle, T. Lewellen and W. Eubank Medical Imaging 2001: Image Processing, 2001, pp. 1609-1620. [2] "PET-CT Image Registration in the Chest Using Free-form Deformations" D. Mattes, D. R. Haynor, H. Vesselle, T. Lewellen and W. Eubank IEEE Transactions in Medical Imaging. Vol.22, No.1, January 2003. pp.120-128. [3] "Optimization of Mutual Information for MultiResolution Image Registration" P. Thevenaz and M. Unser IEEE Transactions in Image Processing, 9(12) December 2000.

Definition at line 111 of file itkMattesMutualInformationImageToImageMetric.h.


Member Typedef Documentation

template<class TFixedImage , class TMovingImage >
typedef Superclass::BSplineTransformIndexArrayType itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::BSplineTransformIndexArrayType
template<class TFixedImage , class TMovingImage >
typedef Superclass::BSplineTransformWeightsType itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::BSplineTransformWeightsType
template<class TFixedImage , class TMovingImage >
typedef SmartPointer<const Self> itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::ConstPointer
template<class TFixedImage , class TMovingImage >
typedef Superclass::CoordinateRepresentationType itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::CoordinateRepresentationType

Type used for representing point components

Reimplemented from itk::ImageToImageMetric< TFixedImage, TMovingImage >.

Definition at line 145 of file itkMattesMutualInformationImageToImageMetric.h.

template<class TFixedImage , class TMovingImage >
typedef BSplineDerivativeKernelFunction<3> itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::CubicBSplineDerivativeFunctionType [private]
template<class TFixedImage , class TMovingImage >
typedef BSplineKernelFunction<3> itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::CubicBSplineFunctionType [private]

Typedefs for BSpline kernel and derivative functions.

Definition at line 274 of file itkMattesMutualInformationImageToImageMetric.h.

template<class TFixedImage , class TMovingImage >
typedef Superclass::DerivativeType itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::DerivativeType

Type of the derivative.

Reimplemented from itk::ImageToImageMetric< TFixedImage, TMovingImage >.

Definition at line 135 of file itkMattesMutualInformationImageToImageMetric.h.

template<class TFixedImage , class TMovingImage >
typedef Superclass::FixedImageConstPointer itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::FixedImageConstPointer
template<class TFixedImage , class TMovingImage >
typedef Superclass::FixedImageSampleContainer itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::FixedImageSampleContainer

FixedImageSamplePoint typedef support.

Reimplemented from itk::ImageToImageMetric< TFixedImage, TMovingImage >.

Definition at line 146 of file itkMattesMutualInformationImageToImageMetric.h.

template<class TFixedImage , class TMovingImage >
typedef Superclass::FixedImageType itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::FixedImageType

Type of the fixed Image.

Reimplemented from itk::ImageToImageMetric< TFixedImage, TMovingImage >.

Definition at line 137 of file itkMattesMutualInformationImageToImageMetric.h.

template<class TFixedImage , class TMovingImage >
typedef Superclass::ImageDerivativesType itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::ImageDerivativesType
template<class TFixedImage , class TMovingImage >
typedef Superclass::IndexValueType itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::IndexValueType
template<class TFixedImage , class TMovingImage >
typedef Superclass::InterpolatorType itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::InterpolatorType

Type of the Interpolator Base class

Reimplemented from itk::ImageToImageMetric< TFixedImage, TMovingImage >.

Definition at line 133 of file itkMattesMutualInformationImageToImageMetric.h.

template<class TFixedImage , class TMovingImage >
typedef JointPDFDerivativesType::IndexType itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::JointPDFDerivativesIndexType [private]
template<class TFixedImage , class TMovingImage >
typedef JointPDFDerivativesType::RegionType itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::JointPDFDerivativesRegionType [private]
template<class TFixedImage , class TMovingImage >
typedef JointPDFDerivativesType::SizeType itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::JointPDFDerivativesSizeType [private]
template<class TFixedImage , class TMovingImage >
typedef Image<PDFValueType, 3> itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::JointPDFDerivativesType
template<class TFixedImage , class TMovingImage >
typedef JointPDFDerivativesType::PixelType itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::JointPDFDerivativesValueType [private]
template<class TFixedImage , class TMovingImage >
typedef JointPDFType::IndexType itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::JointPDFIndexType [private]
template<class TFixedImage , class TMovingImage >
typedef JointPDFType::RegionType itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::JointPDFRegionType [private]
template<class TFixedImage , class TMovingImage >
typedef JointPDFType::SizeType itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::JointPDFSizeType [private]
template<class TFixedImage , class TMovingImage >
typedef Image<PDFValueType, 2> itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::JointPDFType

Typedef for the joint PDF and PDF derivatives are stored as ITK Images.

Definition at line 219 of file itkMattesMutualInformationImageToImageMetric.h.

template<class TFixedImage , class TMovingImage >
typedef JointPDFType::PixelType itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::JointPDFValueType [private]
template<class TFixedImage , class TMovingImage >
typedef Superclass::MeasureType itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::MeasureType

Type of the measure.

Reimplemented from itk::ImageToImageMetric< TFixedImage, TMovingImage >.

Definition at line 134 of file itkMattesMutualInformationImageToImageMetric.h.

template<class TFixedImage , class TMovingImage >
typedef Superclass::MovingImageConstPointer itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::MovingImageConstPointer
template<class TFixedImage , class TMovingImage >
typedef Superclass::MovingImagePointType itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::MovingImagePointType
template<class TFixedImage , class TMovingImage >
typedef Superclass::MovingImageType itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::MovingImageType

Type of the moving Image.

Reimplemented from itk::ImageToImageMetric< TFixedImage, TMovingImage >.

Definition at line 138 of file itkMattesMutualInformationImageToImageMetric.h.

template<class TFixedImage , class TMovingImage >
typedef FixedImageType::OffsetValueType itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::OffsetValueType
template<class TFixedImage , class TMovingImage >
typedef Superclass::ParametersType itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::ParametersType

Type of the parameters.

Reimplemented from itk::ImageToImageMetric< TFixedImage, TMovingImage >.

Definition at line 136 of file itkMattesMutualInformationImageToImageMetric.h.

template<class TFixedImage , class TMovingImage >
typedef float itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::PDFValueType

The marginal PDFs are stored as std::vector.

Definition at line 212 of file itkMattesMutualInformationImageToImageMetric.h.

template<class TFixedImage , class TMovingImage >
typedef SmartPointer<Self> itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::Pointer
template<class TFixedImage , class TMovingImage >
typedef Array2D<PRatioType> itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::PRatioArrayType [private]
template<class TFixedImage , class TMovingImage >
typedef double itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::PRatioType [private]

Helper array for storing the values of the JointPDF ratios.

Definition at line 318 of file itkMattesMutualInformationImageToImageMetric.h.

template<class TFixedImage , class TMovingImage >
typedef MattesMutualInformationImageToImageMetric itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::Self

Standard class typedefs.

Reimplemented from itk::ImageToImageMetric< TFixedImage, TMovingImage >.

Definition at line 117 of file itkMattesMutualInformationImageToImageMetric.h.

template<class TFixedImage , class TMovingImage >
typedef ImageToImageMetric<TFixedImage, TMovingImage> itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::Superclass
template<class TFixedImage , class TMovingImage >
typedef Superclass::TransformJacobianType itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::TransformJacobianType
template<class TFixedImage , class TMovingImage >
typedef Superclass::TransformPointer itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::TransformPointer
template<class TFixedImage , class TMovingImage >
typedef Superclass::TransformType itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::TransformType

Types inherited from Superclass.

Reimplemented from itk::ImageToImageMetric< TFixedImage, TMovingImage >.

Definition at line 127 of file itkMattesMutualInformationImageToImageMetric.h.

template<class TFixedImage , class TMovingImage >
typedef Superclass::WeightsValueType itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::WeightsValueType

Constructor & Destructor Documentation

template<class TFixedImage , class TMovingImage >
itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::MattesMutualInformationImageToImageMetric ( ) [protected]
template<class TFixedImage , class TMovingImage >
virtual itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::~MattesMutualInformationImageToImageMetric ( ) [protected, virtual]
template<class TFixedImage , class TMovingImage >
itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::MattesMutualInformationImageToImageMetric ( const Self ) [private]

Member Function Documentation

template<class TFixedImage , class TMovingImage >
void itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::ComputeFixedImageParzenWindowIndices ( FixedImageSampleContainer samples) [private]

Precompute fixed image parzen window indices.

template<class TFixedImage , class TMovingImage >
void itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::ComputePDFDerivatives ( ThreadIdType  threadID,
unsigned int  sampleNumber,
int  movingImageParzenWindowIndex,
const ImageDerivativesType movingImageGradientValue,
double  cubicBSplineDerivativeValue 
) const [private]

Compute PDF derivative contribution for each parameter.

template<class TFixedImage , class TMovingImage >
virtual::itk::LightObject::Pointer itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::CreateAnother ( void  ) const [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.

template<class TFixedImage , class TMovingImage >
void itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::GetDerivative ( const ParametersType parameters,
DerivativeType Derivative 
) const [virtual]

Get the derivatives of the match measure.

Implements itk::SingleValuedCostFunction.

template<class TFixedImage , class TMovingImage >
const JointPDFType::Pointer itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::GetJointPDF ( ) const [inline]

Get the internal JointPDF image that was used in creating the metric value.

Definition at line 226 of file itkMattesMutualInformationImageToImageMetric.h.

References NULL.

template<class TFixedImage , class TMovingImage >
const JointPDFDerivativesType::Pointer itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::GetJointPDFDerivatives ( ) const [inline]

Get the internal JointPDFDeriviative image that was used in creating the metric derivative value. This is only created when UseExplicitPDFDerivatives is ON, and derivatives are requested.

Definition at line 242 of file itkMattesMutualInformationImageToImageMetric.h.

References NULL.

template<class TFixedImage , class TMovingImage >
virtual const char* itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::GetNameOfClass ( ) const [virtual]

Run-time type information (and related methods).

Reimplemented from itk::ImageToImageMetric< TFixedImage, TMovingImage >.

template<class TFixedImage , class TMovingImage >
virtual const SizeValueType& itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::GetNumberOfHistogramBins ( ) [virtual]

Number of bins to used in the histogram. Typical value is 50. The minimum value is 5 due to the padding required by the Parzen windowing with a cubic-BSpline kernel. Note that even if the metric is used on binary images, the number of bins should at least be equal to five.

template<class TFixedImage , class TMovingImage >
virtual const bool& itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::GetUseExplicitPDFDerivatives ( ) [virtual]

This variable selects the method to be used for computing the Metric derivatives with respect to the Transform parameters. Two modes of computation are available. The choice between one and the other is a trade-off between computation speed and memory allocations. The two modes are described in detail below:

UseExplicitPDFDerivatives = True will compute the Metric derivative by first calculating the derivatives of each one of the Joint PDF bins with respect to each one of the Transform parameters and then accumulating these contributions in the final metric derivative array by using a bin-specific weight. The memory required for storing the intermediate derivatives is a 3D array of doubles with size equals to the product of (number of histogram bins)^2 times number of transform parameters. This method is well suited for Transform with a small number of parameters.

UseExplicitPDFDerivatives = False will compute the Metric derivative by first computing the weights for each one of the Joint PDF bins and caching them into an array. Then it will revisit each one of the PDF bins for computing its weighted contribution to the full derivative array. In this method an extra 2D array is used for storing the weights of each one of the PDF bins. This is an array of doubles with size equals to (number of histogram bins)^2. This method is well suited for Transforms with a large number of parameters, such as, BSplineTransforms.

template<class TFixedImage , class TMovingImage >
MeasureType itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::GetValue ( const ParametersType parameters) const [virtual]

Get the value.

Implements itk::SingleValuedCostFunction.

template<class TFixedImage , class TMovingImage >
void itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::GetValueAndDerivative ( const ParametersType parameters,
MeasureType Value,
DerivativeType Derivative 
) const [virtual]

Get the value and derivatives for single valued optimizers.

Reimplemented from itk::SingleValuedCostFunction.

template<class TFixedImage , class TMovingImage >
void itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::GetValueAndDerivativeThreadPostProcess ( ThreadIdType  threadID,
bool  withinSampleThread 
) const [private, virtual]
template<class TFixedImage , class TMovingImage >
void itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::GetValueAndDerivativeThreadPreProcess ( ThreadIdType  threadID,
bool  withinSampleThread 
) const [private, virtual]
template<class TFixedImage , class TMovingImage >
bool itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::GetValueAndDerivativeThreadProcessSample ( ThreadIdType  threadID,
SizeValueType  fixedImageSample,
const MovingImagePointType mappedPoint,
double  movingImageValue,
const ImageDerivativesType movingImageGradientValue 
) const [private, virtual]
template<class TFixedImage , class TMovingImage >
void itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::GetValueThreadPostProcess ( ThreadIdType  threadID,
bool  withinSampleThread 
) const [private, virtual]
template<class TFixedImage , class TMovingImage >
void itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::GetValueThreadPreProcess ( ThreadIdType  threadID,
bool  withinSampleThread 
) const [private, virtual]
template<class TFixedImage , class TMovingImage >
bool itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::GetValueThreadProcessSample ( ThreadIdType  threadID,
SizeValueType  fixedImageSample,
const MovingImagePointType mappedPoint,
double  movingImageValue 
) const [private, virtual]
template<class TFixedImage , class TMovingImage >
virtual void itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::Initialize ( void  ) throw ( ExceptionObject ) [virtual]

Initialize the Metric by (1) making sure that all the components are present and plugged together correctly, (2) uniformly select NumberOfSpatialSamples within the FixedImageRegion, and (3) allocate memory for pdf data structures.

Reimplemented from itk::ImageToImageMetric< TFixedImage, TMovingImage >.

template<class TFixedImage , class TMovingImage >
static Pointer itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::New ( ) [static]

Method for creation through the object factory.

Reimplemented from itk::Object.

template<class TFixedImage , class TMovingImage >
void itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::operator= ( const Self ) [private]

Mutex lock to protect modification to the reference count

Reimplemented from itk::ImageToImageMetric< TFixedImage, TMovingImage >.

template<class TFixedImage , class TMovingImage >
void itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::PrintSelf ( std::ostream &  os,
Indent  indent 
) const [protected, virtual]

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 >.

template<class TFixedImage , class TMovingImage >
virtual void itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::SetNumberOfHistogramBins ( SizeValueType  _arg) [virtual]

Number of bins to used in the histogram. Typical value is 50. The minimum value is 5 due to the padding required by the Parzen windowing with a cubic-BSpline kernel. Note that even if the metric is used on binary images, the number of bins should at least be equal to five.

template<class TFixedImage , class TMovingImage >
virtual void itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::SetUseExplicitPDFDerivatives ( bool  _arg) [virtual]

This variable selects the method to be used for computing the Metric derivatives with respect to the Transform parameters. Two modes of computation are available. The choice between one and the other is a trade-off between computation speed and memory allocations. The two modes are described in detail below:

UseExplicitPDFDerivatives = True will compute the Metric derivative by first calculating the derivatives of each one of the Joint PDF bins with respect to each one of the Transform parameters and then accumulating these contributions in the final metric derivative array by using a bin-specific weight. The memory required for storing the intermediate derivatives is a 3D array of doubles with size equals to the product of (number of histogram bins)^2 times number of transform parameters. This method is well suited for Transform with a small number of parameters.

UseExplicitPDFDerivatives = False will compute the Metric derivative by first computing the weights for each one of the Joint PDF bins and caching them into an array. Then it will revisit each one of the PDF bins for computing its weighted contribution to the full derivative array. In this method an extra 2D array is used for storing the weights of each one of the PDF bins. This is an array of doubles with size equals to (number of histogram bins)^2. This method is well suited for Transforms with a large number of parameters, such as, BSplineTransforms.

template<class TFixedImage , class TMovingImage >
virtual void itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::UseExplicitPDFDerivativesOff ( ) [virtual]

This variable selects the method to be used for computing the Metric derivatives with respect to the Transform parameters. Two modes of computation are available. The choice between one and the other is a trade-off between computation speed and memory allocations. The two modes are described in detail below:

UseExplicitPDFDerivatives = True will compute the Metric derivative by first calculating the derivatives of each one of the Joint PDF bins with respect to each one of the Transform parameters and then accumulating these contributions in the final metric derivative array by using a bin-specific weight. The memory required for storing the intermediate derivatives is a 3D array of doubles with size equals to the product of (number of histogram bins)^2 times number of transform parameters. This method is well suited for Transform with a small number of parameters.

UseExplicitPDFDerivatives = False will compute the Metric derivative by first computing the weights for each one of the Joint PDF bins and caching them into an array. Then it will revisit each one of the PDF bins for computing its weighted contribution to the full derivative array. In this method an extra 2D array is used for storing the weights of each one of the PDF bins. This is an array of doubles with size equals to (number of histogram bins)^2. This method is well suited for Transforms with a large number of parameters, such as, BSplineTransforms.

template<class TFixedImage , class TMovingImage >
virtual void itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::UseExplicitPDFDerivativesOn ( ) [virtual]

This variable selects the method to be used for computing the Metric derivatives with respect to the Transform parameters. Two modes of computation are available. The choice between one and the other is a trade-off between computation speed and memory allocations. The two modes are described in detail below:

UseExplicitPDFDerivatives = True will compute the Metric derivative by first calculating the derivatives of each one of the Joint PDF bins with respect to each one of the Transform parameters and then accumulating these contributions in the final metric derivative array by using a bin-specific weight. The memory required for storing the intermediate derivatives is a 3D array of doubles with size equals to the product of (number of histogram bins)^2 times number of transform parameters. This method is well suited for Transform with a small number of parameters.

UseExplicitPDFDerivatives = False will compute the Metric derivative by first computing the weights for each one of the Joint PDF bins and caching them into an array. Then it will revisit each one of the PDF bins for computing its weighted contribution to the full derivative array. In this method an extra 2D array is used for storing the weights of each one of the PDF bins. This is an array of doubles with size equals to (number of histogram bins)^2. This method is well suited for Transforms with a large number of parameters, such as, BSplineTransforms.


Member Data Documentation

template<class TFixedImage , class TMovingImage >
CubicBSplineDerivativeFunctionType::Pointer itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::m_CubicBSplineDerivativeKernel [private]
template<class TFixedImage , class TMovingImage >
CubicBSplineFunctionType::Pointer itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::m_CubicBSplineKernel [private]

Cubic BSpline kernel for computing Parzen histograms.

Definition at line 314 of file itkMattesMutualInformationImageToImageMetric.h.

template<class TFixedImage , class TMovingImage >
double itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::m_FixedImageBinSize [private]
template<class TFixedImage , class TMovingImage >
double itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::m_FixedImageNormalizedMin [private]
template<class TFixedImage , class TMovingImage >
double itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::m_FixedImageTrueMax [private]
template<class TFixedImage , class TMovingImage >
double itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::m_FixedImageTrueMin [private]
template<class TFixedImage , class TMovingImage >
bool itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::m_ImplicitDerivativesSecondPass [mutable, private]
template<class TFixedImage , class TMovingImage >
double itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::m_MovingImageBinSize [private]
template<class TFixedImage , class TMovingImage >
std::vector<PDFValueType> itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::m_MovingImageMarginalPDF [mutable, private]

The moving image marginal PDF.

Definition at line 327 of file itkMattesMutualInformationImageToImageMetric.h.

template<class TFixedImage , class TMovingImage >
double itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::m_MovingImageNormalizedMin [private]
template<class TFixedImage , class TMovingImage >
double itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::m_MovingImageTrueMax [private]
template<class TFixedImage , class TMovingImage >
double itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::m_MovingImageTrueMin [private]
template<class TFixedImage , class TMovingImage >
SizeValueType itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::m_NumberOfHistogramBins [private]

Variables to define the marginal and joint histograms.

Definition at line 303 of file itkMattesMutualInformationImageToImageMetric.h.

template<class TFixedImage , class TMovingImage >
PRatioArrayType itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::m_PRatioArray [mutable, private]
template<class TFixedImage , class TMovingImage >
std::vector<std::vector<PDFValueType> > itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::m_ThreaderFixedImageMarginalPDF [mutable, private]
template<class TFixedImage , class TMovingImage >
std::vector<JointPDFType::Pointer> itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::m_ThreaderJointPDF [private]

The joint PDF and PDF derivatives.

Definition at line 331 of file itkMattesMutualInformationImageToImageMetric.h.

template<class TFixedImage , class TMovingImage >
std::vector<JointPDFDerivativesType::Pointer> itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::m_ThreaderJointPDFDerivatives [private]
template<class TFixedImage , class TMovingImage >
std::vector<int> itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::m_ThreaderJointPDFEndBin [private]
template<class TFixedImage , class TMovingImage >
std::vector<int> itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::m_ThreaderJointPDFStartBin [private]
template<class TFixedImage , class TMovingImage >
std::vector<double> itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::m_ThreaderJointPDFSum [mutable, private]
template<class TFixedImage , class TMovingImage >
std::vector<DerivativeType> itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::m_ThreaderMetricDerivative [mutable, private]

Helper variable for accumulating the derivative of the metric.

Definition at line 324 of file itkMattesMutualInformationImageToImageMetric.h.

template<class TFixedImage , class TMovingImage >
bool itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::m_UseExplicitPDFDerivatives [private]
template<class TFixedImage , class TMovingImage >
const unsigned int itk::MattesMutualInformationImageToImageMetric< TFixedImage, TMovingImage >::MovingImageDimension = MovingImageType::ImageDimension [static]

The moving image dimension.

Reimplemented from itk::ImageToImageMetric< TFixedImage, TMovingImage >.

Definition at line 155 of file itkMattesMutualInformationImageToImageMetric.h.


The documentation for this class was generated from the following file: