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

itkMultiScaleHessianBasedMeasureImageFilter.h

Go to the documentation of this file.
00001 /*=========================================================================
00002 
00003   Program:   Insight Segmentation & Registration Toolkit
00004   Module:    $RCSfile: itkMultiScaleHessianBasedMeasureImageFilter.h,v $
00005   Language:  C++
00006   Date:      $Date: 2010-02-04 12:16:43 $
00007   Version:   $Revision: 1.10 $
00008 
00009   Copyright (c) Insight Software 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 __itkMultiScaleHessianBasedMeasureImageFilter_h
00018 #define __itkMultiScaleHessianBasedMeasureImageFilter_h
00019 
00020 #include "itkImageToImageFilter.h"
00021 #include "itkImage.h"
00022 #include "itkHessianRecursiveGaussianImageFilter.h"
00023 
00024 namespace itk
00025 {
00055 template <typename TInputImage,
00056           typename THessianImage, 
00057           typename TOutputImage=TInputImage >
00058 class ITK_EXPORT MultiScaleHessianBasedMeasureImageFilter 
00059 : public
00060 ImageToImageFilter< TInputImage,TOutputImage > 
00061 {
00062 public:
00064   typedef MultiScaleHessianBasedMeasureImageFilter          Self;
00065   typedef ImageToImageFilter<TInputImage,TOutputImage>      Superclass;
00066 
00067   typedef SmartPointer<Self>                                      Pointer;
00068   typedef SmartPointer<const Self>                                ConstPointer;
00069 
00070   typedef TInputImage                                    InputImageType;
00071   typedef TOutputImage                                   OutputImageType;
00072   typedef THessianImage                                  HessianImageType;
00073 
00074   typedef ImageToImageFilter< HessianImageType, OutputImageType  >  HessianToMeasureFilterType;
00075 
00076   typedef typename TInputImage::PixelType                InputPixelType;
00077   typedef typename TOutputImage::PixelType               OutputPixelType;
00078   typedef typename TOutputImage::RegionType              OutputRegionType;
00079  
00081   itkStaticConstMacro(ImageDimension, unsigned int, ::itk::GetImageDimension<InputImageType>::ImageDimension);
00082 
00084   typedef float                                          ScalesPixelType;
00085   typedef Image<ScalesPixelType, itkGetStaticConstMacro(ImageDimension)>  ScalesImageType;
00086 
00088   typedef HessianRecursiveGaussianImageFilter< InputImageType, HessianImageType> HessianFilterType;
00089 
00093   typedef Image< double, itkGetStaticConstMacro(ImageDimension) > UpdateBufferType;
00094   typedef typename UpdateBufferType::ValueType                    BufferValueType;
00095 
00096   typedef typename Superclass::DataObjectPointer          DataObjectPointer;
00097   
00099   itkNewMacro(Self);
00100 
00102   itkTypeMacro(MultiScaleHessianBasedMeasureImageFilter, 
00103                  ImageToImageFilter);
00104 
00106   itkSetMacro(SigmaMinimum, double);
00107   itkGetConstMacro(SigmaMinimum, double);
00109 
00111   itkSetMacro(SigmaMaximum, double);
00112   itkGetConstMacro(SigmaMaximum, double);
00114 
00116   itkSetMacro(NumberOfSigmaSteps, unsigned int);
00117   itkGetConstMacro(NumberOfSigmaSteps, unsigned int);
00119 
00123   itkSetObjectMacro( HessianToMeasureFilter, HessianToMeasureFilterType); 
00124   itkGetObjectMacro( HessianToMeasureFilter, HessianToMeasureFilterType); 
00126 
00133   itkSetMacro(NonNegativeHessianBasedMeasure,bool);
00134   itkGetConstMacro(NonNegativeHessianBasedMeasure,bool);
00135   itkBooleanMacro(NonNegativeHessianBasedMeasure);
00137 
00138   typedef enum { EquispacedSigmaSteps = 0,
00139                  LogarithmicSigmaSteps = 1 } SigmaStepMethodType;
00140 
00143   itkSetMacro(SigmaStepMethod, SigmaStepMethodType);
00144   itkGetConstMacro(SigmaStepMethod, SigmaStepMethodType);
00146 
00148   void SetSigmaStepMethodToEquispaced();
00149 
00151   void SetSigmaStepMethodToLogarithmic();
00152 
00155   const HessianImageType* GetHessianOutput() const;
00156 
00159   const ScalesImageType* GetScalesOutput() const;
00160 
00161   void EnlargeOutputRequestedRegion (DataObject *);
00162 
00165   itkSetMacro(GenerateScalesOutput,bool);
00166   itkGetConstMacro(GenerateScalesOutput,bool);
00167   itkBooleanMacro(GenerateScalesOutput);
00169 
00172   itkSetMacro(GenerateHessianOutput,bool);
00173   itkGetConstMacro(GenerateHessianOutput,bool);
00174   itkBooleanMacro(GenerateHessianOutput);
00176 
00178   virtual DataObjectPointer MakeOutput(unsigned int idx);
00179 
00180 protected:
00181   MultiScaleHessianBasedMeasureImageFilter();
00182   ~MultiScaleHessianBasedMeasureImageFilter() {};
00183   void PrintSelf(std::ostream& os, Indent indent) const;
00184   
00186   void GenerateData( void );
00187 
00188 private:
00189   void UpdateMaximumResponse(double sigma);
00190   double ComputeSigmaValue(int scaleLevel);
00191   
00192   void AllocateUpdateBuffer();
00193 
00194   //purposely not implemented
00195   MultiScaleHessianBasedMeasureImageFilter(const Self&); 
00196   void operator=(const Self&); //purposely not implemented
00197 
00198   bool                        m_NonNegativeHessianBasedMeasure;
00199 
00200   double                      m_SigmaMinimum;
00201   double                      m_SigmaMaximum;
00202 
00203   unsigned int                m_NumberOfSigmaSteps;
00204   SigmaStepMethodType         m_SigmaStepMethod;
00205 
00206   typename HessianToMeasureFilterType::Pointer  m_HessianToMeasureFilter;
00207   typename HessianFilterType::Pointer           m_HessianFilter;
00208   typename UpdateBufferType::Pointer            m_UpdateBuffer;
00209 
00210   bool m_GenerateScalesOutput;
00211   bool m_GenerateHessianOutput;
00212 };
00213 
00214 } // end namespace itk
00215 
00216 #ifndef ITK_MANUAL_INSTANTIATION
00217 #include "itkMultiScaleHessianBasedMeasureImageFilter.txx"
00218 #endif
00219   
00220 #endif
00221 

Generated at Mon Jul 12 2010 19:14:31 for ITK by doxygen 1.7.1 written by Dimitri van Heesch, © 1997-2000