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: 2005/12/03 22:15:25 $
00007   Version:   $Revision: 1.6.10.1 $
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 {
00053 template < class TFixedImage, class TMovingImage > 
00054 class ITK_EXPORT GradientDifferenceImageToImageMetric : 
00055 public ImageToImageMetric< TFixedImage, TMovingImage>
00056 {
00057 public:
00058 
00060   typedef GradientDifferenceImageToImageMetric    Self;
00061   typedef ImageToImageMetric<TFixedImage, TMovingImage >  Superclass;
00062 
00063   typedef SmartPointer<Self>         Pointer;
00064   typedef SmartPointer<const Self>   ConstPointer;
00065 
00067   itkNewMacro(Self);
00068  
00070   itkTypeMacro(GradientDifferenceImageToImageMetric, ImageToImageMetric);
00071 
00072  
00074 // Work around a Visual Studio .NET bug
00075 #if defined(_MSC_VER) && (_MSC_VER == 1300)
00076   typedef double RealType;
00077 #else
00078   typedef typename Superclass::RealType                 RealType;
00079 #endif
00080   typedef typename Superclass::TransformType            TransformType;
00081   typedef typename Superclass::TransformPointer         TransformPointer;
00082   typedef typename Superclass::TransformParametersType  TransformParametersType;
00083   typedef typename Superclass::TransformJacobianType    TransformJacobianType;
00084 
00085   typedef typename Superclass::MeasureType              MeasureType;
00086   typedef typename Superclass::DerivativeType           DerivativeType;
00087   typedef typename Superclass::FixedImageType           FixedImageType;
00088   typedef typename Superclass::MovingImageType          MovingImageType;
00089   typedef typename Superclass::FixedImageConstPointer   FixedImageConstPointer;
00090   typedef typename Superclass::MovingImageConstPointer  MovingImageConstPointer;
00091 
00092   typedef typename TFixedImage::PixelType               FixedImagePixelType;
00093   typedef typename TMovingImage::PixelType              MovedImagePixelType;
00094 
00095   itkStaticConstMacro(FixedImageDimension, unsigned int, TFixedImage::ImageDimension);
00097   typedef itk::Image< FixedImagePixelType, 
00098                       itkGetStaticConstMacro( FixedImageDimension ) >      
00099                     TransformedMovingImageType;
00100 
00101   typedef itk::ResampleImageFilter< MovingImageType, TransformedMovingImageType >  
00102                                                         TransformMovingImageFilterType;
00103 
00106   typedef itk::Image< RealType, 
00107                       itkGetStaticConstMacro( FixedImageDimension ) >
00108                     FixedGradientImageType;
00109 
00110   typedef itk::CastImageFilter< FixedImageType, FixedGradientImageType > 
00111                                                       CastFixedImageFilterType;
00112   typedef typename CastFixedImageFilterType::Pointer CastFixedImageFilterPointer;
00113 
00114   typedef typename FixedGradientImageType::PixelType FixedGradientPixelType;
00115 
00116 
00119   itkStaticConstMacro( MovedImageDimension, unsigned int,
00120                        MovingImageType::ImageDimension );
00121 
00122   typedef itk::Image< RealType,
00123                       itkGetStaticConstMacro( MovedImageDimension ) >
00124                     MovedGradientImageType;
00125 
00126   typedef itk::CastImageFilter< TransformedMovingImageType, MovedGradientImageType > 
00127                                                      CastMovedImageFilterType;
00128   typedef typename CastMovedImageFilterType::Pointer CastMovedImageFilterPointer;
00129 
00130   typedef typename MovedGradientImageType::PixelType MovedGradientPixelType;
00131 
00132 
00134   void GetDerivative( const TransformParametersType & parameters,
00135                             DerivativeType  & derivative ) const;
00136 
00138   MeasureType GetValue( const TransformParametersType & parameters ) const;
00139 
00141   void GetValueAndDerivative( const TransformParametersType & parameters,
00142                               MeasureType& Value, DerivativeType& derivative ) const;
00143 
00146   virtual void Initialize(void) throw ( ExceptionObject );
00147 
00149   void WriteGradientImagesToFiles(void) const;
00150 
00151 protected:
00152   GradientDifferenceImageToImageMetric();
00153   virtual ~GradientDifferenceImageToImageMetric() {};
00154   void PrintSelf(std::ostream& os, Indent indent) const;
00155 
00157   void ComputeMovedGradientRange( void ) const;
00158 
00160   void ComputeVariance( void ) const;
00161 
00163   MeasureType ComputeMeasure( const TransformParametersType &parameters,
00164                               const double *subtractionFactor ) const;
00165 
00166   typedef NeighborhoodOperatorImageFilter<
00167     FixedGradientImageType, FixedGradientImageType > FixedSobelFilter;
00168 
00169   typedef NeighborhoodOperatorImageFilter<
00170     MovedGradientImageType, MovedGradientImageType > MovedSobelFilter;
00171 
00172 private:
00173   GradientDifferenceImageToImageMetric(const Self&); //purposely not implemented
00174   void operator=(const Self&); //purposely not implemented
00175 
00177   mutable MovedGradientPixelType m_Variance[FixedImageDimension];
00178 
00180   mutable MovedGradientPixelType m_MinMovedGradient[MovedImageDimension];
00181   mutable MovedGradientPixelType m_MaxMovedGradient[MovedImageDimension];
00183   mutable FixedGradientPixelType m_MinFixedGradient[FixedImageDimension];
00184   mutable FixedGradientPixelType m_MaxFixedGradient[FixedImageDimension];
00185 
00187   typename TransformMovingImageFilterType::Pointer m_TransformMovingImageFilter;
00188 
00190   CastFixedImageFilterPointer m_CastFixedImageFilter;
00191 
00192   SobelOperator< FixedGradientPixelType, 
00193                  itkGetStaticConstMacro(FixedImageDimension) >
00194                m_FixedSobelOperators[FixedImageDimension];
00195 
00196   typename FixedSobelFilter::Pointer m_FixedSobelFilters[itkGetStaticConstMacro( FixedImageDimension )];
00197 
00198   ZeroFluxNeumannBoundaryCondition< MovedGradientImageType > m_MovedBoundCond;
00199   ZeroFluxNeumannBoundaryCondition< FixedGradientImageType > m_FixedBoundCond;
00200 
00202   CastMovedImageFilterPointer m_CastMovedImageFilter;
00203 
00204   SobelOperator< MovedGradientPixelType, 
00205                  itkGetStaticConstMacro(MovedImageDimension) >
00206                m_MovedSobelOperators[MovedImageDimension];
00207 
00208   typename MovedSobelFilter::Pointer m_MovedSobelFilters[itkGetStaticConstMacro( MovedImageDimension )];
00209 
00210 };
00211 
00212 } // end namespace itk
00213 
00214 #ifndef ITK_MANUAL_INSTANTIATION
00215 #include "itkGradientDifferenceImageToImageMetric.txx"
00216 #endif
00217 
00218 #endif
00219 
00220 
00221 

Generated at Wed May 24 23:12:40 2006 for ITK by doxygen 1.3.5 written by Dimitri van Heesch, © 1997-2000