/*========================================================================= * * Copyright NumFOCUS * * Licensed under the Apache License, Version 2.0 (the "License"); * you may not use this file except in compliance with the License. * You may obtain a copy of the License at * * https://www.apache.org/licenses/LICENSE-2.0.txt * * Unless required by applicable law or agreed to in writing, software * distributed under the License is distributed on an "AS IS" BASIS, * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied. * See the License for the specific language governing permissions and * limitations under the License. * *=========================================================================*/ #ifndef itkJointHistogramMutualInformationImageToImageMetricv4_h #define itkJointHistogramMutualInformationImageToImageMetricv4_h #include "itkImageToImageMetricv4.h" #include "itkImage.h" #include "itkBSplineDerivativeKernelFunction.h" #include "itkJointHistogramMutualInformationComputeJointPDFThreader.h" #include "itkJointHistogramMutualInformationGetValueAndDerivativeThreader.h" namespace itk { /** \class JointHistogramMutualInformationImageToImageMetricv4 * \brief Computes the mutual information between two images to be * registered using the method referenced below. * * References: * [1] "Optimization of Mutual Information for MultiResolution Image * Registration" * P. Thevenaz and M. Unser * IEEE Transactions in Image Processing, 9(12) December 2000. * * \ingroup ITKMetricsv4 */ template > class ITK_TEMPLATE_EXPORT JointHistogramMutualInformationImageToImageMetricv4 : public ImageToImageMetricv4 { public: ITK_DISALLOW_COPY_AND_MOVE(JointHistogramMutualInformationImageToImageMetricv4); /** Standard class type aliases. */ using Self = JointHistogramMutualInformationImageToImageMetricv4; using Superclass = ImageToImageMetricv4; using Pointer = SmartPointer; using ConstPointer = SmartPointer; /** Method for creation through the object factory. */ itkNewMacro(Self); /** \see LightObject::GetNameOfClass() */ itkOverrideGetNameOfClassMacro(JointHistogramMutualInformationImageToImageMetricv4); /** Type used for representing parameter values */ using typename Superclass::CoordinateRepresentationType; /** Type used internally for computations */ /** It should be possible to derive the internal computation type from the class object. */ using InternalComputationValueType = TInternalComputationValueType; /** Type of the parameters. */ using typename Superclass::ParametersType; using typename Superclass::ParametersValueType; using typename Superclass::NumberOfParametersType; /** Superclass type alias */ using typename Superclass::MeasureType; using typename Superclass::DerivativeType; using typename Superclass::FixedImagePointType; using typename Superclass::FixedImagePixelType; using FixedImageGradientType = typename Superclass::FixedGradientPixelType; using typename Superclass::MovingImagePointType; using typename Superclass::MovingImagePixelType; using MovingImageGradientType = typename Superclass::MovingGradientPixelType; using FixedTransformJacobianType = typename Superclass::FixedTransformType::JacobianType; using MovingTransformJacobianType = typename Superclass::MovingTransformType::JacobianType; using VirtualImageType = typename Superclass::VirtualImageType; using typename Superclass::VirtualIndexType; using typename Superclass::VirtualPointType; using typename Superclass::VirtualPointSetType; /* Image dimension accessors */ static constexpr typename TVirtualImage::ImageDimensionType VirtualImageDimension = TVirtualImage::ImageDimension; static constexpr typename TMovingImage::ImageDimensionType MovingImageDimension = TMovingImage::ImageDimension; /** Value type of the PDF */ using PDFValueType = TInternalComputationValueType; /** Typedef for the joint PDF and marginal PDF are stored as ITK Images. */ using MarginalPDFType = Image; using MarginalPDFIndexType = typename MarginalPDFType::IndexType; using MarginalPDFPointType = typename MarginalPDFType::PointType; using JointPDFType = Image; using JointPDFIndexType = typename JointPDFType::IndexType; using JointPDFPointType = typename JointPDFType::PointType; using JointPDFIndexValueType = typename JointPDFType::IndexValueType; /** Get the JointPDF. Valid after GetValueAndDerivative has been called. */ itkGetModifiableObjectMacro(JointPDF, JointPDFType); // Declare the type for the derivative calculation using JPDFGradientFilterType = itk::GradientRecursiveGaussianImageFilter; using JPDFGradientImageType = typename JPDFGradientFilterType::OutputImageType; using JPDFGradientImagePointer = typename JPDFGradientImageType::Pointer; using MarginalGradientFilterType = itk::GradientRecursiveGaussianImageFilter; using MarginalGradientImageType = typename MarginalGradientFilterType::OutputImageType; using MarginalGradientImagePointer = typename MarginalGradientImageType::Pointer; /** pdf interpolator */ using JointPDFInterpolatorType = LinearInterpolateImageFunction; using JointPDFInterpolatorPointer = typename JointPDFInterpolatorType::Pointer; using MarginalPDFInterpolatorType = LinearInterpolateImageFunction; using MarginalPDFInterpolatorPointer = typename MarginalPDFInterpolatorType::Pointer; /** Joint PDF types */ using JointPDFValueType = typename JointPDFType::PixelType; using JointPDFRegionType = typename JointPDFType::RegionType; using JointPDFSizeType = typename JointPDFType::SizeType; using JointPDFSpacingType = typename JointPDFType::SpacingType; /** Get/Set the number of histogram bins */ itkSetClampMacro(NumberOfHistogramBins, SizeValueType, 5, NumericTraits::max()); itkGetConstReferenceMacro(NumberOfHistogramBins, SizeValueType); /** Get/Set option to smooth the joint pdf after it's updated */ itkSetMacro(VarianceForJointPDFSmoothing, TInternalComputationValueType); itkGetMacro(VarianceForJointPDFSmoothing, TInternalComputationValueType); /** Initialize the metric. Make sure all essential inputs are plugged in. */ void Initialize() override; MeasureType GetValue() const override; protected: JointHistogramMutualInformationImageToImageMetricv4(); ~JointHistogramMutualInformationImageToImageMetricv4() override = default; /** Update the histograms for use in GetValueAndDerivative * Results are returned in \c value and \c derivative. */ void InitializeForIteration() const override; /** Compute the metric value. For internal use. */ MeasureType ComputeValue() const; /** Compute the point location with the JointPDF image. Returns false if the * point is not inside the image. */ inline void ComputeJointPDFPoint(const FixedImagePixelType fixedImageValue, const MovingImagePixelType movingImageValue, JointPDFPointType & jointPDFpoint) const; friend class JointHistogramMutualInformationComputeJointPDFThreaderBase< ThreadedImageRegionPartitioner, Self>; friend class JointHistogramMutualInformationComputeJointPDFThreaderBase; friend class JointHistogramMutualInformationComputeJointPDFThreader< ThreadedImageRegionPartitioner, Self>; friend class JointHistogramMutualInformationComputeJointPDFThreader; using JointHistogramMutualInformationDenseComputeJointPDFThreaderType = JointHistogramMutualInformationComputeJointPDFThreader, Self>; using JointHistogramMutualInformationSparseComputeJointPDFThreaderType = JointHistogramMutualInformationComputeJointPDFThreader; typename JointHistogramMutualInformationDenseComputeJointPDFThreaderType::Pointer m_JointHistogramMutualInformationDenseComputeJointPDFThreader; typename JointHistogramMutualInformationSparseComputeJointPDFThreaderType::Pointer m_JointHistogramMutualInformationSparseComputeJointPDFThreader; friend class JointHistogramMutualInformationGetValueAndDerivativeThreader< ThreadedImageRegionPartitioner, Superclass, Self>; friend class JointHistogramMutualInformationGetValueAndDerivativeThreader; using JointHistogramMutualInformationDenseGetValueAndDerivativeThreaderType = JointHistogramMutualInformationGetValueAndDerivativeThreader< ThreadedImageRegionPartitioner, Superclass, Self>; using JointHistogramMutualInformationSparseGetValueAndDerivativeThreaderType = JointHistogramMutualInformationGetValueAndDerivativeThreader; /** Standard PrintSelf method. */ void PrintSelf(std::ostream & os, Indent indent) const override; /** Count of the number of valid histogram points. */ SizeValueType m_JointHistogramTotalCount{ 0 }; private: /** The fixed image marginal PDF */ typename MarginalPDFType::Pointer m_FixedImageMarginalPDF{}; /** The moving image marginal PDF. */ typename MarginalPDFType::Pointer m_MovingImageMarginalPDF{}; /** The joint PDF and PDF derivatives. */ mutable typename JointPDFType::Pointer m_JointPDF{}; /** Flag to control smoothing of joint pdf */ TInternalComputationValueType m_VarianceForJointPDFSmoothing{}; /** Variables to define the marginal and joint histograms. */ SizeValueType m_NumberOfHistogramBins{}; TInternalComputationValueType m_FixedImageTrueMin{}; TInternalComputationValueType m_FixedImageTrueMax{}; TInternalComputationValueType m_MovingImageTrueMin{}; TInternalComputationValueType m_MovingImageTrueMax{}; TInternalComputationValueType m_FixedImageBinSize{}; TInternalComputationValueType m_MovingImageBinSize{}; TInternalComputationValueType m_JointPDFSum{}; JointPDFSpacingType m_JointPDFSpacing{}; TInternalComputationValueType m_Log2{}; JointPDFIndexValueType m_Padding{}; }; } // end namespace itk #ifndef ITK_MANUAL_INSTANTIATION # include "itkJointHistogramMutualInformationImageToImageMetricv4.hxx" #endif #endif