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

itkFastMarchingUpwindGradientImageFilter.h

Go to the documentation of this file.
00001 /*=========================================================================
00002 
00003   Program:   Insight Segmentation & Registration Toolkit
00004   Module:    $RCSfile: itkFastMarchingUpwindGradientImageFilter.h,v $
00005   Language:  C++
00006   Date:      $Date: 2007/04/20 11:52:27 $
00007   Version:   $Revision: 1.7 $
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 
00018 #ifndef __itkFastMarchingUpwindGradientImageFilter_h
00019 #define __itkFastMarchingUpwindGradientImageFilter_h
00020 
00021 #include "itkFastMarchingImageFilter.h"
00022 #include "itkImage.h"
00023 
00024 namespace itk
00025 {
00057 template <
00058   class TLevelSet, 
00059   class TSpeedImage = Image<float,::itk::GetImageDimension<TLevelSet>::ImageDimension> >
00060 class ITK_EXPORT FastMarchingUpwindGradientImageFilter :
00061     public FastMarchingImageFilter<TLevelSet,TSpeedImage>
00062 {
00063 
00064 public:
00066   typedef FastMarchingUpwindGradientImageFilter          Self;
00067   typedef FastMarchingImageFilter<TLevelSet,TSpeedImage> Superclass;
00068   typedef SmartPointer<Self>                             Pointer;
00069   typedef SmartPointer<const Self>                       ConstPointer;
00070 
00072   itkNewMacro(Self);
00073 
00075   itkTypeMacro(FastMarchingUpwindGradientImageFilter, FastMarchingImageFilter);
00076 
00078   typedef typename Superclass::LevelSetType            LevelSetType;
00079   typedef typename Superclass::SpeedImageType          SpeedImageType;
00080   typedef typename Superclass::LevelSetImageType       LevelSetImageType;
00081   typedef typename Superclass::LevelSetPointer         LevelSetPointer;
00082   typedef typename Superclass::SpeedImageConstPointer  SpeedImageConstPointer;
00083   typedef typename Superclass::LabelImageType          LabelImageType;
00084   typedef typename Superclass::PixelType               PixelType;
00085   typedef typename Superclass::AxisNodeType            AxisNodeType;
00086   typedef typename Superclass::NodeType                NodeType;
00087   typedef typename Superclass::NodeContainer           NodeContainer;
00088   typedef typename Superclass::NodeContainerPointer    NodeContainerPointer;
00089 
00090   typedef typename Superclass::IndexType          IndexType;
00091   typedef typename Superclass::OutputSpacingType  OutputSpacingType;
00092   typedef typename Superclass::LevelSetIndexType  LevelSetIndexType;
00093 
00094   typedef typename Superclass::OutputPointType   PointType;
00095 
00096 
00098   itkStaticConstMacro(SetDimension, unsigned int,Superclass::SetDimension); 
00099 
00103   void SetTargetPoints( NodeContainer * points )
00104     { 
00105     m_TargetPoints = points;
00106     this->Modified();
00107     };
00109 
00111   NodeContainerPointer GetTargetPoints( )
00112     { return m_TargetPoints; }
00113 
00115   NodeContainerPointer GetReachedTargetPoints( )
00116     { return m_ReachedTargetPoints; }
00117 
00119   typedef CovariantVector<PixelType, 
00120           itkGetStaticConstMacro(SetDimension)> GradientPixelType;
00121 
00123   typedef Image<GradientPixelType, 
00124           itkGetStaticConstMacro(SetDimension)> GradientImageType;
00125 
00127   typedef typename GradientImageType::Pointer GradientImagePointer;
00128 
00130   GradientImagePointer GetGradientImage() const
00131     { return m_GradientImage; }
00132 
00135   itkSetMacro( GenerateGradientImage, bool );
00136 
00138   itkGetConstReferenceMacro( GenerateGradientImage, bool );
00139   itkBooleanMacro( GenerateGradientImage );
00141 
00145   itkSetMacro( TargetOffset, double );
00146 
00148   itkGetConstReferenceMacro( TargetOffset, double );
00149 
00153   itkSetMacro( TargetReachedMode, int );
00154   itkGetConstReferenceMacro( TargetReachedMode, int );
00155   void SetTargetReachedModeToNoTargets()
00156     { this->SetTargetReachedMode(NoTargets); }
00157   void SetTargetReachedModeToOneTarget() 
00158     { this->SetTargetReachedMode(OneTarget); }
00159   void SetTargetReachedModeToSomeTargets( long numberOfTargets )
00160     {
00161     this->SetTargetReachedMode(SomeTargets);
00162     m_NumberOfTargets = numberOfTargets;
00163     }
00164   void SetTargetReachedModeToAllTargets() 
00165     { this->SetTargetReachedMode(AllTargets); }
00167 
00169   itkGetConstReferenceMacro( NumberOfTargets, long );
00170 
00172   itkGetConstReferenceMacro( TargetValue, double );
00173 
00174   enum
00175     {
00176       NoTargets,
00177       OneTarget,
00178       SomeTargets,
00179       AllTargets
00180     };
00181 
00182 #ifdef ITK_USE_CONCEPT_CHECKING
00183 
00184   itkConceptMacro(LevelSetDoubleDivisionOperatorsCheck,
00185     (Concept::DivisionOperators<typename TLevelSet::PixelType, double>));
00186 
00188 #endif
00189 
00190 protected:
00191   FastMarchingUpwindGradientImageFilter();
00192   ~FastMarchingUpwindGradientImageFilter(){};
00193   void PrintSelf( std::ostream& os, Indent indent ) const;
00194 
00195   virtual void Initialize( LevelSetImageType * );
00196 
00197   void GenerateData();
00198 
00199   virtual void UpdateNeighbors( const IndexType& index, 
00200                                 const SpeedImageType *, LevelSetImageType * );
00201 
00202   virtual void ComputeGradient( const IndexType& index ,
00203                               const LevelSetImageType * output, 
00204                               const LabelImageType * labelImage,
00205                               GradientImageType * gradientImage);
00206 
00207 private:
00208   FastMarchingUpwindGradientImageFilter(const Self&); //purposely not implemented
00209   void operator=(const Self&); //purposely not implemented
00210   
00211   NodeContainerPointer m_TargetPoints;
00212   NodeContainerPointer m_ReachedTargetPoints;
00213 
00214   GradientImagePointer m_GradientImage;
00215  
00216   bool m_GenerateGradientImage;
00217 
00218   double m_TargetOffset;
00219 
00220   int m_TargetReachedMode;
00221 
00222   double m_TargetValue;
00223 
00224   long m_NumberOfTargets;
00225 
00226 };
00227 
00228 } // namespace itk
00229 
00230 
00231 #ifndef ITK_MANUAL_INSTANTIATION
00232 #include "itkFastMarchingUpwindGradientImageFilter.txx"
00233 #endif
00234 
00235 #endif
00236 

Generated at Sun Sep 23 12:28:22 2007 for ITK by doxygen 1.5.1 written by Dimitri van Heesch, © 1997-2000