go home Home | Main Page | Modules | Namespace List | Class Hierarchy | Alphabetical List | Data Structures | File List | Namespace Members | Data Fields | Globals | Related Pages
itkParzenWindowHistogramImageToImageMetric.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 itkParzenWindowHistogramImageToImageMetric_h
19#define itkParzenWindowHistogramImageToImageMetric_h
20
23#include <vector>
24
25
26namespace itk
27{
74template <class TFixedImage, class TMovingImage>
76 : public AdvancedImageToImageMetric<TFixedImage, TMovingImage>
77{
78public:
80
84 using Pointer = SmartPointer<Self>;
85 using ConstPointer = SmartPointer<const Self>;
86
89
91 using typename Superclass::CoordinateRepresentationType;
92 using typename Superclass::MovingImageType;
93 using typename Superclass::MovingImagePixelType;
94 using typename Superclass::MovingImageConstPointer;
95 using typename Superclass::FixedImageType;
96 using typename Superclass::FixedImageConstPointer;
97 using typename Superclass::FixedImageRegionType;
98 using typename Superclass::TransformType;
99 using typename Superclass::TransformPointer;
100 using typename Superclass::InputPointType;
101 using typename Superclass::OutputPointType;
102 using typename Superclass::TransformParametersType;
103 using typename Superclass::TransformJacobianType;
104 using typename Superclass::InterpolatorType;
105 using typename Superclass::InterpolatorPointer;
106 using typename Superclass::RealType;
107 using typename Superclass::GradientPixelType;
108 using typename Superclass::GradientImageType;
109 using typename Superclass::GradientImagePointer;
110 using typename Superclass::FixedImageMaskType;
112 using typename Superclass::MovingImageMaskType;
114 using typename Superclass::MeasureType;
115 using typename Superclass::DerivativeType;
116 using typename Superclass::DerivativeValueType;
117 using typename Superclass::ParametersType;
118 using typename Superclass::FixedImagePixelType;
120 using typename Superclass::ImageSamplerType;
121 using typename Superclass::ImageSamplerPointer;
129 using typename Superclass::ThreadInfoType;
130
132 itkStaticConstMacro(FixedImageDimension, unsigned int, FixedImageType::ImageDimension);
133
135 itkStaticConstMacro(MovingImageDimension, unsigned int, MovingImageType::ImageDimension);
136
143 void
144 Initialize() override;
145
150 void
151 GetDerivative(const ParametersType & parameters, DerivativeType & Derivative) const override;
152
158 void
159 GetValueAndDerivative(const ParametersType & parameters,
160 MeasureType & value,
161 DerivativeType & derivative) const override;
162
169 itkSetClampMacro(NumberOfFixedHistogramBins, unsigned long, 4, NumericTraits<unsigned long>::max());
170 itkGetConstMacro(NumberOfFixedHistogramBins, unsigned long);
178 itkSetClampMacro(NumberOfMovingHistogramBins, unsigned long, 4, NumericTraits<unsigned long>::max());
179 itkGetConstMacro(NumberOfMovingHistogramBins, unsigned long);
182 itkSetClampMacro(FixedKernelBSplineOrder, unsigned int, 0, 3);
183 itkGetConstMacro(FixedKernelBSplineOrder, unsigned int);
184
186 itkSetClampMacro(MovingKernelBSplineOrder, unsigned int, 0, 3);
187 itkGetConstMacro(MovingKernelBSplineOrder, unsigned int);
188
192 itkSetMacro(UseExplicitPDFDerivatives, bool);
193 itkGetConstReferenceMacro(UseExplicitPDFDerivatives, bool);
194 itkBooleanMacro(UseExplicitPDFDerivatives);
195
199 itkSetMacro(UseDerivative, bool);
200 itkGetConstMacro(UseDerivative, bool);
201
205 itkSetMacro(UseFiniteDifferenceDerivative, bool);
206 itkGetConstMacro(UseFiniteDifferenceDerivative, bool);
207
212 itkSetMacro(FiniteDifferencePerturbation, double);
213 itkGetConstMacro(FiniteDifferencePerturbation, double);
214
215protected:
218
221
223 void
224 PrintSelf(std::ostream & os, Indent indent) const override;
225
229 using typename Superclass::FixedImageIndexType;
231 using OffsetValueType = typename FixedImageType::OffsetValueType;
233 using typename Superclass::FixedImagePointType;
239
243 using MarginalPDFType = Array<PDFValueType>;
245 using JointPDFPointer = typename JointPDFType::Pointer;
247 using JointPDFDerivativesPointer = typename JointPDFDerivativesType::Pointer;
249 using IncrementalMarginalPDFPointer = typename IncrementalMarginalPDFType::Pointer;
250 using JointPDFIndexType = JointPDFType::IndexType;
251 using JointPDFRegionType = JointPDFType::RegionType;
252 using JointPDFSizeType = JointPDFType::SizeType;
253 using JointPDFDerivativesIndexType = JointPDFDerivativesType::IndexType;
254 using JointPDFDerivativesRegionType = JointPDFDerivativesType::RegionType;
255 using JointPDFDerivativesSizeType = JointPDFDerivativesType::SizeType;
256 using IncrementalMarginalPDFIndexType = IncrementalMarginalPDFType::IndexType;
257 using IncrementalMarginalPDFRegionType = IncrementalMarginalPDFType::RegionType;
258 using IncrementalMarginalPDFSizeType = IncrementalMarginalPDFType::SizeType;
259 using ParzenValueContainerType = Array<PDFValueType>;
260
264
268 mutable double m_Alpha{ 0.0 };
269 mutable DerivativeType m_PerturbedAlphaRight{};
270 mutable DerivativeType m_PerturbedAlphaLeft{};
271
273 mutable MarginalPDFType m_FixedImageMarginalPDF{};
274 mutable MarginalPDFType m_MovingImageMarginalPDF{};
275 JointPDFPointer m_JointPDF{ nullptr };
276 JointPDFDerivativesPointer m_JointPDFDerivatives{ nullptr };
277 JointPDFDerivativesPointer m_IncrementalJointPDFRight{};
278 JointPDFDerivativesPointer m_IncrementalJointPDFLeft{};
279 IncrementalMarginalPDFPointer m_FixedIncrementalMarginalPDFRight{ nullptr };
280 IncrementalMarginalPDFPointer m_MovingIncrementalMarginalPDFRight{ nullptr };
281 IncrementalMarginalPDFPointer m_FixedIncrementalMarginalPDFLeft{ nullptr };
282 IncrementalMarginalPDFPointer m_MovingIncrementalMarginalPDFLeft{ nullptr };
283 mutable JointPDFRegionType m_JointPDFWindow{}; // no need for mutable anymore?
284 double m_MovingImageNormalizedMin{ 0.0 };
285 double m_FixedImageNormalizedMin{ 0.0 };
286 double m_FixedImageBinSize{ 0.0 };
287 double m_MovingImageBinSize{ 0.0 };
288 double m_FixedParzenTermToIndexOffset{ 0.5 };
289 double m_MovingParzenTermToIndexOffset{ -1.0 };
290
292 KernelFunctionPointer m_FixedKernel{ nullptr };
293 KernelFunctionPointer m_MovingKernel{ nullptr };
294 KernelFunctionPointer m_DerivativeMovingKernel{ nullptr };
295
297 void
299
301 void
302 ThreadedComputePDFs(ThreadIdType threadId);
303
305 void
307
309 static ITK_THREAD_RETURN_FUNCTION_CALL_CONVENTION
311
313 void
315
321 static void
322 EvaluateParzenValues(double parzenWindowTerm,
323 OffsetValueType parzenWindowIndex,
324 const KernelFunctionType & kernel,
325 PDFValueType * parzenValues);
326
330 virtual void
331 UpdateJointPDFAndDerivatives(const RealType fixedImageValue,
332 const RealType movingImageValue,
333 const DerivativeType * imageJacobian,
334 const NonZeroJacobianIndicesType * nzji,
335 JointPDFType * jointPDF) const;
336
347 virtual void
348 UpdateJointPDFAndIncrementalPDFs(RealType fixedImageValue,
349 RealType movingImageValue,
350 RealType movingMaskValue,
351 const DerivativeType & movingImageValuesRight,
352 const DerivativeType & movingImageValuesLeft,
353 const DerivativeType & movingMaskValuesRight,
354 const DerivativeType & movingMaskValuesLeft,
355 const NonZeroJacobianIndicesType & nzji) const;
356
362 void
364 double factor,
365 const DerivativeType & imageJacobian,
366 const NonZeroJacobianIndicesType & nzji) const;
367
369 void
370 NormalizeJointPDF(JointPDFType * pdf, const double factor) const;
371
373 void
374 NormalizeJointPDFDerivatives(JointPDFDerivativesType * pdf, const double factor) const;
375
380 void
381 ComputeMarginalPDF(const JointPDFType * jointPDF, MarginalPDFType & marginalPDF, const unsigned int direction) const;
382
386 virtual void
388 IncrementalMarginalPDFType * fixedIncrementalMarginalPDF,
389 IncrementalMarginalPDFType * movingIncrementalMarginalPDF) const;
390
400 virtual void
401 ComputePDFsAndPDFDerivatives(const ParametersType & parameters) const;
402
426 virtual void
427 ComputePDFsAndIncrementalPDFs(const ParametersType & parameters) const;
428
437 virtual void
438 ComputePDFsSingleThreaded(const ParametersType & parameters) const;
439
440 virtual void
441 ComputePDFs(const ParametersType & parameters) const;
442
444 virtual void
446
447 virtual void
449
454 virtual void
455 GetValueAndAnalyticDerivative(const ParametersType & itkNotUsed(parameters),
456 MeasureType & itkNotUsed(value),
457 DerivativeType & itkNotUsed(derivative)) const
458 {}
459
464 virtual void
465 GetValueAndFiniteDifferenceDerivative(const ParametersType & itkNotUsed(parameters),
466 MeasureType & itkNotUsed(value),
467 DerivativeType & itkNotUsed(derivative)) const
468 {}
469
470private:
472 mutable std::vector<JointPDFPointer> m_ThreaderJointPDFs{};
473
477 struct ParzenWindowHistogramMultiThreaderParameterType // can't we use the one from AdvancedImageToImageMetric ?
478 {
480 };
481 ParzenWindowHistogramMultiThreaderParameterType m_ParzenWindowHistogramThreaderParameters{};
482
488 itkPadStruct(ITK_CACHE_LINE_ALIGNMENT,
490 PaddedParzenWindowHistogramGetValueAndDerivativePerThreadStruct);
491 itkAlignedTypedef(ITK_CACHE_LINE_ALIGNMENT,
492 PaddedParzenWindowHistogramGetValueAndDerivativePerThreadStruct,
493 AlignedParzenWindowHistogramGetValueAndDerivativePerThreadStruct);
494 mutable std::vector<AlignedParzenWindowHistogramGetValueAndDerivativePerThreadStruct>
496
498 unsigned long m_NumberOfFixedHistogramBins{ 32 };
499 unsigned long m_NumberOfMovingHistogramBins{ 32 };
500 unsigned int m_FixedKernelBSplineOrder{ 0 };
501 unsigned int m_MovingKernelBSplineOrder{ 3 };
502 bool m_UseDerivative{ false };
503 bool m_UseExplicitPDFDerivatives{ true };
504 bool m_UseFiniteDifferenceDerivative{ false };
505 double m_FiniteDifferencePerturbation{ 1.0 };
506};
507
508} // end namespace itk
509
510#ifndef ITK_MANUAL_INSTANTIATION
511# include "itkParzenWindowHistogramImageToImageMetric.hxx"
512#endif
513
514#endif // end #ifndef itkParzenWindowHistogramImageToImageMetric_h
An extension of the ITK ImageToImageMetric. It is the intended base class for all elastix metrics.
typename TransformType::OutputPointType MovingImagePointType
typename ImageSamplerType::OutputVectorContainerPointer ImageSampleContainerPointer
typename MovingImageType::RegionType MovingImageRegionType
typename AdvancedTransformType::NonZeroJacobianIndicesType NonZeroJacobianIndicesType
typename FixedImageType::PixelType FixedImagePixelType
typename DerivativeType::ValueType DerivativeValueType
typename FixedImageIndexType::IndexValueType FixedImageIndexValueType
typename MovingImageType::IndexType MovingImageIndexType
typename ImageSamplerType::OutputVectorContainerType ImageSampleContainerType
ImageMaskSpatialObject< Self::FixedImageDimension > FixedImageMaskType
typename FixedImageType::IndexType FixedImageIndexType
SmartPointer< MovingImageMaskType > MovingImageMaskPointer
MultiThreaderBase::WorkUnitInfo ThreadInfoType
typename BSplineInterpolatorType::CovariantVectorType MovingImageDerivativeType
typename MovingImageLimiterType::OutputType MovingImageLimiterOutputType
typename TransformType::InputPointType FixedImagePointType
typename FixedImageLimiterType::OutputType FixedImageLimiterOutputType
SmartPointer< FixedImageMaskType > FixedImageMaskPointer
BSplineInterpolateImageFunction< MovingImageType, CoordinateRepresentationType, double > BSplineInterpolatorType
typename ImageSamplerType::Pointer ImageSamplerPointer
typename InterpolatorType::ContinuousIndexType MovingImageContinuousIndexType
ImageMaskSpatialObject< Self::MovingImageDimension > MovingImageMaskType
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.
void InitializeThreadingParameters() const override
ITK_DISALLOW_COPY_AND_MOVE(ParzenWindowHistogramImageToImageMetric)
~ParzenWindowHistogramImageToImageMetric() override=default
virtual void ComputePDFsAndIncrementalPDFs(const ParametersType &parameters) const
typename AdvancedTransformType::NonZeroJacobianIndicesType NonZeroJacobianIndicesType
void NormalizeJointPDFDerivatives(JointPDFDerivativesType *pdf, const double factor) const
std::vector< AlignedParzenWindowHistogramGetValueAndDerivativePerThreadStruct > m_ParzenWindowHistogramGetValueAndDerivativePerThreadVariables
void ThreadedComputePDFs(ThreadIdType threadId)
virtual void ComputePDFsSingleThreaded(const ParametersType &parameters) const
itkPadStruct(ITK_CACHE_LINE_ALIGNMENT, ParzenWindowHistogramGetValueAndDerivativePerThreadStruct, PaddedParzenWindowHistogramGetValueAndDerivativePerThreadStruct)
static ITK_THREAD_RETURN_FUNCTION_CALL_CONVENTION ComputePDFsThreaderCallback(void *arg)
virtual void ComputePDFsAndPDFDerivatives(const ParametersType &parameters) const
void ComputeMarginalPDF(const JointPDFType *jointPDF, MarginalPDFType &marginalPDF, const unsigned int direction) const
void GetDerivative(const ParametersType &parameters, DerivativeType &Derivative) const override
itkAlignedTypedef(ITK_CACHE_LINE_ALIGNMENT, PaddedParzenWindowHistogramGetValueAndDerivativePerThreadStruct, AlignedParzenWindowHistogramGetValueAndDerivativePerThreadStruct)
void UpdateJointPDFDerivatives(const JointPDFIndexType &pdfIndex, double factor, const DerivativeType &imageJacobian, const NonZeroJacobianIndicesType &nzji) const
virtual void GetValueAndAnalyticDerivative(const ParametersType &, MeasureType &, DerivativeType &) const
virtual void ComputePDFs(const ParametersType &parameters) const
typename IncrementalMarginalPDFType::Pointer IncrementalMarginalPDFPointer
virtual void UpdateJointPDFAndIncrementalPDFs(RealType fixedImageValue, RealType movingImageValue, RealType movingMaskValue, const DerivativeType &movingImageValuesRight, const DerivativeType &movingImageValuesLeft, const DerivativeType &movingMaskValuesRight, const DerivativeType &movingMaskValuesLeft, const NonZeroJacobianIndicesType &nzji) const
static void EvaluateParzenValues(double parzenWindowTerm, OffsetValueType parzenWindowIndex, const KernelFunctionType &kernel, PDFValueType *parzenValues)
void NormalizeJointPDF(JointPDFType *pdf, const double factor) const
virtual void UpdateJointPDFAndDerivatives(const RealType fixedImageValue, const RealType movingImageValue, const DerivativeType *imageJacobian, const NonZeroJacobianIndicesType *nzji, JointPDFType *jointPDF) const
itkStaticConstMacro(FixedImageDimension, unsigned int, FixedImageType::ImageDimension)
void GetValueAndDerivative(const ParametersType &parameters, MeasureType &value, DerivativeType &derivative) const override
void PrintSelf(std::ostream &os, Indent indent) const override
virtual void ComputeIncrementalMarginalPDFs(const JointPDFDerivativesType *incrementalPDF, IncrementalMarginalPDFType *fixedIncrementalMarginalPDF, IncrementalMarginalPDFType *movingIncrementalMarginalPDF) const
virtual void GetValueAndFiniteDifferenceDerivative(const ParametersType &, MeasureType &, DerivativeType &) const
itkStaticConstMacro(MovingImageDimension, unsigned int, MovingImageType::ImageDimension)


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