ITK  6.0.0
Insight Toolkit
itkGradientDifferenceImageToImageMetric.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 itkGradientDifferenceImageToImageMetric_h
19#define itkGradientDifferenceImageToImageMetric_h
20
22
23#include "itkSobelOperator.h"
25#include "itkPoint.h"
26#include "itkCastImageFilter.h"
28
29namespace itk
30{
57template <typename TFixedImage, typename TMovingImage>
58class ITK_TEMPLATE_EXPORT GradientDifferenceImageToImageMetric : public ImageToImageMetric<TFixedImage, TMovingImage>
59{
60public:
61 ITK_DISALLOW_COPY_AND_MOVE(GradientDifferenceImageToImageMetric);
62
66
69
71 itkNewMacro(Self);
72
74 itkOverrideGetNameOfClassMacro(GradientDifferenceImageToImageMetric);
75
77 using typename Superclass::RealType;
78 using typename Superclass::TransformType;
79 using typename Superclass::TransformPointer;
80 using typename Superclass::TransformParametersType;
81 using typename Superclass::TransformJacobianType;
82
83 using typename Superclass::MeasureType;
84 using typename Superclass::DerivativeType;
85 using typename Superclass::FixedImageType;
86 using typename Superclass::MovingImageType;
87 using typename Superclass::FixedImageConstPointer;
88 using typename Superclass::MovingImageConstPointer;
89
90 using FixedImagePixelType = typename TFixedImage::PixelType;
91 using MovedImagePixelType = typename TMovingImage::PixelType;
92
93 static constexpr unsigned int FixedImageDimension = TFixedImage::ImageDimension;
96
98
102
105
107
110 static constexpr unsigned int MovedImageDimension = MovingImageType::ImageDimension;
111
113
116
118
120 void
121 GetDerivative(const TransformParametersType & parameters, DerivativeType & derivative) const override;
122
125 GetValue(const TransformParametersType & parameters) const override;
126
128 void
130 MeasureType & Value,
131 DerivativeType & Derivative) const override;
132
134 void
135 Initialize() override;
136
139 itkSetMacro(DerivativeDelta, double);
140 itkGetConstReferenceMacro(DerivativeDelta, double);
143protected:
146 void
147 PrintSelf(std::ostream & os, Indent indent) const override;
148
150 void
152
154 void
156
159 ComputeMeasure(const TransformParametersType & parameters, const double * subtractionFactor) const;
160
162
164
165private:
167 mutable MovedGradientPixelType m_Variance[FixedImageDimension]{};
168
170 mutable MovedGradientPixelType m_MinMovedGradient[MovedImageDimension]{};
171 mutable MovedGradientPixelType m_MaxMovedGradient[MovedImageDimension]{};
172
174 mutable FixedGradientPixelType m_MinFixedGradient[FixedImageDimension]{};
175 mutable FixedGradientPixelType m_MaxFixedGradient[FixedImageDimension]{};
176
178 typename TransformMovingImageFilterType::Pointer m_TransformMovingImageFilter{};
179
181 CastFixedImageFilterPointer m_CastFixedImageFilter{};
182
183 SobelOperator<FixedGradientPixelType, Self::FixedImageDimension> m_FixedSobelOperators[FixedImageDimension]{};
184
185 typename FixedSobelFilter::Pointer m_FixedSobelFilters[Self::FixedImageDimension]{};
186
189
191 CastMovedImageFilterPointer m_CastMovedImageFilter{};
192
193 SobelOperator<MovedGradientPixelType, Self::MovedImageDimension> m_MovedSobelOperators[MovedImageDimension]{};
194
195 typename MovedSobelFilter::Pointer m_MovedSobelFilters[Self::MovedImageDimension]{};
196
197 double m_DerivativeDelta{};
198};
199} // end namespace itk
200
201#ifndef ITK_MANUAL_INSTANTIATION
202# include "itkGradientDifferenceImageToImageMetric.hxx"
203#endif
204
205#endif
Array class with size defined at construction time.
Definition: itkArray.h:48
Casts input pixels to output pixel type.
Computes similarity between two objects to be registered.
~GradientDifferenceImageToImageMetric() override=default
typename CastFixedImageFilterType::Pointer CastFixedImageFilterPointer
void PrintSelf(std::ostream &os, Indent indent) const override
typename MovedGradientImageType::PixelType MovedGradientPixelType
void GetDerivative(const TransformParametersType &parameters, DerivativeType &derivative) const override
MeasureType GetValue(const TransformParametersType &parameters) const override
MeasureType ComputeMeasure(const TransformParametersType &parameters, const double *subtractionFactor) const
typename CastMovedImageFilterType::Pointer CastMovedImageFilterPointer
typename FixedGradientImageType::PixelType FixedGradientPixelType
void GetValueAndDerivative(const TransformParametersType &parameters, MeasureType &Value, DerivativeType &Derivative) const override
Computes similarity between regions of two images.
typename TransformType::ParametersType TransformParametersType
Templated n-dimensional image class.
Definition: itkImage.h:89
TPixel PixelType
Definition: itkImage.h:108
Control indentation during Print() invocation.
Definition: itkIndent.h:50
Light weight base class for most itk classes.
Applies a single NeighborhoodOperator to an image region.
Resample an image via a coordinate transform.
A function object that determines a neighborhood of values at an image boundary according to a Neuman...
SmartPointer< Self > Pointer
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....