ITK  4.9.0
Insight Segmentation and Registration Toolkit
itkRegionCompetitionImageFilter.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Insight Segmentation & Registration Toolkit
4  Module: itkRegionCompetitionImageFilter.h
5  Language: C++
6  Date: $Date$
7  Version: $Revision$
8 
9  Copyright (c) Insight Software Consortium. All rights reserved.
10  See ITKCopyright.txt or http://www.itk.org/HTML/Copyright.htm for details.
11 
12  This software is distributed WITHOUT ANY WARRANTY; without even
13  the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
14  PURPOSE. See the above copyright notices for more information.
15 
16 =========================================================================*/
17 #ifndef itkRegionCompetitionImageFilter_h
18 #define itkRegionCompetitionImageFilter_h
19 
20 #include "itkImage.h"
21 #include "itkImageToImageFilter.h"
22 
23 #include <vector>
24 
25 namespace itk
26 {
27 
40 template <class TInputImage, class TOutputImage>
42  public ImageToImageFilter<TInputImage,TOutputImage>
43 {
44 public:
50 
52  itkNewMacro(Self);
53 
56 
57  typedef TInputImage InputImageType;
58  typedef typename InputImageType::Pointer InputImagePointer;
59  typedef typename InputImageType::ConstPointer InputImageConstPointer;
60  typedef typename InputImageType::RegionType InputImageRegionType;
61  typedef typename InputImageType::SizeType InputSizeType;
62  typedef typename InputImageType::PixelType InputImagePixelType;
63  typedef typename InputImageType::IndexType IndexType;
65 
66  typedef TOutputImage OutputImageType;
67  typedef typename OutputImageType::Pointer OutputImagePointer;
68  typedef typename OutputImageType::RegionType OutputImageRegionType;
69  typedef typename OutputImageType::PixelType OutputImagePixelType;
70 
71 
73  itkStaticConstMacro(InputImageDimension, unsigned int, TInputImage::ImageDimension);
74  itkStaticConstMacro(OutputImageDimension, unsigned int, TOutputImage::ImageDimension);
76 
79  itkSetMacro( MaximumNumberOfIterations, unsigned int );
80  itkGetMacro( MaximumNumberOfIterations, unsigned int );
82 
84  itkGetMacro( CurrentIterationNumber, unsigned int );
85 
87  itkGetMacro( TotalNumberOfPixelsChanged, unsigned int );
88 
90  void SetInputLabels( const TOutputImage * inputLabelImage );
91 
92 #ifdef ITK_USE_CONCEPT_CHECKING
93 
100 
102 #endif
103 
104 protected:
107 
108  void GenerateData();
109 
110  void PrintSelf ( std::ostream& os, Indent indent ) const;
111 
112 private:
113  RegionCompetitionImageFilter(const Self&); //purposely not implemented
114  void operator=(const Self&); //purposely not implemented
115 
116 
118 
120 
122 
123  void InitializeNeighborhood();
124 
126 
128 
130 
132 
133  void SwapSeedArrays();
134 
135  void ClearSecondSeedArray();
136 
138 
140 
142 
143  void ComputeBirthThreshold();
144 
145  itkSetMacro( CurrentPixelIndex, IndexType );
146  itkGetConstReferenceMacro( CurrentPixelIndex, IndexType );
147 
148  typedef std::vector<IndexType> SeedArrayType;
149 
152 
154 
155  typedef std::vector<OutputImagePixelType> SeedNewValuesArrayType;
156 
158 
163 
165 
166  //
167  // Variables used for addressing the Neighbors.
168  // This could be factorized into a helper class.
169  //
171 
172  typedef std::vector< OffsetValueType > NeighborOffsetArrayType;
173 
175 
176 
177  //
178  // Helper cache variables
179  //
183 
186 
188 
190 
192 
193  mutable unsigned int m_NumberOfLabels;
194 };
195 
196 } // end namespace itk
197 
198 #ifndef ITK_MANUAL_INSTANTIATION
199 #include "itkRegionCompetitionImageFilter.hxx"
200 #endif
201 
202 #endif
itk::Image< unsigned char, InputImageDimension > SeedMaskImageType
itk::Neighborhood< InputImagePixelType, InputImageDimension > NeighborhoodType
OffsetValueType m_OffsetTable[InputImageDimension+1]
void SetInputLabels(const TOutputImage *inputLabelImage)
signed long OffsetValueType
Definition: itkIntTypes.h:154
std::vector< OutputImagePixelType > SeedNewValuesArrayType
InputImageType::ConstPointer InputImageConstPointer
Base class for all process objects that output image data.
InputImageType::OffsetValueType OffsetValueType
InputImageType::RegionType InputImageRegionType
void PrintSelf(std::ostream &os, Indent indent) const
std::vector< OffsetValueType > NeighborOffsetArrayType
Base class for filters that take an image as input and produce an image as output.
Control indentation during Print() invocation.
Definition: itkIndent.h:49
Perform front-propagation from different starting labeled regions.
#define itkConceptMacro(name, concept)
Templated n-dimensional image class.
Definition: itkImage.h:75
bool TestForAvailabilityAtCurrentPixel() const
ImageToImageFilter< TInputImage, TOutputImage > Superclass