ITK  6.0.0
Insight Toolkit
itkMattesMutualInformationImageToImageMetric.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright NumFOCUS
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  * https://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 <memory> // For unique_ptr.
28 #include <mutex>
29 
30 
31 namespace itk
32 {
116 template <typename TFixedImage, typename TMovingImage>
118  : public ImageToImageMetric<TFixedImage, TMovingImage>
119 {
120 public:
121  ITK_DISALLOW_COPY_AND_MOVE(MattesMutualInformationImageToImageMetric);
122 
128 
130  itkNewMacro(Self);
131 
133  itkOverrideGetNameOfClassMacro(MattesMutualInformationImageToImageMetric);
134 
136  using typename Superclass::TransformType;
137  using typename Superclass::TransformPointer;
138  using typename Superclass::TransformJacobianType;
139  using typename Superclass::InterpolatorType;
140  using typename Superclass::MeasureType;
141  using typename Superclass::DerivativeType;
142  using typename Superclass::ParametersType;
143  using typename Superclass::FixedImageType;
144  using typename Superclass::MovingImageType;
145  using typename Superclass::MovingImagePointType;
146  using typename Superclass::FixedImageConstPointer;
147  using typename Superclass::MovingImageConstPointer;
148  using typename Superclass::BSplineTransformWeightsType;
149  using typename Superclass::BSplineTransformIndexArrayType;
150 
151  using typename Superclass::CoordinateRepresentationType;
152  using typename Superclass::FixedImageSampleContainer;
153  using typename Superclass::ImageDerivativesType;
154  using typename Superclass::WeightsValueType;
155  using typename Superclass::IndexValueType;
156 
158 
160  static constexpr unsigned int MovingImageDimension = MovingImageType::ImageDimension;
161 
169  void
170  Initialize() override;
171 
174  GetValue(const ParametersType & parameters) const override;
175 
177  void
178  GetDerivative(const ParametersType & parameters, DerivativeType & derivative) const override;
179 
181  void
182  GetValueAndDerivative(const ParametersType & parameters,
183  MeasureType & value,
184  DerivativeType & derivative) const override;
185 
192  itkSetClampMacro(NumberOfHistogramBins, SizeValueType, 5, NumericTraits<SizeValueType>::max());
193  itkGetConstReferenceMacro(NumberOfHistogramBins, SizeValueType);
220  itkSetMacro(UseExplicitPDFDerivatives, bool);
221  itkGetConstReferenceMacro(UseExplicitPDFDerivatives, bool);
222  itkBooleanMacro(UseExplicitPDFDerivatives);
226  using PDFValueType = double; // NOTE: floating point precision is not as stable. Double precision proves faster and
227  // more robust in real-world testing.
228 
232 
237  const typename JointPDFType::Pointer
238  GetJointPDF() const
239  {
240  if (this->m_MMIMetricPerThreadVariables == nullptr)
241  {
242  return JointPDFType::Pointer(nullptr);
243  }
244  return this->m_MMIMetricPerThreadVariables[0].JointPDF;
245  }
254  const typename JointPDFDerivativesType::Pointer
256  {
257  if (this->m_MMIMetricPerThreadVariables == nullptr)
258  {
259  return JointPDFDerivativesType::Pointer(nullptr);
260  }
261  return this->m_MMIMetricPerThreadVariables[0].JointPDFDerivatives;
262  }
265 protected:
267  ~MattesMutualInformationImageToImageMetric() override = default;
268  void
269  PrintSelf(std::ostream & os, Indent indent) const override;
270 
271 private:
280 
284 
286  void
287  CommonGetValueProcessing() const;
288 
290  void
291  ComputeFixedImageParzenWindowIndices(FixedImageSampleContainer & samples);
292 
294  void
295  ComputePDFDerivatives(ThreadIdType threadId,
296  unsigned int sampleNumber,
297  int pdfMovingIndex,
298  const ImageDerivativesType & movingImageGradientValue,
299  PDFValueType cubicBSplineDerivativeValue) const;
300 
301  void
302  GetValueThreadPreProcess(ThreadIdType threadId, bool withinSampleThread) const override;
303  void
304  GetValueThreadPostProcess(ThreadIdType threadId, bool withinSampleThread) const override;
305  // NOTE: The signature in base class requires that movingImageValue is of type double
306  bool
307  GetValueThreadProcessSample(ThreadIdType threadId,
308  SizeValueType fixedImageSample,
309  const MovingImagePointType & mappedPoint,
310  double movingImageValue) const override;
311 
312  void
313  GetValueAndDerivativeThreadPreProcess(ThreadIdType threadId, bool withinSampleThread) const override;
314  void
315  GetValueAndDerivativeThreadPostProcess(ThreadIdType threadId, bool withinSampleThread) const override;
316  // NOTE: The signature in base class requires that movingImageValue is of type double
317  bool
318  GetValueAndDerivativeThreadProcessSample(ThreadIdType threadId,
319  SizeValueType fixedImageSample,
320  const MovingImagePointType & mappedPoint,
321  double movingImageValue,
322  const ImageDerivativesType & movingImageGradientValue) const override;
323 
325  SizeValueType m_NumberOfHistogramBins{ 50 };
326  PDFValueType m_MovingImageNormalizedMin{ 0.0 };
327  PDFValueType m_FixedImageNormalizedMin{ 0.0 };
328  PDFValueType m_FixedImageTrueMin{ 0.0 };
329  PDFValueType m_FixedImageTrueMax{ 0.0 };
330  PDFValueType m_MovingImageTrueMin{ 0.0 };
331  PDFValueType m_MovingImageTrueMax{ 0.0 };
332  PDFValueType m_FixedImageBinSize{ 0.0 };
333  PDFValueType m_MovingImageBinSize{ 0.0 };
334 
338 
339  mutable PRatioArrayType m_PRatioArray{};
340 
342  using MarginalPDFType = std::vector<PDFValueType>;
343  mutable MarginalPDFType m_MovingImageMarginalPDF{};
344 
346  {
349 
351 
354 
358 
360 
362  };
363 
364 #if !defined(ITK_WRAPPING_PARSER)
365  itkPadStruct(ITK_CACHE_LINE_ALIGNMENT, MMIMetricPerThreadStruct, PaddedMMIMetricPerThreadStruct);
366  itkAlignedTypedef(ITK_CACHE_LINE_ALIGNMENT, PaddedMMIMetricPerThreadStruct, AlignedMMIMetricPerThreadStruct);
367  // Due to a bug in older version of Visual Studio where std::vector resize
368  // uses a value instead of a const reference, this must be a pointer.
369  // See
370  // https://thetweaker.wordpress.com/2010/05/05/stdvector-of-aligned-elements/
371  // https://connect.microsoft.com/VisualStudio/feedback/details/692988
372  std::unique_ptr<AlignedMMIMetricPerThreadStruct[]> m_MMIMetricPerThreadVariables;
373 #endif
374 
375  bool m_UseExplicitPDFDerivatives{ true };
376  mutable bool m_ImplicitDerivativesSecondPass{ false };
377 };
378 } // end namespace itk
379 
380 #ifndef ITK_MANUAL_INSTANTIATION
381 # include "itkMattesMutualInformationImageToImageMetric.hxx"
382 #endif
383 
384 #endif
Pointer
SmartPointer< Self > Pointer
Definition: itkAddImageFilter.h:93
itk::ImageToImageMetric
Computes similarity between regions of two images.
Definition: itkImageToImageMetric.h:54
itk::BSplineKernelFunction
BSpline kernel used for density estimation and nonparametric regression.
Definition: itkBSplineKernelFunction.h:43
itk::SingleValuedCostFunction::MeasureType
double MeasureType
Definition: itkSingleValuedCostFunction.h:50
itk::OptimizerParameters< TInternalComputationValueType >
itk::MattesMutualInformationImageToImageMetric::MMIMetricPerThreadStruct::JointPDF
JointPDFType::Pointer JointPDF
Definition: itkMattesMutualInformationImageToImageMetric.h:356
itk::Index< VImageDimension >
itk::MattesMutualInformationImageToImageMetric::GetJointPDFDerivatives
const JointPDFDerivativesType::Pointer GetJointPDFDerivatives() const
Definition: itkMattesMutualInformationImageToImageMetric.h:255
itk::MattesMutualInformationImageToImageMetric::MarginalPDFType
std::vector< PDFValueType > MarginalPDFType
Definition: itkMattesMutualInformationImageToImageMetric.h:342
itk::MattesMutualInformationImageToImageMetric::JointPDFValueType
JointPDFType::PixelType JointPDFValueType
Definition: itkMattesMutualInformationImageToImageMetric.h:273
itk::Size< VImageDimension >
itk::MattesMutualInformationImageToImageMetric::MMIMetricPerThreadStruct::JointPDFStartBin
int JointPDFStartBin
Definition: itkMattesMutualInformationImageToImageMetric.h:347
itk::MattesMutualInformationImageToImageMetric::MMIMetricPerThreadStruct::FixedImageMarginalPDF
MarginalPDFType FixedImageMarginalPDF
Definition: itkMattesMutualInformationImageToImageMetric.h:361
itkPoint.h
itk::ImageRegion
An image region represents a structured region of data.
Definition: itkImageRegion.h:80
itk::ImageToImageMetric::MovingImagePointType
typename TransformType::OutputPointType MovingImagePointType
Definition: itkImageToImageMetric.h:100
itk::MattesMutualInformationImageToImageMetric::m_MMIMetricPerThreadVariables
std::unique_ptr< AlignedMMIMetricPerThreadStruct[]> m_MMIMetricPerThreadVariables
Definition: itkMattesMutualInformationImageToImageMetric.h:372
itk::GTest::TypedefsAndConstructors::Dimension2::SizeType
ImageBaseType::SizeType SizeType
Definition: itkGTestTypedefsAndConstructors.h:49
itk::MattesMutualInformationImageToImageMetric::OffsetValueType
typename FixedImageType::OffsetValueType OffsetValueType
Definition: itkMattesMutualInformationImageToImageMetric.h:157
itk::MattesMutualInformationImageToImageMetric::MMIMetricPerThreadStruct::JointPDFDerivatives
JointPDFDerivativesType::Pointer JointPDFDerivatives
Definition: itkMattesMutualInformationImageToImageMetric.h:357
itk::MattesMutualInformationImageToImageMetric::MMIMetricPerThreadStruct::JointPDFSum
PDFValueType JointPDFSum
Definition: itkMattesMutualInformationImageToImageMetric.h:350
itk::SmartPointer< Self >
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::IndexValueType
long IndexValueType
Definition: itkIntTypes.h:93
itk::ThreadIdType
unsigned int ThreadIdType
Definition: itkIntTypes.h:102
itk::GTest::TypedefsAndConstructors::Dimension2::IndexType
ImageBaseType::IndexType IndexType
Definition: itkGTestTypedefsAndConstructors.h:50
itk::LightObject
Light weight base class for most itk classes.
Definition: itkLightObject.h:55
itkImageToImageMetric.h
itk::GTest::TypedefsAndConstructors::Dimension2::RegionType
ImageBaseType::RegionType RegionType
Definition: itkGTestTypedefsAndConstructors.h:54
itk::OffsetValueType
long OffsetValueType
Definition: itkIntTypes.h:97
itk::ImageToImageMetric::FixedImageSampleContainer
std::vector< FixedImageSamplePoint > FixedImageSampleContainer
Definition: itkImageToImageMetric.h:409
itkIndex.h
itk::MattesMutualInformationImageToImageMetric::MMIMetricPerThreadStruct::JointPDFEndBin
int JointPDFEndBin
Definition: itkMattesMutualInformationImageToImageMetric.h:348
itk::Image::PixelType
TPixel PixelType
Definition: itkImage.h:108
itk::NumericTraits
Define additional traits for native types such as int or float.
Definition: itkNumericTraits.h:60
itk::CovariantVector
A templated class holding a n-Dimensional covariant vector.
Definition: itkCovariantVector.h:70
itk::MattesMutualInformationImageToImageMetric::MMIMetricPerThreadStruct::MetricDerivative
DerivativeType MetricDerivative
Definition: itkMattesMutualInformationImageToImageMetric.h:353
itk::MattesMutualInformationImageToImageMetric::JointPDFDerivativesValueType
JointPDFDerivativesType::PixelType JointPDFDerivativesValueType
Definition: itkMattesMutualInformationImageToImageMetric.h:277
itkArray2D.h
itk::MattesMutualInformationImageToImageMetric::MMIMetricPerThreadStruct::Jacobian
TransformType::JacobianType Jacobian
Definition: itkMattesMutualInformationImageToImageMetric.h:359
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnatomicalOrientation.h:29
itk::BSplineDerivativeKernelFunction
Derivative of a BSpline kernel used for density estimation and nonparametric regression.
Definition: itkBSplineDerivativeKernelFunction.h:43
itk::Array
Array class with size defined at construction time.
Definition: itkArray.h:47
itk::MattesMutualInformationImageToImageMetric::MMIMetricPerThreadStruct
Definition: itkMattesMutualInformationImageToImageMetric.h:345
itk::Image
Templated n-dimensional image class.
Definition: itkImage.h:88
itk::MattesMutualInformationImageToImageMetric::PRatioType
PDFValueType PRatioType
Definition: itkMattesMutualInformationImageToImageMetric.h:336
itk::Array2D< PRatioType >
itk::MattesMutualInformationImageToImageMetric::GetJointPDF
const JointPDFType::Pointer GetJointPDF() const
Definition: itkMattesMutualInformationImageToImageMetric.h:238
itk::MattesMutualInformationImageToImageMetric::PDFValueType
double PDFValueType
Definition: itkMattesMutualInformationImageToImageMetric.h:226
itk::SizeValueType
unsigned long SizeValueType
Definition: itkIntTypes.h:86
itkBSplineDerivativeKernelFunction.h
itk::MattesMutualInformationImageToImageMetric
Computes the mutual information between two images to be registered using the method of Mattes et al.
Definition: itkMattesMutualInformationImageToImageMetric.h:117