go home Home | Main Page | Modules | Namespace List | Class Hierarchy | Alphabetical List | Data Structures | File List | Namespace Members | Data Fields | Globals | Related Pages
itkGradientDifferenceImageToImageMetric2.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 /*=========================================================================
19 
20  Program: Insight Segmentation & Registration Toolkit
21  Module: $RCSfile: itkGradientDifferenceImageToImageMetric2.h,v $
22  Language: C++
23  Date: $Date: 2011-29-04 14:33 $
24  Version: $Revision: 2.0 $
25 
26  Copyright (c) Insight Software Consortium. All rights reserved.
27  See ITKCopyright.txt or http://www.itk.org/HTML/Copyright.htm for details.
28 
29  This software is distributed WITHOUT ANY WARRANTY; without even
30  the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
31  PURPOSE. See the above copyright notices for more information.
32 
33 =========================================================================*/
34 #ifndef __itkGradientDifferenceImageToImageMetric2_h
35 #define __itkGradientDifferenceImageToImageMetric2_h
36 
38 
39 #include "itkSobelOperator.h"
40 #include "itkNeighborhoodOperatorImageFilter.h"
41 #include "itkPoint.h"
42 #include "itkCastImageFilter.h"
43 #include "itkResampleImageFilter.h"
44 #include "itkOptimizer.h"
47 
48 namespace itk
49 {
74 template< class TFixedImage, class TMovingImage >
76  public AdvancedImageToImageMetric< TFixedImage, TMovingImage >
77 {
78 public:
79 
83 
85  typedef SmartPointer< const Self > ConstPointer;
86 
88  itkNewMacro( Self );
89 
91  itkTypeMacro( GradientDifferenceImageToImageMetric, ImageToImageMetric );
92 
95  #if defined( _MSC_VER ) && ( _MSC_VER == 1300 )
96  typedef double RealType;
97  #else
98  typedef typename Superclass::RealType RealType;
99  #endif
100 
102  typedef typename TransformType::ScalarType ScalarType;
107  typedef typename InterpolatorType::Pointer InterpolatorPointer;
114  typedef typename TFixedImage::PixelType FixedImagePixelType;
115  typedef typename TMovingImage::PixelType MovedImagePixelType;
116  typedef typename MovingImageType::RegionType MovingImageRegionType;
117  typedef typename itk::Optimizer OptimizerType;
118  typedef typename OptimizerType::ScalesType ScalesType;
119 
120  itkStaticConstMacro( FixedImageDimension, unsigned int,
121  FixedImageType::ImageDimension );
122  itkStaticConstMacro( MovedImageDimension, unsigned int,
123  MovingImageType::ImageDimension );
124 
128  typedef itk::Image< FixedImagePixelType, itkGetStaticConstMacro( FixedImageDimension ) >
130  typedef itk::ResampleImageFilter< MovingImageType, TransformedMovingImageType >
135  typedef itk::Image< RealType, itkGetStaticConstMacro( FixedImageDimension ) >
137  typedef itk::CastImageFilter< FixedImageType, FixedGradientImageType >
139  typedef typename CastFixedImageFilterType::Pointer CastFixedImageFilterPointer;
140  typedef typename FixedGradientImageType::PixelType FixedGradientPixelType;
141  typedef itk::Image< RealType, itkGetStaticConstMacro( MovedImageDimension ) >
143  typedef itk::CastImageFilter< TransformedMovingImageType, MovedGradientImageType >
145  typedef typename CastMovedImageFilterType::Pointer CastMovedImageFilterPointer;
146  typedef typename MovedGradientImageType::PixelType MovedGradientPixelType;
147 
149  void GetDerivative( const TransformParametersType & parameters,
150  DerivativeType & derivative ) const;
151 
153  MeasureType GetValue( const TransformParametersType & parameters ) const;
154 
157  MeasureType & Value, DerivativeType & derivative ) const;
158 
159  virtual void Initialize( void ) throw ( ExceptionObject );
160 
162  void WriteGradientImagesToFiles( void ) const;
163 
165  itkSetMacro( Scales, ScalesType );
166  itkGetConstReferenceMacro( Scales, ScalesType );
167 
170  itkSetMacro( DerivativeDelta, double );
171  itkGetConstReferenceMacro( DerivativeDelta, double );
172 
173 protected:
174 
177  void PrintSelf( std::ostream & os, Indent indent ) const;
178 
180  void ComputeMovedGradientRange( void ) const;
181 
183  void ComputeVariance( void ) const;
184 
187  const double * subtractionFactor ) const;
188 
189  typedef NeighborhoodOperatorImageFilter<
191 
192  typedef NeighborhoodOperatorImageFilter<
194 
195 private:
196 
197  GradientDifferenceImageToImageMetric( const Self & ); // purposely not implemented
198  void operator=( const Self & ); // purposely not implemented
199 
201  mutable MovedGradientPixelType m_Variance[ FixedImageDimension ];
202 
204  mutable MovedGradientPixelType m_MinMovedGradient[ MovedImageDimension ];
205  mutable MovedGradientPixelType m_MaxMovedGradient[ MovedImageDimension ];
206 
208  mutable FixedGradientPixelType m_MinFixedGradient[ FixedImageDimension ];
209  mutable FixedGradientPixelType m_MaxFixedGradient[ FixedImageDimension ];
210 
212  typename TransformMovingImageFilterType::Pointer m_TransformMovingImageFilter;
213 
216 
217  SobelOperator< FixedGradientPixelType,
218  itkGetStaticConstMacro( FixedImageDimension ) >
219  m_FixedSobelOperators[ FixedImageDimension ];
220 
221  typename FixedSobelFilter::Pointer m_FixedSobelFilters[ itkGetStaticConstMacro
222  ( FixedImageDimension ) ];
223 
226 
229 
230  SobelOperator< MovedGradientPixelType,
231  itkGetStaticConstMacro( MovedImageDimension ) >
232  m_MovedSobelOperators[ MovedImageDimension ];
233 
234  typename MovedSobelFilter::Pointer m_MovedSobelFilters[ itkGetStaticConstMacro
235  ( MovedImageDimension ) ];
236 
241 
242 };
243 
244 } // end namespace itk
245 
246 #ifndef ITK_MANUAL_INSTANTIATION
247 #include "itkGradientDifferenceImageToImageMetric2.hxx"
248 #endif
249 
250 #endif
itk::GradientDifferenceImageToImageMetric::MovingImageRegionType
MovingImageType::RegionType MovingImageRegionType
Definition: itkGradientDifferenceImageToImageMetric2.h:116
itk::AdvancedImageToImageMetric::TransformParametersType
Superclass::TransformParametersType TransformParametersType
Definition: itkAdvancedImageToImageMetric.h:113
itk::GradientDifferenceImageToImageMetric::FixedSobelFilter
NeighborhoodOperatorImageFilter< FixedGradientImageType, FixedGradientImageType > FixedSobelFilter
Definition: itkGradientDifferenceImageToImageMetric2.h:190
itk::GradientDifferenceImageToImageMetric::m_Scales
ScalesType m_Scales
Definition: itkGradientDifferenceImageToImageMetric2.h:237
itk::GradientDifferenceImageToImageMetric::ScalarType
TransformType::ScalarType ScalarType
Definition: itkGradientDifferenceImageToImageMetric2.h:102
itk::GradientDifferenceImageToImageMetric::FixedImageType
Superclass::FixedImageType FixedImageType
Definition: itkGradientDifferenceImageToImageMetric2.h:110
itk::GradientDifferenceImageToImageMetric::GetDerivative
void GetDerivative(const TransformParametersType &parameters, DerivativeType &derivative) const
itk::GradientDifferenceImageToImageMetric::FixedImagePixelType
TFixedImage::PixelType FixedImagePixelType
Definition: itkGradientDifferenceImageToImageMetric2.h:114
itk::GradientDifferenceImageToImageMetric::GetValueAndDerivative
void GetValueAndDerivative(const TransformParametersType &parameters, MeasureType &Value, DerivativeType &derivative) const
itk::GradientDifferenceImageToImageMetric::TransformType
Superclass::TransformType TransformType
Definition: itkGradientDifferenceImageToImageMetric2.h:101
itk::GradientDifferenceImageToImageMetric::m_CombinationTransform
CombinationTransformPointer m_CombinationTransform
Definition: itkGradientDifferenceImageToImageMetric2.h:240
SmartPointer< Self >
itk::AdvancedImageToImageMetric
An extension of the ITK ImageToImageMetric. It is the intended base class for all elastix metrics.
Definition: itkAdvancedImageToImageMetric.h:81
itk::GradientDifferenceImageToImageMetric::CastFixedImageFilterType
itk::CastImageFilter< FixedImageType, FixedGradientImageType > CastFixedImageFilterType
Definition: itkGradientDifferenceImageToImageMetric2.h:138
itk::GradientDifferenceImageToImageMetric::DerivativeType
Superclass::DerivativeType DerivativeType
Definition: itkGradientDifferenceImageToImageMetric2.h:109
itk::GradientDifferenceImageToImageMetric::m_MinFixedGradient
FixedGradientPixelType m_MinFixedGradient[FixedImageDimension]
Definition: itkGradientDifferenceImageToImageMetric2.h:208
itk::AdvancedImageToImageMetric::DerivativeType
Superclass::DerivativeType DerivativeType
Definition: itkAdvancedImageToImageMetric.h:128
itk::GradientDifferenceImageToImageMetric::ScalesType
OptimizerType::ScalesType ScalesType
Definition: itkGradientDifferenceImageToImageMetric2.h:118
itkAdvancedCombinationTransform.h
itk::GradientDifferenceImageToImageMetric::WriteGradientImagesToFiles
void WriteGradientImagesToFiles(void) const
itk::GradientDifferenceImageToImageMetric::CastMovedImageFilterPointer
CastMovedImageFilterType::Pointer CastMovedImageFilterPointer
Definition: itkGradientDifferenceImageToImageMetric2.h:145
itk::GradientDifferenceImageToImageMetric::GradientDifferenceImageToImageMetric
GradientDifferenceImageToImageMetric()
itk::GradientDifferenceImageToImageMetric::FixedGradientImageType
itk::Image< RealType, itkGetStaticConstMacro(FixedImageDimension) > FixedGradientImageType
Definition: itkGradientDifferenceImageToImageMetric2.h:136
itk::GradientDifferenceImageToImageMetric
Computes similarity between two objects to be registered.
Definition: itkGradientDifferenceImageToImageMetric2.h:77
itk::GradientDifferenceImageToImageMetric::MovedImagePixelType
TMovingImage::PixelType MovedImagePixelType
Definition: itkGradientDifferenceImageToImageMetric2.h:115
itk::GradientDifferenceImageToImageMetric::m_CastFixedImageFilter
CastFixedImageFilterPointer m_CastFixedImageFilter
Definition: itkGradientDifferenceImageToImageMetric2.h:215
itk::GradientDifferenceImageToImageMetric::InterpolatorPointer
InterpolatorType::Pointer InterpolatorPointer
Definition: itkGradientDifferenceImageToImageMetric2.h:107
itk::GradientDifferenceImageToImageMetric::ComputeVariance
void ComputeVariance(void) const
itk::AdvancedCombinationTransform
This class combines two transforms: an 'initial transform' with a 'current transform'.
Definition: itkAdvancedCombinationTransform.h:58
itk::GradientDifferenceImageToImageMetric::m_Variance
MovedGradientPixelType m_Variance[FixedImageDimension]
Definition: itkGradientDifferenceImageToImageMetric2.h:201
itk::GradientDifferenceImageToImageMetric::~GradientDifferenceImageToImageMetric
virtual ~GradientDifferenceImageToImageMetric()
Definition: itkGradientDifferenceImageToImageMetric2.h:176
itk::AdvancedImageToImageMetric::InterpolatorType
Superclass::InterpolatorType InterpolatorType
Definition: itkAdvancedImageToImageMetric.h:115
itk::AdvancedImageToImageMetric::TransformJacobianType
Superclass::TransformJacobianType TransformJacobianType
Definition: itkAdvancedImageToImageMetric.h:114
itk::GradientDifferenceImageToImageMetric::OptimizerType
itk::Optimizer OptimizerType
Definition: itkGradientDifferenceImageToImageMetric2.h:117
itk::GradientDifferenceImageToImageMetric::MovingImageConstPointer
Superclass::MovingImageConstPointer MovingImageConstPointer
Definition: itkGradientDifferenceImageToImageMetric2.h:113
itk::GradientDifferenceImageToImageMetric::FixedGradientPixelType
FixedGradientImageType::PixelType FixedGradientPixelType
Definition: itkGradientDifferenceImageToImageMetric2.h:140
itk::GradientDifferenceImageToImageMetric::RealType
Superclass::RealType RealType
Definition: itkGradientDifferenceImageToImageMetric2.h:91
itk::GradientDifferenceImageToImageMetric::m_FixedSobelFilters
FixedSobelFilter::Pointer m_FixedSobelFilters[itkGetStaticConstMacro(FixedImageDimension)]
Definition: itkGradientDifferenceImageToImageMetric2.h:222
itk::GradientDifferenceImageToImageMetric::TransformedMovingImageType
itk::Image< FixedImagePixelType, itkGetStaticConstMacro(FixedImageDimension) > TransformedMovingImageType
Definition: itkGradientDifferenceImageToImageMetric2.h:129
itk::GradientDifferenceImageToImageMetric::CastFixedImageFilterPointer
CastFixedImageFilterType::Pointer CastFixedImageFilterPointer
Definition: itkGradientDifferenceImageToImageMetric2.h:139
itk::AdvancedRayCastInterpolateImageFunction
Projective interpolation of an image at specified positions.
Definition: itkAdvancedRayCastInterpolateImageFunction.h:58
itk::GradientDifferenceImageToImageMetric::RayCastInterpolatorPointer
RayCastInterpolatorType::Pointer RayCastInterpolatorPointer
Definition: itkGradientDifferenceImageToImageMetric2.h:134
itk::GradientDifferenceImageToImageMetric::Superclass
AdvancedImageToImageMetric< TFixedImage, TMovingImage > Superclass
Definition: itkGradientDifferenceImageToImageMetric2.h:82
itk::GradientDifferenceImageToImageMetric::Initialize
virtual void Initialize(void)
itk::GradientDifferenceImageToImageMetric::FixedImageConstPointer
Superclass::FixedImageConstPointer FixedImageConstPointer
Definition: itkGradientDifferenceImageToImageMetric2.h:112
itk::GradientDifferenceImageToImageMetric::TransformParametersType
Superclass::TransformParametersType TransformParametersType
Definition: itkGradientDifferenceImageToImageMetric2.h:104
itk::GradientDifferenceImageToImageMetric::TransformMovingImageFilterType
itk::ResampleImageFilter< MovingImageType, TransformedMovingImageType > TransformMovingImageFilterType
Definition: itkGradientDifferenceImageToImageMetric2.h:131
itk::GradientDifferenceImageToImageMetric::m_TransformMovingImageFilter
TransformMovingImageFilterType::Pointer m_TransformMovingImageFilter
Definition: itkGradientDifferenceImageToImageMetric2.h:212
itk::AdvancedImageToImageMetric::FixedImageType
Superclass::FixedImageType FixedImageType
Definition: itkAdvancedImageToImageMetric.h:105
itk::AdvancedImageToImageMetric< MetricBase< TElastix >::FixedImageType, MetricBase< TElastix >::MovingImageType >::ScalarType
TransformType::ScalarType ScalarType
Definition: itkAdvancedImageToImageMetric.h:152
itk::GradientDifferenceImageToImageMetric::m_MovedSobelOperators
SobelOperator< MovedGradientPixelType, itkGetStaticConstMacro(MovedImageDimension) > m_MovedSobelOperators[MovedImageDimension]
Definition: itkGradientDifferenceImageToImageMetric2.h:232
itk::GradientDifferenceImageToImageMetric::ComputeMeasure
MeasureType ComputeMeasure(const TransformParametersType &parameters, const double *subtractionFactor) const
itk::GradientDifferenceImageToImageMetric::m_MaxFixedGradient
FixedGradientPixelType m_MaxFixedGradient[FixedImageDimension]
Definition: itkGradientDifferenceImageToImageMetric2.h:209
itkAdvancedRayCastInterpolateImageFunction.h
itk::GradientDifferenceImageToImageMetric::MovedGradientPixelType
MovedGradientImageType::PixelType MovedGradientPixelType
Definition: itkGradientDifferenceImageToImageMetric2.h:146
itk::GradientDifferenceImageToImageMetric::CombinationTransformType
itk::AdvancedCombinationTransform< ScalarType, FixedImageDimension > CombinationTransformType
Definition: itkGradientDifferenceImageToImageMetric2.h:126
itk::GradientDifferenceImageToImageMetric::itkStaticConstMacro
itkStaticConstMacro(FixedImageDimension, unsigned int, FixedImageType::ImageDimension)
itk::GradientDifferenceImageToImageMetric::m_MovedSobelFilters
MovedSobelFilter::Pointer m_MovedSobelFilters[itkGetStaticConstMacro(MovedImageDimension)]
Definition: itkGradientDifferenceImageToImageMetric2.h:235
itk::GradientDifferenceImageToImageMetric::PrintSelf
void PrintSelf(std::ostream &os, Indent indent) const
itk::GradientDifferenceImageToImageMetric::ConstPointer
SmartPointer< const Self > ConstPointer
Definition: itkGradientDifferenceImageToImageMetric2.h:85
itk::AdvancedImageToImageMetric::RealType
Superclass::RealType RealType
Definition: itkAdvancedImageToImageMetric.h:117
itk::GradientDifferenceImageToImageMetric::MeasureType
Superclass::MeasureType MeasureType
Definition: itkGradientDifferenceImageToImageMetric2.h:108
ZeroFluxNeumannBoundaryCondition< FixedGradientImageType >
itk::GradientDifferenceImageToImageMetric::itkStaticConstMacro
itkStaticConstMacro(MovedImageDimension, unsigned int, MovingImageType::ImageDimension)
itk::GradientDifferenceImageToImageMetric::Pointer
SmartPointer< Self > Pointer
Definition: itkGradientDifferenceImageToImageMetric2.h:84
itk::GradientDifferenceImageToImageMetric::m_MinMovedGradient
MovedGradientPixelType m_MinMovedGradient[MovedImageDimension]
Definition: itkGradientDifferenceImageToImageMetric2.h:204
itk::GradientDifferenceImageToImageMetric::MovingImageType
Superclass::MovingImageType MovingImageType
Definition: itkGradientDifferenceImageToImageMetric2.h:111
itk::AdvancedImageToImageMetric::MeasureType
Superclass::MeasureType MeasureType
Definition: itkAdvancedImageToImageMetric.h:127
itk::GradientDifferenceImageToImageMetric::m_MaxMovedGradient
MovedGradientPixelType m_MaxMovedGradient[MovedImageDimension]
Definition: itkGradientDifferenceImageToImageMetric2.h:205
itk
Definition: itkAdvancedImageToImageMetric.h:40
itk::GradientDifferenceImageToImageMetric::m_FixedBoundCond
ZeroFluxNeumannBoundaryCondition< FixedGradientImageType > m_FixedBoundCond
Definition: itkGradientDifferenceImageToImageMetric2.h:225
itk::GradientDifferenceImageToImageMetric::MovedGradientImageType
itk::Image< RealType, itkGetStaticConstMacro(MovedImageDimension) > MovedGradientImageType
Definition: itkGradientDifferenceImageToImageMetric2.h:142
ZeroFluxNeumannBoundaryCondition< MovedGradientImageType >
itk::GradientDifferenceImageToImageMetric::MovedSobelFilter
NeighborhoodOperatorImageFilter< MovedGradientImageType, MovedGradientImageType > MovedSobelFilter
Definition: itkGradientDifferenceImageToImageMetric2.h:193
itk::GradientDifferenceImageToImageMetric::TransformPointer
Superclass::TransformPointer TransformPointer
Definition: itkGradientDifferenceImageToImageMetric2.h:103
itk::GradientDifferenceImageToImageMetric::CastMovedImageFilterType
itk::CastImageFilter< TransformedMovingImageType, MovedGradientImageType > CastMovedImageFilterType
Definition: itkGradientDifferenceImageToImageMetric2.h:144
itk::GradientDifferenceImageToImageMetric::m_MovedBoundCond
ZeroFluxNeumannBoundaryCondition< MovedGradientImageType > m_MovedBoundCond
Definition: itkGradientDifferenceImageToImageMetric2.h:224
itk::AdvancedImageToImageMetric::MovingImageType
Superclass::MovingImageType MovingImageType
Definition: itkAdvancedImageToImageMetric.h:101
itk::GradientDifferenceImageToImageMetric::m_CastMovedImageFilter
CastMovedImageFilterPointer m_CastMovedImageFilter
Definition: itkGradientDifferenceImageToImageMetric2.h:228
itk::GradientDifferenceImageToImageMetric::m_DerivativeDelta
double m_DerivativeDelta
Definition: itkGradientDifferenceImageToImageMetric2.h:238
itk::GradientDifferenceImageToImageMetric::ComputeMovedGradientRange
void ComputeMovedGradientRange(void) const
itk::GradientDifferenceImageToImageMetric::InterpolatorType
Superclass::InterpolatorType InterpolatorType
Definition: itkGradientDifferenceImageToImageMetric2.h:106
itk::GradientDifferenceImageToImageMetric::TransformJacobianType
Superclass::TransformJacobianType TransformJacobianType
Definition: itkGradientDifferenceImageToImageMetric2.h:105
itk::GradientDifferenceImageToImageMetric::operator=
void operator=(const Self &)
itk::GradientDifferenceImageToImageMetric::GradientDifferenceImageToImageMetric
GradientDifferenceImageToImageMetric(const Self &)
itk::AdvancedImageToImageMetric::TransformPointer
Superclass::TransformPointer TransformPointer
Definition: itkAdvancedImageToImageMetric.h:110
itk::AdvancedImageToImageMetric::TransformType
Superclass::TransformType TransformType
Definition: itkAdvancedImageToImageMetric.h:109
itk::GradientDifferenceImageToImageMetric::RayCastInterpolatorType
itk::AdvancedRayCastInterpolateImageFunction< MovingImageType, ScalarType > RayCastInterpolatorType
Definition: itkGradientDifferenceImageToImageMetric2.h:133
itk::GradientDifferenceImageToImageMetric::m_FixedSobelOperators
SobelOperator< FixedGradientPixelType, itkGetStaticConstMacro(FixedImageDimension) > m_FixedSobelOperators[FixedImageDimension]
Definition: itkGradientDifferenceImageToImageMetric2.h:219
itkAdvancedImageToImageMetric.h
itk::GradientDifferenceImageToImageMetric::CombinationTransformPointer
CombinationTransformType::Pointer CombinationTransformPointer
Definition: itkGradientDifferenceImageToImageMetric2.h:127
itk::AdvancedImageToImageMetric::FixedImageConstPointer
Superclass::FixedImageConstPointer FixedImageConstPointer
Definition: itkAdvancedImageToImageMetric.h:107
itk::AdvancedCombinationTransform::Pointer
SmartPointer< Self > Pointer
Definition: itkAdvancedCombinationTransform.h:65
itk::GradientDifferenceImageToImageMetric::Self
GradientDifferenceImageToImageMetric Self
Definition: itkGradientDifferenceImageToImageMetric2.h:81
itk::AdvancedImageToImageMetric::MovingImageConstPointer
Superclass::MovingImageConstPointer MovingImageConstPointer
Definition: itkAdvancedImageToImageMetric.h:104
itk::GradientDifferenceImageToImageMetric::m_Rescalingfactor
double m_Rescalingfactor
Definition: itkGradientDifferenceImageToImageMetric2.h:239
itk::GradientDifferenceImageToImageMetric::GetValue
MeasureType GetValue(const TransformParametersType &parameters) const


Generated on OURCE_DATE_EPOCH for elastix by doxygen 1.8.18 elastix logo