Main Page   Groups   Namespace List   Class Hierarchy   Alphabetical List   Compound List   File List   Namespace Members   Compound Members   File Members   Concepts

itkGradientDifferenceImageToImageMetric.h

Go to the documentation of this file.
00001 /*=========================================================================
00002 
00003   Program:   Insight Segmentation & Registration Toolkit
00004   Module:    $RCSfile: itkGradientDifferenceImageToImageMetric.h,v $
00005   Language:  C++
00006   Date:      $Date: 2004/07/12 16:17:30 $
00007   Version:   $Revision: 1.6 $
00008 
00009   Copyright (c) 2002 Insight Consortium. All rights reserved.
00010   See ITKCopyright.txt or http://www.itk.org/HTML/Copyright.htm for details.
00011 
00012      This software is distributed WITHOUT ANY WARRANTY; without even 
00013      the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR 
00014      PURPOSE.  See the above copyright notices for more information.
00015 
00016 =========================================================================*/
00017 #ifndef __itkGradientDifferenceImageToImageMetric_h
00018 #define __itkGradientDifferenceImageToImageMetric_h
00019 
00020 #include "itkImageToImageMetric.h"
00021 
00022 #include "itkSobelOperator.h"
00023 #include "itkNeighborhoodOperatorImageFilter.h"
00024 #include "itkPoint.h"
00025 #include "itkCastImageFilter.h"
00026 #include "itkResampleImageFilter.h"
00027 
00028 namespace itk
00029 {
00048 template < class TFixedImage, class TMovingImage > 
00049 class ITK_EXPORT GradientDifferenceImageToImageMetric : 
00050 public ImageToImageMetric< TFixedImage, TMovingImage>
00051 {
00052 public:
00053 
00055   typedef GradientDifferenceImageToImageMetric    Self;
00056   typedef ImageToImageMetric<TFixedImage, TMovingImage >  Superclass;
00057 
00058   typedef SmartPointer<Self>         Pointer;
00059   typedef SmartPointer<const Self>   ConstPointer;
00060 
00062   itkNewMacro(Self);
00063  
00065   itkTypeMacro(GradientDifferenceImageToImageMetric, ImageToImageMetric);
00066 
00067  
00069 // Work around a Visual Studio .NET bug
00070 #if defined(_MSC_VER) && (_MSC_VER == 1300)
00071   typedef double RealType;
00072 #else
00073   typedef typename Superclass::RealType                 RealType;
00074 #endif
00075   typedef typename Superclass::TransformType            TransformType;
00076   typedef typename Superclass::TransformPointer         TransformPointer;
00077   typedef typename Superclass::TransformParametersType  TransformParametersType;
00078   typedef typename Superclass::TransformJacobianType    TransformJacobianType;
00079 
00080   typedef typename Superclass::MeasureType              MeasureType;
00081   typedef typename Superclass::DerivativeType           DerivativeType;
00082   typedef typename Superclass::FixedImageType           FixedImageType;
00083   typedef typename Superclass::MovingImageType          MovingImageType;
00084   typedef typename Superclass::FixedImageConstPointer   FixedImageConstPointer;
00085   typedef typename Superclass::MovingImageConstPointer  MovingImageConstPointer;
00086 
00087   typedef typename TFixedImage::PixelType               FixedImagePixelType;
00088   typedef typename TMovingImage::PixelType              MovedImagePixelType;
00089 
00090   itkStaticConstMacro(FixedImageDimension, unsigned int, TFixedImage::ImageDimension);
00092   typedef itk::Image< FixedImagePixelType, 
00093                       itkGetStaticConstMacro( FixedImageDimension ) >      
00094                     TransformedMovingImageType;
00095 
00096   typedef itk::ResampleImageFilter< MovingImageType, TransformedMovingImageType >  
00097                                                         TransformMovingImageFilterType;
00098 
00101   typedef itk::Image< RealType, 
00102                       itkGetStaticConstMacro( FixedImageDimension ) >
00103                     FixedGradientImageType;
00104 
00105   typedef itk::CastImageFilter< FixedImageType, FixedGradientImageType > 
00106                                                       CastFixedImageFilterType;
00107   typedef typename CastFixedImageFilterType::Pointer CastFixedImageFilterPointer;
00108 
00109   typedef typename FixedGradientImageType::PixelType FixedGradientPixelType;
00110 
00111 
00114   itkStaticConstMacro( MovedImageDimension, unsigned int,
00115                        MovingImageType::ImageDimension );
00116 
00117   typedef itk::Image< RealType,
00118                       itkGetStaticConstMacro( MovedImageDimension ) >
00119                     MovedGradientImageType;
00120 
00121   typedef itk::CastImageFilter< TransformedMovingImageType, MovedGradientImageType > 
00122                                                      CastMovedImageFilterType;
00123   typedef typename CastMovedImageFilterType::Pointer CastMovedImageFilterPointer;
00124 
00125   typedef typename MovedGradientImageType::PixelType MovedGradientPixelType;
00126 
00127 
00129   void GetDerivative( const TransformParametersType & parameters,
00130                             DerivativeType  & derivative ) const;
00131 
00133   MeasureType GetValue( const TransformParametersType & parameters ) const;
00134 
00136   void GetValueAndDerivative( const TransformParametersType & parameters,
00137                               MeasureType& Value, DerivativeType& derivative ) const;
00138 
00141   virtual void Initialize(void) throw ( ExceptionObject );
00142 
00144   void WriteGradientImagesToFiles(void) const;
00145 
00146 protected:
00147   GradientDifferenceImageToImageMetric();
00148   virtual ~GradientDifferenceImageToImageMetric() {};
00149   void PrintSelf(std::ostream& os, Indent indent) const;
00150 
00152   void ComputeMovedGradientRange( void ) const;
00153 
00155   void ComputeVariance( void ) const;
00156 
00158   MeasureType ComputeMeasure( const TransformParametersType &parameters,
00159                               const double *subtractionFactor ) const;
00160 
00161   typedef NeighborhoodOperatorImageFilter<
00162     FixedGradientImageType, FixedGradientImageType > FixedSobelFilter;
00163 
00164   typedef NeighborhoodOperatorImageFilter<
00165     MovedGradientImageType, MovedGradientImageType > MovedSobelFilter;
00166 
00167 private:
00168   GradientDifferenceImageToImageMetric(const Self&); //purposely not implemented
00169   void operator=(const Self&); //purposely not implemented
00170 
00172   mutable MovedGradientPixelType m_Variance[FixedImageDimension];
00173 
00175   mutable MovedGradientPixelType m_MinMovedGradient[MovedImageDimension];
00176   mutable MovedGradientPixelType m_MaxMovedGradient[MovedImageDimension];
00178   mutable FixedGradientPixelType m_MinFixedGradient[FixedImageDimension];
00179   mutable FixedGradientPixelType m_MaxFixedGradient[FixedImageDimension];
00180 
00182   typename TransformMovingImageFilterType::Pointer m_TransformMovingImageFilter;
00183 
00185   CastFixedImageFilterPointer m_CastFixedImageFilter;
00186 
00187   SobelOperator< FixedGradientPixelType, 
00188                  itkGetStaticConstMacro(FixedImageDimension) >
00189                m_FixedSobelOperators[FixedImageDimension];
00190 
00191   typename FixedSobelFilter::Pointer m_FixedSobelFilters[itkGetStaticConstMacro( FixedImageDimension )];
00192 
00193 
00195   CastMovedImageFilterPointer m_CastMovedImageFilter;
00196 
00197   SobelOperator< MovedGradientPixelType, 
00198                  itkGetStaticConstMacro(MovedImageDimension) >
00199                m_MovedSobelOperators[MovedImageDimension];
00200 
00201   typename MovedSobelFilter::Pointer m_MovedSobelFilters[itkGetStaticConstMacro( MovedImageDimension )];
00202 
00203 };
00204 
00205 } // end namespace itk
00206 
00207 #ifndef ITK_MANUAL_INSTANTIATION
00208 #include "itkGradientDifferenceImageToImageMetric.txx"
00209 #endif
00210 
00211 #endif
00212 
00213 
00214 

Generated at Tue Mar 29 23:52:43 2005 for ITK by doxygen 1.3.9.1 written by Dimitri van Heesch, © 1997-2000