go home Home | Main Page | Modules | Namespace List | Class Hierarchy | Alphabetical List | Data Structures | File List | Namespace Members | Data Fields | Globals | Related Pages
itkParzenWindowMutualInformationImageToImageMetric.h
Go to the documentation of this file.
1/*=========================================================================
2 *
3 * Copyright UMC Utrecht and contributors
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 itkParzenWindowMutualInformationImageToImageMetric_h
19#define itkParzenWindowMutualInformationImageToImageMetric_h
20
22
23#include "itkArray2D.h"
24
25namespace itk
26{
27
74template <class TFixedImage, class TMovingImage>
76 : public ParzenWindowHistogramImageToImageMetric<TFixedImage, TMovingImage>
77{
78public:
80
84 using Pointer = SmartPointer<Self>;
85 using ConstPointer = SmartPointer<const Self>;
86
88 itkNewMacro(Self);
89
92
94 using typename Superclass::CoordinateRepresentationType;
95 using typename Superclass::MovingImageType;
96 using typename Superclass::MovingImagePixelType;
97 using typename Superclass::MovingImageConstPointer;
98 using typename Superclass::FixedImageType;
99 using typename Superclass::FixedImageConstPointer;
100 using typename Superclass::FixedImageRegionType;
101 using typename Superclass::TransformType;
102 using typename Superclass::TransformPointer;
103 using typename Superclass::InputPointType;
104 using typename Superclass::OutputPointType;
105 using typename Superclass::TransformParametersType;
106 using typename Superclass::TransformJacobianType;
108 using typename Superclass::InterpolatorType;
109 using typename Superclass::InterpolatorPointer;
110 using typename Superclass::RealType;
111 using typename Superclass::GradientPixelType;
112 using typename Superclass::GradientImageType;
113 using typename Superclass::GradientImagePointer;
114 using typename Superclass::FixedImageMaskType;
116 using typename Superclass::MovingImageMaskType;
118 using typename Superclass::MeasureType;
119 using typename Superclass::DerivativeType;
120 using typename Superclass::DerivativeValueType;
121 using typename Superclass::ParametersType;
122 using typename Superclass::FixedImagePixelType;
124 using typename Superclass::ImageSamplerType;
125 using typename Superclass::ImageSamplerPointer;
133 using typename Superclass::ThreadInfoType;
134
136 itkStaticConstMacro(FixedImageDimension, unsigned int, FixedImageType::ImageDimension);
137
139 itkStaticConstMacro(MovingImageDimension, unsigned int, MovingImageType::ImageDimension);
140
142 MeasureType
143 GetValue(const ParametersType & parameters) const override;
144
146 itkGetConstMacro(UseJacobianPreconditioning, bool);
147 itkSetMacro(UseJacobianPreconditioning, bool);
148
149protected:
152
155
159 using typename Superclass::FixedImageIndexType;
162 using typename Superclass::FixedImagePointType;
167 using typename Superclass::PDFValueType;
169 using typename Superclass::MarginalPDFType;
170 using typename Superclass::JointPDFType;
173 using typename Superclass::JointPDFIndexType;
174 using typename Superclass::JointPDFRegionType;
175 using typename Superclass::JointPDFSizeType;
180 using typename Superclass::KernelFunctionType;
182
190 void
191 GetValueAndAnalyticDerivative(const ParametersType & parameters,
192 MeasureType & value,
193 DerivativeType & derivative) const override;
194
204 virtual void
205 GetValueAndAnalyticDerivativeLowMemory(const ParametersType & parameters,
206 MeasureType & value,
207 DerivativeType & derivative) const;
208
214 void
215 GetValueAndFiniteDifferenceDerivative(const ParametersType & parameters,
216 MeasureType & value,
217 DerivativeType & derivative) const override;
218
220 virtual void
221 ComputeJacobianPreconditioner(const TransformJacobianType & jac,
222 const NonZeroJacobianIndicesType & nzji,
223 DerivativeType & preconditioner,
224 DerivativeType & divisor) const;
225
227 void
229
235 ParzenWindowMutualInformationMultiThreaderParameterType m_ParzenWindowMutualInformationThreaderParameters{};
236
238 void
240
242 void
243 AfterThreadedComputeDerivativeLowMemory(DerivativeType & derivative) const;
244
246 static ITK_THREAD_RETURN_FUNCTION_CALL_CONVENTION
248
250 void
252
253private:
256 using PRatioArrayType = vnl_matrix<PRatioType>;
257 mutable PRatioArrayType m_PRatioArray{};
258
260 bool m_UseJacobianPreconditioning{ false };
261
263 void
264 ComputeDerivativeLowMemorySingleThreaded(DerivativeType & derivative) const;
265
266 void
267 ComputeDerivativeLowMemory(DerivativeType & derivative) const;
268
270 void
271 UpdateDerivativeLowMemory(const RealType fixedImageValue,
272 const RealType movingImageValue,
273 const DerivativeType & imageJacobian,
274 const NonZeroJacobianIndicesType & nzji,
275 DerivativeType & derivative) const;
276
278 void
279 ComputeValueAndPRatioArray(double & MI) const;
280};
281
282} // end namespace itk
283
284#ifndef ITK_MANUAL_INSTANTIATION
285# include "itkParzenWindowMutualInformationImageToImageMetric.hxx"
286#endif
287
288#endif // end #ifndef itkParzenWindowMutualInformationImageToImageMetric_h
typename AdvancedTransformType::NumberOfParametersType NumberOfParametersType
This class is a base class for any image sampler.
Kernel used for density estimation and nonparameteric regression.
Base class for all ITK limiter function objects.
A base class for image metrics based on a joint histogram computed using Parzen Windowing.
typename TransformType::OutputPointType MovingImagePointType
typename ImageSamplerType::OutputVectorContainerPointer ImageSampleContainerPointer
typename MovingImageType::RegionType MovingImageRegionType
typename AdvancedTransformType::NonZeroJacobianIndicesType NonZeroJacobianIndicesType
typename FixedImageIndexType::IndexValueType FixedImageIndexValueType
typename ImageSamplerType::OutputVectorContainerType ImageSampleContainerType
ImageMaskSpatialObject< Self::FixedImageDimension > FixedImageMaskType
SmartPointer< MovingImageMaskType > MovingImageMaskPointer
typename BSplineInterpolatorType::CovariantVectorType MovingImageDerivativeType
typename MovingImageLimiterType::OutputType MovingImageLimiterOutputType
typename TransformType::InputPointType FixedImagePointType
typename FixedImageLimiterType::OutputType FixedImageLimiterOutputType
BSplineInterpolateImageFunction< MovingImageType, CoordinateRepresentationType, double > BSplineInterpolatorType
typename InterpolatorType::ContinuousIndexType MovingImageContinuousIndexType
ImageMaskSpatialObject< Self::MovingImageDimension > MovingImageMaskType
Computes the mutual information between two images to be registered using the method of Mattes et al.
virtual void ComputeJacobianPreconditioner(const TransformJacobianType &jac, const NonZeroJacobianIndicesType &nzji, DerivativeType &preconditioner, DerivativeType &divisor) const
void ComputeDerivativeLowMemorySingleThreaded(DerivativeType &derivative) const
typename AdvancedTransformType::NonZeroJacobianIndicesType NonZeroJacobianIndicesType
void GetValueAndAnalyticDerivative(const ParametersType &parameters, MeasureType &value, DerivativeType &derivative) const override
itkStaticConstMacro(MovingImageDimension, unsigned int, MovingImageType::ImageDimension)
ITK_DISALLOW_COPY_AND_MOVE(ParzenWindowMutualInformationImageToImageMetric)
void ThreadedComputeDerivativeLowMemory(ThreadIdType threadId)
void GetValueAndFiniteDifferenceDerivative(const ParametersType &parameters, MeasureType &value, DerivativeType &derivative) const override
void UpdateDerivativeLowMemory(const RealType fixedImageValue, const RealType movingImageValue, const DerivativeType &imageJacobian, const NonZeroJacobianIndicesType &nzji, DerivativeType &derivative) const
void ComputeDerivativeLowMemory(DerivativeType &derivative) const
MeasureType GetValue(const ParametersType &parameters) const override
itkStaticConstMacro(FixedImageDimension, unsigned int, FixedImageType::ImageDimension)
virtual void GetValueAndAnalyticDerivativeLowMemory(const ParametersType &parameters, MeasureType &value, DerivativeType &derivative) const
static ITK_THREAD_RETURN_FUNCTION_CALL_CONVENTION ComputeDerivativeLowMemoryThreaderCallback(void *arg)
void AfterThreadedComputeDerivativeLowMemory(DerivativeType &derivative) const


Generated on 2024-07-17 for elastix by doxygen 1.11.0 (9b424b03c9833626cd435af22a444888fbbb192d) elastix logo