ITK  4.13.0
Insight Segmentation and Registration Toolkit
itkGradientDifferenceImageToImageMetric.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright Insight Software Consortium
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 itkGradientDifferenceImageToImageMetric_h
19 #define itkGradientDifferenceImageToImageMetric_h
20 
21 #include "itkImageToImageMetric.h"
22 
23 #include "itkSobelOperator.h"
25 #include "itkPoint.h"
26 #include "itkCastImageFilter.h"
27 #include "itkResampleImageFilter.h"
28 
29 namespace itk
30 {
57 template< typename TFixedImage, typename TMovingImage >
58 class ITK_TEMPLATE_EXPORT GradientDifferenceImageToImageMetric:
59  public ImageToImageMetric< TFixedImage, TMovingImage >
60 {
61 public:
62 
66 
69 
71  itkNewMacro(Self);
72 
75 
77  typedef typename Superclass::RealType RealType;
78  typedef typename Superclass::TransformType TransformType;
79  typedef typename Superclass::TransformPointer TransformPointer;
80  typedef typename Superclass::TransformParametersType TransformParametersType;
81  typedef typename Superclass::TransformJacobianType TransformJacobianType;
82 
83  typedef typename Superclass::MeasureType MeasureType;
84  typedef typename Superclass::DerivativeType DerivativeType;
85  typedef typename Superclass::FixedImageType FixedImageType;
86  typedef typename Superclass::MovingImageType MovingImageType;
87  typedef typename Superclass::FixedImageConstPointer FixedImageConstPointer;
88  typedef typename Superclass::MovingImageConstPointer MovingImageConstPointer;
89 
90  typedef typename TFixedImage::PixelType FixedImagePixelType;
91  typedef typename TMovingImage::PixelType MovedImagePixelType;
92 
93  itkStaticConstMacro(FixedImageDimension, unsigned int, TFixedImage::ImageDimension);
96  itkGetStaticConstMacro(FixedImageDimension) >
98 
101 
105 
108 
110 
113  itkStaticConstMacro(MovedImageDimension, unsigned int, MovingImageType::ImageDimension);
114 
116 
119 
121 
123  void GetDerivative(const TransformParametersType & parameters,
124  DerivativeType & derivative) const ITK_OVERRIDE;
125 
127  MeasureType GetValue(const TransformParametersType & parameters) const ITK_OVERRIDE;
128 
130  void GetValueAndDerivative(const TransformParametersType & parameters,
131  MeasureType & Value, DerivativeType & derivative) const ITK_OVERRIDE;
132 
135  virtual void Initialize(void) ITK_OVERRIDE;
136 
138  void WriteGradientImagesToFiles() const;
139 
142  itkSetMacro(DerivativeDelta, double);
143  itkGetConstReferenceMacro(DerivativeDelta, double);
145 
146 protected:
148  virtual ~GradientDifferenceImageToImageMetric() ITK_OVERRIDE {}
149  void PrintSelf(std::ostream & os, Indent indent) const ITK_OVERRIDE;
150 
152  void ComputeMovedGradientRange() const;
153 
155  void ComputeVariance() const;
156 
158  MeasureType ComputeMeasure(const TransformParametersType & parameters,
159  const double *subtractionFactor) const;
160 
163 
166 
167 private:
168  ITK_DISALLOW_COPY_AND_ASSIGN(GradientDifferenceImageToImageMetric);
169 
171  mutable MovedGradientPixelType m_Variance[FixedImageDimension];
172 
174  mutable MovedGradientPixelType m_MinMovedGradient[MovedImageDimension];
175  mutable MovedGradientPixelType m_MaxMovedGradient[MovedImageDimension];
176 
178  mutable FixedGradientPixelType m_MinFixedGradient[FixedImageDimension];
179  mutable FixedGradientPixelType m_MaxFixedGradient[FixedImageDimension];
180 
183 
186 
188  itkGetStaticConstMacro(FixedImageDimension) >
189  m_FixedSobelOperators[FixedImageDimension];
190 
191  typename FixedSobelFilter::Pointer m_FixedSobelFilters[itkGetStaticConstMacro(FixedImageDimension)];
192 
195 
198 
200  itkGetStaticConstMacro(MovedImageDimension) >
201  m_MovedSobelOperators[MovedImageDimension];
202 
203  typename MovedSobelFilter::Pointer m_MovedSobelFilters[itkGetStaticConstMacro(MovedImageDimension)];
204 
206 };
207 } // end namespace itk
208 
209 #ifndef ITK_MANUAL_INSTANTIATION
210 #include "itkGradientDifferenceImageToImageMetric.hxx"
211 #endif
212 
213 #endif
Array class with size defined at construction time.
Definition: itkArray.h:50
A function object that determines a neighborhood of values at an image boundary according to a Neuman...
Light weight base class for most itk classes.
NeighborhoodOperatorImageFilter< FixedGradientImageType, FixedGradientImageType > FixedSobelFilter
ImageToImageMetric< TFixedImage, TMovingImage > Superclass
Resample an image via a coordinate transform.
TPixel PixelType
Definition: itkImage.h:89
itk::Image< RealType, itkGetStaticConstMacro(MovedImageDimension) > MovedGradientImageType
A NeighborhoodOperator for performing a directional Sobel edge-detection operation at a pixel locatio...
itk::Image< RealType, itkGetStaticConstMacro(FixedImageDimension) > FixedGradientImageType
NeighborhoodOperatorImageFilter< MovedGradientImageType, MovedGradientImageType > MovedSobelFilter
itk::ResampleImageFilter< MovingImageType, TransformedMovingImageType > TransformMovingImageFilterType
ZeroFluxNeumannBoundaryCondition< FixedGradientImageType > m_FixedBoundCond
TransformMovingImageFilterType::Pointer m_TransformMovingImageFilter
Computes similarity between two objects to be registered.
itk::CastImageFilter< TransformedMovingImageType, MovedGradientImageType > CastMovedImageFilterType
itk::Image< FixedImagePixelType, itkGetStaticConstMacro(FixedImageDimension) > TransformedMovingImageType
Control indentation during Print() invocation.
Definition: itkIndent.h:49
Applies a single NeighborhoodOperator to an image region.
Computes similarity between regions of two images.
itk::CastImageFilter< FixedImageType, FixedGradientImageType > CastFixedImageFilterType
Templated n-dimensional image class.
Definition: itkImage.h:75
ZeroFluxNeumannBoundaryCondition< MovedGradientImageType > m_MovedBoundCond
Casts input pixels to output pixel type.