ITK  4.13.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;
132  typedef typename Superclass::TransformPointer TransformPointer;
133  typedef typename Superclass::TransformJacobianType TransformJacobianType;
134  typedef typename Superclass::InterpolatorType InterpolatorType;
135  typedef typename Superclass::MeasureType MeasureType;
136  typedef typename Superclass::DerivativeType DerivativeType;
137  typedef typename Superclass::ParametersType ParametersType;
138  typedef typename Superclass::FixedImageType FixedImageType;
139  typedef typename Superclass::MovingImageType MovingImageType;
140  typedef typename Superclass::MovingImagePointType MovingImagePointType;
141  typedef typename Superclass::FixedImageConstPointer FixedImageConstPointer;
142  typedef typename Superclass::MovingImageConstPointer MovingImageConstPointer;
143  typedef typename Superclass::BSplineTransformWeightsType BSplineTransformWeightsType;
144  typedef typename Superclass::BSplineTransformIndexArrayType BSplineTransformIndexArrayType;
145 
146  typedef typename Superclass::CoordinateRepresentationType CoordinateRepresentationType;
147  typedef typename Superclass::FixedImageSampleContainer FixedImageSampleContainer;
148  typedef typename Superclass::ImageDerivativesType ImageDerivativesType;
149  typedef typename Superclass::WeightsValueType WeightsValueType;
151 
153 
155  itkStaticConstMacro(MovingImageDimension, unsigned int,
156  MovingImageType::ImageDimension);
157 
165  virtual void Initialize(void) 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 
256  virtual ~MattesMutualInformationImageToImageMetric() ITK_OVERRIDE;
257  void PrintSelf(std::ostream & os, Indent indent) const ITK_OVERRIDE;
258 
259 private:
260  ITK_DISALLOW_COPY_AND_ASSIGN(MattesMutualInformationImageToImageMetric);
261 
263  typedef JointPDFType::PixelType JointPDFValueType;
264  typedef JointPDFType::RegionType JointPDFRegionType;
270 
274 
276  void ComputeFixedImageParzenWindowIndices( FixedImageSampleContainer & samples);
277 
279  void ComputePDFDerivatives(ThreadIdType threadId, unsigned int sampleNumber, int movingImageParzenWindowIndex,
281  & movingImageGradientValue,
282  PDFValueType cubicBSplineDerivativeValue) const;
283 
284  virtual void GetValueThreadPreProcess(ThreadIdType threadId, bool withinSampleThread) const ITK_OVERRIDE;
285  virtual void GetValueThreadPostProcess(ThreadIdType threadId, bool withinSampleThread) const ITK_OVERRIDE;
286  //NOTE: The signature in base class requires that movingImageValue is of type double
287  virtual bool GetValueThreadProcessSample(ThreadIdType threadId, SizeValueType fixedImageSample,
288  const MovingImagePointType & mappedPoint,
289  double movingImageValue) const ITK_OVERRIDE;
290 
291  virtual void GetValueAndDerivativeThreadPreProcess( ThreadIdType threadId, bool withinSampleThread) const ITK_OVERRIDE;
292  virtual void GetValueAndDerivativeThreadPostProcess( ThreadIdType threadId, bool withinSampleThread) const ITK_OVERRIDE;
293  //NOTE: The signature in base class requires that movingImageValue is of type double
294  virtual bool GetValueAndDerivativeThreadProcessSample(ThreadIdType threadId, SizeValueType fixedImageSample,
295  const MovingImagePointType & mappedPoint,
296  double movingImageValue, const ImageDerivativesType &
297  movingImageGradientValue) const ITK_OVERRIDE;
298 
300  SizeValueType m_NumberOfHistogramBins;
301  PDFValueType m_MovingImageNormalizedMin;
302  PDFValueType m_FixedImageNormalizedMin;
303  PDFValueType m_FixedImageTrueMin;
304  PDFValueType m_FixedImageTrueMax;
305  PDFValueType m_MovingImageTrueMin;
306  PDFValueType m_MovingImageTrueMax;
307  PDFValueType m_FixedImageBinSize;
308  PDFValueType m_MovingImageBinSize;
309 
311  typename CubicBSplineFunctionType::Pointer m_CubicBSplineKernel;
312  typename CubicBSplineDerivativeFunctionType::Pointer m_CubicBSplineDerivativeKernel;
313 
316  typedef Array2D<PRatioType> PRatioArrayType;
317 
318  mutable PRatioArrayType m_PRatioArray;
319 
321  typedef std::vector< PDFValueType > MarginalPDFType;
322  mutable MarginalPDFType m_MovingImageMarginalPDF;
323 
325  {
328 
330 
333 
337 
339 
340  MarginalPDFType FixedImageMarginalPDF;
341  };
342 
343 #if !defined(ITK_WRAPPING_PARSER)
344  itkPadStruct( ITK_CACHE_LINE_ALIGNMENT, MMIMetricPerThreadStruct,
345  PaddedMMIMetricPerThreadStruct);
346  itkAlignedTypedef( ITK_CACHE_LINE_ALIGNMENT, PaddedMMIMetricPerThreadStruct,
347  AlignedMMIMetricPerThreadStruct );
348  // Due to a bug in older version of Visual Studio where std::vector resize
349  // uses a value instead of a const reference, this must be a pointer.
350  // See
351  // http://thetweaker.wordpress.com/2010/05/05/stdvector-of-aligned-elements/
352  // http://connect.microsoft.com/VisualStudio/feedback/details/692988
353  mutable AlignedMMIMetricPerThreadStruct * m_MMIMetricPerThreadVariables;
354 #endif
355 
358 };
359 } // end namespace itk
360 
361 #ifndef ITK_MANUAL_INSTANTIATION
362 #include "itkMattesMutualInformationImageToImageMetric.hxx"
363 #endif
364 
365 #endif
Array class with size defined at construction time.
Definition: itkArray.h:50
Light weight base class for most itk classes.
const JointPDFDerivativesType::Pointer GetJointPDFDerivatives() const
Represent the size (bounds) of a n-dimensional image.
Definition: itkSize.h:52
signed long OffsetValueType
Definition: itkIntTypes.h:154
signed long IndexValueType
Definition: itkIntTypes.h:150
unsigned long SizeValueType
Definition: itkIntTypes.h:143
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
unsigned int ThreadIdType
Definition: itkIntTypes.h:159
std::vector< FixedImageSamplePoint > FixedImageSampleContainer
Computes the mutual information between two images to be registered using the method of Mattes et al...
Control indentation during Print() invocation.
Definition: itkIndent.h:49
Define additional traits for native types such as int or float.
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.
Templated n-dimensional image class.
Definition: itkImage.h:75