ITK  4.8.0
Insight Segmentation and Registration Toolkit
itkVoronoiSegmentationImageFilterBase.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright Insight Software Consortium
4  *
5  * Licensed under the Apache License, Version 2.0 (the "License");
6  * you may not use this file except in compliance with the License.
7  * You may obtain a copy of the License at
8  *
9  * http://www.apache.org/licenses/LICENSE-2.0.txt
10  *
11  * Unless required by applicable law or agreed to in writing, software
12  * distributed under the License is distributed on an "AS IS" BASIS,
13  * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
14  * See the License for the specific language governing permissions and
15  * limitations under the License.
16  *
17  *=========================================================================*/
18 #ifndef itkVoronoiSegmentationImageFilterBase_h
19 #define itkVoronoiSegmentationImageFilterBase_h
20 
21 #include "itkImageToImageFilter.h"
23 #include "itkImage.h"
24 
25 namespace itk
26 {
58 template< typename TInputImage, typename TOutputImage, typename TBinaryPriorImage = Image< unsigned char, 2 > >
60  public ImageToImageFilter< TInputImage, TOutputImage >
61 {
62 public:
68 
70  itkNewMacro(Self);
71 
74 
76  itkStaticConstMacro(ImageDimension, unsigned int,
77  TInputImage::ImageDimension);
78 
80  typedef TInputImage InputImageType;
81  typedef typename TInputImage::Pointer InputImagePointer;
82  typedef typename TInputImage::ConstPointer InputImageConstPointer;
83  typedef typename TInputImage::IndexType IndexType;
84  typedef typename TInputImage::SizeType SizeType;
85  typedef typename TInputImage::RegionType RegionType;
86  typedef typename TInputImage::PixelType PixelType;
87 
88  typedef TOutputImage OutputImageType;
89  typedef typename TOutputImage::PixelType OutputPixelType;
90 
97  typedef typename CellType::PointIdIterator PointIdIterator;
103  typedef std::vector< PointType > PointTypeVector;
104  typedef std::deque< PointType > PointTypeDeque;
105  typedef TBinaryPriorImage BinaryObjectImage;
107  typedef std::vector< IndexType > IndexList;
108 
112 
114  itkSetMacro(NumberOfSeeds, int);
115  itkGetConstMacro(NumberOfSeeds, int);
117 
119  itkSetMacro(MinRegion, int);
120  itkGetConstMacro(MinRegion, int);
122 
125  itkSetMacro(Steps, int);
126  itkGetConstMacro(Steps, int);
128 
130  itkGetConstMacro(LastStepSeeds, int);
131 
133  itkGetConstMacro(NumberOfSeedsToAdded, int);
134 
136  itkSetMacro(UseBackgroundInAPrior, bool);
137  itkGetConstMacro(UseBackgroundInAPrior, bool);
139 
141  itkSetMacro(OutputBoundary, bool);
142  itkGetConstMacro(OutputBoundary, bool);
144 
147  itkSetMacro(InteractiveSegmentation, bool);
148  itkGetConstMacro(InteractiveSegmentation, bool);
149  itkBooleanMacro(InteractiveSegmentation);
151 
153  itkSetMacro(MeanDeviation, double);
154  itkGetConstMacro(MeanDeviation, double);
156 
158  itkSetMacro(Size, SizeType);
159  itkGetConstMacro(Size, SizeType);
161 
164  virtual void TakeAPrior(const BinaryObjectImage *){}
165 
167  void RunSegment();
168 
170  void RunSegmentOneStep();
171 
173  virtual void MakeSegmentBoundary();
174 
175  virtual void MakeSegmentObject();
176 
179  { return m_WorkingVD; }
180 
181 #if !defined( CABLE_CONFIGURATION ) // generates invalid iterator instantiation
182  // with msvc
186  void SetSeeds(int num, SeedsIterator begin)
187  {
188  m_NumberOfSeeds = num;
189  m_WorkingVD->SetSeeds(num, begin);
190  }
192 
193 #endif
194 
198  void SetSeeds(SeedsType & seeds)
199  {
200  m_NumberOfSeeds = seeds.size();
201  typename SeedsType::iterator it = seeds.begin();
202  m_WorkingVD->SetSeeds(m_NumberOfSeeds, it);
203  }
205 
207  PointType GetSeed(int SeedID)
208  { return m_WorkingVD->GetSeed(SeedID); }
209 
211  void DrawDiagram(VDImagePointer result, unsigned char incolor,
212  unsigned char outcolor, unsigned char boundcolor);
213 
214  void BeforeNextStep();
215 
218  virtual void GenerateInputRequestedRegion() ITK_OVERRIDE;
219 
222  virtual void EnlargeOutputRequestedRegion(DataObject *output) ITK_OVERRIDE;
223 
224 protected:
226  ~VoronoiSegmentationImageFilterBase();
227  virtual void PrintSelf(std::ostream & os, Indent indent) const ITK_OVERRIDE;
228 
229  void GenerateData() ITK_OVERRIDE; //general pipeline function.
230 
234  int m_Steps;
238 
239  std::vector< int > m_NumberOfPixels;
240  std::vector< unsigned char > m_Label;
241 
244  bool m_OutputBoundary; //if =1 then output the boundaries, if = 0 then
245  // output the object.
247 
249 
251 
252  std::vector< PointType > m_SeedsToAdded;
253 
254  // private methods:
255  // Classify all the voronoi cells as interior , exterior or boundary.
256  virtual void ClassifyDiagram();
257 
258  // Generate the seeds to be added by dividing the boundary cells.
259  virtual void GenerateAddingSeeds();
260 
261  // Compute the statistics of the pixels inside the cell.
262  void GetPixelIndexFromPolygon(PointTypeDeque VertList, IndexList *PixelPool);
263 
264  virtual bool TestHomogeneity(IndexList &)
265  { return 1; }
266 
267  void FillPolygon(PointTypeDeque vertlist, OutputPixelType color = 1);
268 
269  // Draw a straight line to the output image.
270  void drawLine(PointType p1, PointType p2);
271 
272  // Draw the intermedia Voronoi Diagram structure.
273  void drawVDline(VDImagePointer result, PointType p1, PointType p2, unsigned char color);
274 
275 private:
276  VoronoiSegmentationImageFilterBase(const Self &); //purposely not implemented
277  void operator=(const Self &); //purposely not implemented
278 };
279 } //end namespace
280 
281 #ifndef ITK_MANUAL_INSTANTIATION
282 #include "itkVoronoiSegmentationImageFilterBase.hxx"
283 #endif
284 
285 #endif
Superclass::CellType CellType
void drawVDline(VDImagePointer result, PointType p1, PointType p2, unsigned char color)
Base class for VoronoiSegmentationImageFilter.
void GetPixelIndexFromPolygon(PointTypeDeque VertList, IndexList *PixelPool)
std::vector< VoronoiEdge >::iterator VoronoiEdgeIterator
Represent the size (bounds) of a n-dimensional image.
Definition: itkSize.h:52
void FillPolygon(PointTypeDeque vertlist, OutputPixelType color=1)
SeedsType::iterator SeedsIterator
virtual void PrintSelf(std::ostream &os, Indent indent) const override
VoronoiDiagram2DGenerator< double > VoronoiDiagramGenerator
virtual void GenerateInputRequestedRegion() override
Implements the 2-Dimensional Voronoi Diagram.
ImageToImageFilter< TInputImage, TOutputImage > Superclass
Base class for all process objects that output image data.
virtual void EnlargeOutputRequestedRegion(DataObject *output) override
INTvector::iterator NeighborIdIterator
Base class for filters that take an image as input and produce an image as output.
Implement the Sweep Line Algorithm for the construction of the 2D Voronoi Diagram.
Control indentation during Print() invocation.
Definition: itkIndent.h:49
std::vector< PointType > SeedsType
Superclass::CellAutoPointer CellAutoPointer
void DrawDiagram(VDImagePointer result, unsigned char incolor, unsigned char outcolor, unsigned char boundcolor)
A templated class holding a geometric point in n-Dimensional space.
Definition: itkPoint.h:51
void drawLine(PointType p1, PointType p2)
Base class for all data objects in ITK.
Templated n-dimensional image class.
Definition: itkImage.h:75