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

itkWatershedSegmenter.h

Go to the documentation of this file.
00001 /*=========================================================================
00002 
00003   Program:   Insight Segmentation & Registration Toolkit
00004   Module:    $RCSfile: itkWatershedSegmenter.h,v $
00005   Language:  C++
00006   Date:      $Date: 2002/09/11 19:39:01 $
00007   Version:   $Revision: 1.8 $
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 __itkWatershedSegmenter_h
00018 #define __itkWatershedSegmenter_h
00019 
00020 #if defined(_MSC_VER)
00021 #pragma warning ( disable : 4786 )
00022 #endif
00023 
00024 #include "itk_hash_map.h"
00025 #include "itkWatershedBoundary.h"
00026 #include "itkWatershedSegmentTable.h"
00027 #include "itkWatershedEquivalencyTable.h"
00028 #include "itkImage.h"
00029 
00030 namespace itk
00031 {
00032 namespace watershed
00033 {
00089 template <class TInputImage>
00090 class Segmenter
00091   : public ProcessObject
00092 {
00093 public:
00095   typedef Segmenter Self;
00096 
00098   typedef TInputImage InputImageType;
00099   itkStaticConstMacro(ImageDimension, unsigned int,
00100                       TInputImage::ImageDimension);
00101   typedef Image<unsigned long, itkGetStaticConstMacro(ImageDimension)> OutputImageType;
00102   typedef typename InputImageType::RegionType ImageRegionType;
00103   typedef typename InputImageType::PixelType  InputPixelType;
00104   typedef Boundary<InputPixelType, itkGetStaticConstMacro(ImageDimension)> BoundaryType;
00105   typedef typename BoundaryType::IndexType BoundaryIndexType;
00106   typedef typename BoundaryType::FlatHashValueType BoundaryFlatHashValueType;
00107   typedef SegmentTable<InputPixelType> SegmentTableType;
00108   typedef DataObject::Pointer DataObjectPointer;
00109   
00112   typedef ProcessObject Superclass;
00113   typedef SmartPointer<Self> Pointer;
00114   typedef SmartPointer<const Self> ConstPointer;
00115   itkNewMacro(Self);
00116   itkTypeMacro(Segmenter, ProcessObject);
00117 
00119   typedef typename InputImageType::Pointer InputImageTypePointer;
00120   typedef typename OutputImageType::Pointer OutputImageTypePointer;
00121   typedef typename SegmentTableType::Pointer SegmentTableTypePointer;
00122   typedef typename BoundaryType::Pointer BoundaryTypePointer;
00123     
00125   static unsigned long NULL_LABEL;
00126 
00128   static short NULL_FLOW;
00129  
00131   InputImageType * GetInputImage(void)
00132     { return static_cast<InputImageType *>
00133         (this->ProcessObject::GetInput(0));    }
00134   void SetInputImage(InputImageType *img)
00135     {  this->ProcessObject::SetNthInput(0, img); }
00136 
00139   OutputImageType * GetOutputImage(void)
00140     { return static_cast<OutputImageType *>
00141         (this->ProcessObject::GetOutput(0)); }
00142   void SetOutputImage(OutputImageType *img)
00143     { this->ProcessObject::SetNthOutput(0, img);    }
00144   
00147   SegmentTableType * GetSegmentTable(void)
00148     { return static_cast<SegmentTableType *>
00149         (this->ProcessObject::GetOutput(1)); }
00150   void SetSegmentTable(SegmentTableType *s)
00151     { this->ProcessObject::SetNthOutput(1, s); }
00152   
00155   BoundaryType * GetBoundary(void)
00156     { return static_cast<BoundaryType *>
00157         (this->ProcessObject::GetOutput(2)); }
00158   void SetBoundary(BoundaryType *b)
00159     { this->ProcessObject::SetNthOutput(2,b); }
00160   
00162   void GenerateData();
00163 
00170   void SetLargestPossibleRegion(ImageRegionType reg)
00171   {
00172     if (reg == m_LargestPossibleRegion) return;
00173     m_LargestPossibleRegion = reg;
00174     this->Modified();
00175   }
00176   ImageRegionType GetLargestPossibleRegion() const
00177     {      return m_LargestPossibleRegion; }
00178 
00181   static void RelabelImage(OutputImageTypePointer,
00182                            ImageRegionType,
00183                            EquivalencyTable::Pointer);
00184   
00186   virtual DataObjectPointer MakeOutput(unsigned int idx);
00187 
00190   itkSetMacro(CurrentLabel, unsigned long);
00191   itkGetMacro(CurrentLabel, unsigned long);
00192 
00202   itkSetClampMacro(Threshold, double, 0.0, 1.0);
00203   itkGetMacro(Threshold, double);
00204 
00208   itkSetMacro(DoBoundaryAnalysis, bool);
00209   itkGetMacro(DoBoundaryAnalysis, bool);
00210 
00214   itkGetMacro(Minimum, InputPixelType);
00215   itkSetMacro(Minimum, InputPixelType);
00216   itkGetMacro(Maximum, InputPixelType);
00217   itkSetMacro(Maximum, InputPixelType);
00218   
00223   itkGetMacro(SortEdgeLists, bool);
00224   itkSetMacro(SortEdgeLists, bool);
00225 
00226 protected:
00229   struct flat_region_t
00230   {
00231     unsigned long   *min_label_ptr;
00232     InputPixelType  bounds_min;
00233     //    InputPixelType  bounds_max; // <-- may not be necc.
00234     InputPixelType  value;
00235     bool            is_on_boundary;
00236     flat_region_t() : is_on_boundary(false) {}
00237   };
00238 
00240   typedef itk::hash_map<unsigned long, flat_region_t, itk::hash<unsigned long> >
00241     flat_region_table_t;
00242 
00243   struct connectivity_t
00244   {
00245     unsigned int size;
00246     unsigned int *index;
00247     typename InputImageType::OffsetType *direction;
00248   };
00249 
00254   typedef itk::hash_map<unsigned long, InputPixelType, itk::hash<unsigned long> 
00255   > edge_table_t;
00256   
00257   typedef itk::hash_map<unsigned long, edge_table_t, itk::hash<unsigned long>
00258   > edge_table_hash_t;
00259   
00260   Segmenter();
00261   Segmenter(const Self&) {}
00262   virtual ~Segmenter()
00263     {
00264       if (m_Connectivity.index != 0)     delete[] m_Connectivity.index;
00265       if (m_Connectivity.direction !=0 ) delete[] m_Connectivity.direction;
00266     }
00267   void PrintSelf(std::ostream& os, Indent indent) const;
00268   void operator=(const Self&) {}
00269   
00272   virtual void GenerateConnectivity();
00273   
00277   void GenerateInputRequestedRegion();
00278   void GenerateOutputRequestedRegion(DataObject *output);
00279   void UpdateOutputInformation();
00280 
00283   void InitializeBoundary();
00284 
00288   void AnalyzeBoundaryFlow(InputImageTypePointer,
00289                            flat_region_table_t &,
00290                            InputPixelType);
00291 
00295   void BuildRetainingWall(InputImageTypePointer,
00296                           ImageRegionType, InputPixelType);
00297 
00300   void LabelMinima(InputImageTypePointer,
00301                    ImageRegionType, flat_region_table_t &,
00302                    InputPixelType);
00303 
00307   void GradientDescent(InputImageTypePointer, ImageRegionType);
00308 
00311   void DescendFlatRegions(flat_region_table_t &, ImageRegionType);
00312 
00315   void UpdateSegmentTable(InputImageTypePointer, ImageRegionType);
00316 
00320   void CollectBoundaryInformation(flat_region_table_t &);
00321   
00325   static void Threshold(InputImageTypePointer destination,
00326                         InputImageTypePointer source,
00327                         const ImageRegionType source_region,
00328                         const ImageRegionType destination_region,
00329                         InputPixelType threshold);
00330 
00332   static void MinMax(InputImageTypePointer img,
00333                      ImageRegionType region,
00334                      InputPixelType &min,
00335                      InputPixelType &max);
00336 
00338   static void MergeFlatRegions(flat_region_table_t &, EquivalencyTable::Pointer);
00339     
00341   static void SetImageValues(InputImageTypePointer img,
00342                              const ImageRegionType region,
00343                              InputPixelType value);
00344 
00345   static void SetImageValues(OutputImageTypePointer img,
00346                              const ImageRegionType region,
00347                              unsigned long value);
00348 
00350   //  bool CheckLabeledBoundaries();
00351   
00354   connectivity_t m_Connectivity;
00355 
00356 private:
00358   //  void PrintFlatRegions(flat_region_table_t &t);
00359 
00363   ImageRegionType m_LargestPossibleRegion;
00364 
00365   bool m_SortEdgeLists;
00366   bool m_DoBoundaryAnalysis;
00367   double m_Threshold;
00368   double m_MaximumFloodLevel;
00369   InputPixelType m_Minimum;
00370   InputPixelType m_Maximum;
00371   unsigned long m_CurrentLabel;
00372 };
00373   
00374 }// end namespace watershed
00375 }// end namespace itk
00376 
00377 #ifndef ITK_MANUAL_INSTANTIATION
00378 #include "itkWatershedSegmenter.txx"
00379 #endif
00380 
00381 #endif

Generated at Fri May 21 01:15:42 2004 for ITK by doxygen 1.2.15 written by Dimitri van Heesch, © 1997-2000