ITK  4.9.0
Insight Segmentation and Registration Toolkit
itkMattesMutualInformationImageToImageMetric.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright Insight Software Consortium
4  *
5  * Licensed under the Apache License, Version 2.0 (the "License");
6  * you may not use this file except in compliance with the License.
7  * You may obtain a copy of the License at
8  *
9  * http://www.apache.org/licenses/LICENSE-2.0.txt
10  *
11  * Unless required by applicable law or agreed to in writing, software
12  * distributed under the License is distributed on an "AS IS" BASIS,
13  * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
14  * See the License for the specific language governing permissions and
15  * limitations under the License.
16  *
17  *=========================================================================*/
18 #ifndef itkMattesMutualInformationImageToImageMetric_h
19 #define itkMattesMutualInformationImageToImageMetric_h
20 
21 #include "itkImageToImageMetric.h"
22 #include "itkPoint.h"
23 #include "itkIndex.h"
25 #include "itkArray2D.h"
26 
27 #include "itkSimpleFastMutexLock.h"
28 
29 
30 namespace itk
31 {
111 template <typename TFixedImage, typename TMovingImage>
113  public ImageToImageMetric<TFixedImage, TMovingImage>
114 {
115 public:
116 
122 
124  itkNewMacro(Self);
125 
129 
131  typedef typename Superclass::TransformType TransformType;
145 
151 
153 
155  itkStaticConstMacro(MovingImageDimension, unsigned int,
156  MovingImageType::ImageDimension);
157 
165  virtual void Initialize(void) throw ( ExceptionObject ) ITK_OVERRIDE;
166 
168  MeasureType GetValue(const ParametersType & parameters) const ITK_OVERRIDE;
169 
171  void GetDerivative(const ParametersType & parameters, DerivativeType & Derivative) const ITK_OVERRIDE;
172 
174  void GetValueAndDerivative(const ParametersType & parameters, MeasureType & Value, DerivativeType & Derivative) const ITK_OVERRIDE;
175 
182  itkSetClampMacro( NumberOfHistogramBins, SizeValueType,
183  5, NumericTraits<SizeValueType>::max() );
184  itkGetConstReferenceMacro(NumberOfHistogramBins, SizeValueType);
186 
211  itkSetMacro(UseExplicitPDFDerivatives, bool);
212  itkGetConstReferenceMacro(UseExplicitPDFDerivatives, bool);
213  itkBooleanMacro(UseExplicitPDFDerivatives);
215 
217  typedef double PDFValueType; //NOTE: floating point precision is not as stable. Double precision proves faster and more robust in real-world testing.
218 
220  typedef Image<PDFValueType, 2> JointPDFType;
221  typedef Image<PDFValueType, 3> JointPDFDerivativesType;
222 
227  const typename JointPDFType::Pointer GetJointPDF () const
228  {
229  if( this->m_MMIMetricPerThreadVariables == ITK_NULLPTR )
230  {
231  return JointPDFType::Pointer(ITK_NULLPTR);
232  }
233  return this->m_MMIMetricPerThreadVariables[0].JointPDF;
234  }
236 
244  {
245  if( this->m_MMIMetricPerThreadVariables == ITK_NULLPTR )
246  {
247  return JointPDFDerivativesType::Pointer(ITK_NULLPTR);
248  }
249  return this->m_MMIMetricPerThreadVariables[0].JointPDFDerivatives;
250  }
252 
253 protected:
254 
257  void PrintSelf(std::ostream & os, Indent indent) const ITK_OVERRIDE;
258 
259 private:
260  MattesMutualInformationImageToImageMetric(const Self &) ITK_DELETE_FUNCTION;
261  void operator=(const Self &) ITK_DELETE_FUNCTION;
262 
263  typedef JointPDFType::IndexType JointPDFIndexType;
264  typedef JointPDFType::PixelType JointPDFValueType;
265  typedef JointPDFType::RegionType JointPDFRegionType;
266  typedef JointPDFType::SizeType JointPDFSizeType;
271 
275 
278 
280  void ComputePDFDerivatives(ThreadIdType threadId, unsigned int sampleNumber, int movingImageParzenWindowIndex,
282  & movingImageGradientValue,
283  PDFValueType cubicBSplineDerivativeValue) const;
284 
285  virtual void GetValueThreadPreProcess(ThreadIdType threadId, bool withinSampleThread) const ITK_OVERRIDE;
286  virtual void GetValueThreadPostProcess(ThreadIdType threadId, bool withinSampleThread) const ITK_OVERRIDE;
287  //NOTE: The signature in base class requires that movingImageValue is of type double
288  virtual bool GetValueThreadProcessSample(ThreadIdType threadId, SizeValueType fixedImageSample,
289  const MovingImagePointType & mappedPoint,
290  double movingImageValue) const ITK_OVERRIDE;
291 
292  virtual void GetValueAndDerivativeThreadPreProcess( ThreadIdType threadId, bool withinSampleThread) const ITK_OVERRIDE;
293  virtual void GetValueAndDerivativeThreadPostProcess( ThreadIdType threadId, bool withinSampleThread) const ITK_OVERRIDE;
294  //NOTE: The signature in base class requires that movingImageValue is of type double
295  virtual bool GetValueAndDerivativeThreadProcessSample(ThreadIdType threadId, SizeValueType fixedImageSample,
296  const MovingImagePointType & mappedPoint,
297  double movingImageValue, const ImageDerivativesType &
298  movingImageGradientValue) const ITK_OVERRIDE;
299 
310 
312  typename CubicBSplineFunctionType::Pointer m_CubicBSplineKernel;
313  typename CubicBSplineDerivativeFunctionType::Pointer m_CubicBSplineDerivativeKernel;
314 
317  typedef Array2D<PRatioType> PRatioArrayType;
318 
319  mutable PRatioArrayType m_PRatioArray;
320 
322  typedef std::vector< PDFValueType > MarginalPDFType;
323  mutable MarginalPDFType m_MovingImageMarginalPDF;
324 
326  {
329 
331 
334 
338 
340 
341  MarginalPDFType FixedImageMarginalPDF;
342  };
343 
344 #if !defined(ITK_WRAPPING_PARSER)
345  itkPadStruct( ITK_CACHE_LINE_ALIGNMENT, MMIMetricPerThreadStruct,
346  PaddedMMIMetricPerThreadStruct);
347  itkAlignedTypedef( ITK_CACHE_LINE_ALIGNMENT, PaddedMMIMetricPerThreadStruct,
348  AlignedMMIMetricPerThreadStruct );
349  // Due to a bug in older version of Visual Studio where std::vector resize
350  // uses a value instead of a const reference, this must be a pointer.
351  // See
352  // http://thetweaker.wordpress.com/2010/05/05/stdvector-of-aligned-elements/
353  // http://connect.microsoft.com/VisualStudio/feedback/details/692988
354  mutable AlignedMMIMetricPerThreadStruct * m_MMIMetricPerThreadVariables;
355 #endif
356 
359 };
360 } // end namespace itk
361 
362 #ifndef ITK_MANUAL_INSTANTIATION
363 #include "itkMattesMutualInformationImageToImageMetric.hxx"
364 #endif
365 
366 #endif
Array class with size defined at construction time.
Definition: itkArray.h:50
itkAlignedTypedef(ITK_CACHE_LINE_ALIGNMENT, PaddedMMIMetricPerThreadStruct, AlignedMMIMetricPerThreadStruct)
Light weight base class for most itk classes.
const JointPDFDerivativesType::Pointer GetJointPDFDerivatives() const
BSplineTransformType::WeightsType BSplineTransformWeightsType
SmartPointer< Self > Pointer
Definition: itkImage.h:81
void ComputePDFDerivatives(ThreadIdType threadId, unsigned int sampleNumber, int movingImageParzenWindowIndex, const ImageDerivativesType &movingImageGradientValue, PDFValueType cubicBSplineDerivativeValue) const
Superclass::ParametersValueType CoordinateRepresentationType
Represent the size (bounds) of a n-dimensional image.
Definition: itkSize.h:52
signed long OffsetValueType
Definition: itkIntTypes.h:154
void PrintSelf(std::ostream &os, Indent indent) const override
TransformType::Pointer TransformPointer
CovariantVector< double, itkGetStaticConstMacro(MovingImageDimension) > ImageDerivativesType
virtual bool GetValueThreadProcessSample(ThreadIdType threadId, SizeValueType fixedImageSample, const MovingImagePointType &mappedPoint, double movingImageValue) const override
Superclass::DerivativeType DerivativeType
MovingImageType::ConstPointer MovingImageConstPointer
TransformType::OutputPointType MovingImagePointType
CubicBSplineDerivativeFunctionType::Pointer m_CubicBSplineDerivativeKernel
virtual void GetValueAndDerivativeThreadPostProcess(ThreadIdType threadId, bool withinSampleThread) const override
unsigned long SizeValueType
Definition: itkIntTypes.h:143
MeasureType GetValue(const ParametersType &parameters) const override
Derivative of a BSpline kernel used for density estimation and nonparameteric regression.
BSpline kernel used for density estimation and nonparameteric regression.
Transform points and vectors from an input space to an output space.
Definition: itkTransform.h:82
Array2D class representing a 2D array with size defined at construction time.
Definition: itkArray2D.h:45
virtual void GetValueThreadPreProcess(ThreadIdType threadId, bool withinSampleThread) const override
Superclass::ParametersType ParametersType
virtual bool GetValueAndDerivativeThreadProcessSample(ThreadIdType threadId, SizeValueType fixedImageSample, const MovingImagePointType &mappedPoint, double movingImageValue, const ImageDerivativesType &movingImageGradientValue) const override
Standard exception handling object.
BSplineTransformIndexArrayType::ValueType IndexValueType
void ComputeFixedImageParzenWindowIndices(FixedImageSampleContainer &samples)
unsigned int ThreadIdType
Definition: itkIntTypes.h:159
std::vector< FixedImageSamplePoint > FixedImageSampleContainer
TransformType::JacobianType TransformJacobianType
virtual void GetValueAndDerivativeThreadPreProcess(ThreadIdType threadId, bool withinSampleThread) const override
Computes the mutual information between two images to be registered using the method of Mattes et al...
BSplineTransformWeightsType::ValueType WeightsValueType
InterpolateImageFunction< MovingImageType, CoordinateRepresentationType > InterpolatorType
Control indentation during Print() invocation.
Definition: itkIndent.h:49
FixedImageType::ConstPointer FixedImageConstPointer
Define additional traits for native types such as int or float.
Superclass::MeasureType MeasureType
Computes similarity between regions of two images.
A templated class holding a geometric point in n-Dimensional space.
Definition: itkPoint.h:52
A templated class holding a n-Dimensional covariant vector.
virtual void GetValueThreadPostProcess(ThreadIdType threadId, bool withinSampleThread) const override
void GetDerivative(const ParametersType &parameters, DerivativeType &Derivative) const override
Templated n-dimensional image class.
Definition: itkImage.h:75
itkPadStruct(ITK_CACHE_LINE_ALIGNMENT, MMIMetricPerThreadStruct, PaddedMMIMetricPerThreadStruct)
void GetValueAndDerivative(const ParametersType &parameters, MeasureType &Value, DerivativeType &Derivative) const override
BSplineTransformType::ParameterIndexArrayType BSplineTransformIndexArrayType