ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
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{
54template <typename TFixedImage, typename TMovingImage>
55class ITK_TEMPLATE_EXPORT GradientDifferenceImageToImageMetric : public ImageToImageMetric<TFixedImage, TMovingImage>
56{
57public:
58 ITK_DISALLOW_COPY_AND_MOVE(GradientDifferenceImageToImageMetric);
59
63
66
68 itkNewMacro(Self);
69
71 itkOverrideGetNameOfClassMacro(GradientDifferenceImageToImageMetric);
72
74 using typename Superclass::RealType;
75 using typename Superclass::TransformType;
76 using typename Superclass::TransformPointer;
79
80 using typename Superclass::MeasureType;
81 using typename Superclass::DerivativeType;
82 using typename Superclass::FixedImageType;
83 using typename Superclass::MovingImageType;
86
87 using FixedImagePixelType = typename TFixedImage::PixelType;
88 using MovedImagePixelType = typename TMovingImage::PixelType;
89
90 static constexpr unsigned int FixedImageDimension = TFixedImage::ImageDimension;
93
95
97
99
102
104
106
107 static constexpr unsigned int MovedImageDimension = MovingImageType::ImageDimension;
108
110
113
115
117 void
118 GetDerivative(const TransformParametersType & parameters, DerivativeType & derivative) const override;
119
122 GetValue(const TransformParametersType & parameters) const override;
123
125 void
127 MeasureType & Value,
128 DerivativeType & Derivative) const override;
129
131 void
132 Initialize() override;
133
136 itkSetMacro(DerivativeDelta, double);
137 itkGetConstReferenceMacro(DerivativeDelta, double);
139
140protected:
143 void
144 PrintSelf(std::ostream & os, Indent indent) const override;
145
147 void
149
151 void
153
156 ComputeMeasure(const TransformParametersType & parameters, const double * subtractionFactor) const;
157
159
161
162private:
165
169
173
176
179
181
183
186
189
191
193
195};
196} // end namespace itk
197
198#ifndef ITK_MANUAL_INSTANTIATION
199# include "itkGradientDifferenceImageToImageMetric.hxx"
200#endif
201
202#endif
Casts input pixels to output pixel type.
SobelOperator< FixedGradientPixelType, Self::FixedImageDimension > m_FixedSobelOperators[FixedImageDimension]
ImageToImageMetric< TFixedImage, TMovingImage > Superclass
TransformMovingImageFilterType::Pointer m_TransformMovingImageFilter
~GradientDifferenceImageToImageMetric() override=default
typename CastFixedImageFilterType::Pointer CastFixedImageFilterPointer
FixedGradientPixelType m_MaxFixedGradient[FixedImageDimension]
typename TransformType::ParametersType TransformParametersType
itk::CastImageFilter< TransformedMovingImageType, MovedGradientImageType > CastMovedImageFilterType
NeighborhoodOperatorImageFilter< FixedGradientImageType, FixedGradientImageType > FixedSobelFilter
MovedGradientPixelType m_MinMovedGradient[MovedImageDimension]
void PrintSelf(std::ostream &os, Indent indent) const override
MovedGradientPixelType m_MaxMovedGradient[MovedImageDimension]
ZeroFluxNeumannBoundaryCondition< FixedGradientImageType > m_FixedBoundCond
MovedSobelFilter::Pointer m_MovedSobelFilters[Self::MovedImageDimension]
typename MovedGradientImageType::PixelType MovedGradientPixelType
FixedGradientPixelType m_MinFixedGradient[FixedImageDimension]
void GetDerivative(const TransformParametersType &parameters, DerivativeType &derivative) const override
itk::Image< RealType, Self::MovedImageDimension > MovedGradientImageType
MeasureType GetValue(const TransformParametersType &parameters) const override
NeighborhoodOperatorImageFilter< MovedGradientImageType, MovedGradientImageType > MovedSobelFilter
MeasureType ComputeMeasure(const TransformParametersType &parameters, const double *subtractionFactor) const
ZeroFluxNeumannBoundaryCondition< MovedGradientImageType > m_MovedBoundCond
FixedSobelFilter::Pointer m_FixedSobelFilters[Self::FixedImageDimension]
typename CastMovedImageFilterType::Pointer CastMovedImageFilterPointer
SobelOperator< MovedGradientPixelType, Self::MovedImageDimension > m_MovedSobelOperators[MovedImageDimension]
typename FixedGradientImageType::PixelType FixedGradientPixelType
itk::CastImageFilter< FixedImageType, FixedGradientImageType > CastFixedImageFilterType
itk::Image< RealType, Self::FixedImageDimension > FixedGradientImageType
itk::ResampleImageFilter< MovingImageType, TransformedMovingImageType > TransformMovingImageFilterType
void GetValueAndDerivative(const TransformParametersType &parameters, MeasureType &Value, DerivativeType &Derivative) const override
itk::Image< FixedImagePixelType, Self::FixedImageDimension > TransformedMovingImageType
typename FixedImageType::ConstPointer FixedImageConstPointer
typename NumericTraits< MovingImagePixelType >::RealType RealType
typename TransformType::Pointer TransformPointer
Array< ParametersValueType > DerivativeType
typename MovingImageType::ConstPointer MovingImageConstPointer
typename TransformType::ParametersType TransformParametersType
Transform< CoordinateRepresentationType, Self::MovingImageDimension, Self::FixedImageDimension > TransformType
typename TransformType::JacobianType TransformJacobianType
Templated n-dimensional image class.
Definition itkImage.h:89
Control indentation during Print() invocation.
Definition itkIndent.h:50
Applies a single NeighborhoodOperator to an image region.
Resample an image via a coordinate transform.
Implements transparent reference counting.
A NeighborhoodOperator for performing a directional Sobel edge-detection operation at a pixel locatio...
A function object that determines a neighborhood of values at an image boundary according to a Neuman...
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....