ITK  5.0.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 > >
59 class ITK_TEMPLATE_EXPORT VoronoiSegmentationImageFilterBase:
60  public ImageToImageFilter< TInputImage, TOutputImage >
61 {
62 public:
63  ITK_DISALLOW_COPY_AND_ASSIGN(VoronoiSegmentationImageFilterBase);
64 
70 
72  itkNewMacro(Self);
73 
76 
78  static constexpr unsigned int ImageDimension = TInputImage::ImageDimension;
79 
81  using InputImageType = TInputImage;
82  using InputImagePointer = typename TInputImage::Pointer;
83  using InputImageConstPointer = typename TInputImage::ConstPointer;
85  using SizeType = typename TInputImage::SizeType;
87  using PixelType = typename TInputImage::PixelType;
88 
89  using OutputImageType = TOutputImage;
90  using OutputPixelType = typename TOutputImage::PixelType;
91 
98  using PointIdIterator = typename CellType::PointIdIterator;
104  using PointTypeVector = std::vector< PointType >;
105  using PointTypeDeque = std::deque< PointType >;
106  using BinaryObjectImage = TBinaryPriorImage;
108  using IndexList = std::vector< IndexType >;
109 
113 
115  itkSetMacro(NumberOfSeeds, int);
116  itkGetConstMacro(NumberOfSeeds, int);
118 
120  itkSetMacro(MinRegion, SizeValueType);
121  itkGetConstMacro(MinRegion, SizeValueType);
123 
126  itkSetMacro(Steps, int);
127  itkGetConstMacro(Steps, int);
129 
131  itkGetConstMacro(LastStepSeeds, int);
132 
134  itkGetConstMacro(NumberOfSeedsToAdded, int);
135 
137  itkSetMacro(UseBackgroundInAPrior, bool);
138  itkGetConstMacro(UseBackgroundInAPrior, bool);
140 
142  itkSetMacro(OutputBoundary, bool);
143  itkGetConstMacro(OutputBoundary, bool);
145 
148  itkSetMacro(InteractiveSegmentation, bool);
149  itkGetConstMacro(InteractiveSegmentation, bool);
150  itkBooleanMacro(InteractiveSegmentation);
152 
154  itkSetMacro(MeanDeviation, double);
155  itkGetConstMacro(MeanDeviation, double);
157 
159  itkSetMacro(Size, SizeType);
160  itkGetConstMacro(Size, SizeType);
162 
165  virtual void TakeAPrior(const BinaryObjectImage *){}
166 
168  void RunSegment();
169 
171  void RunSegmentOneStep();
172 
174  virtual void MakeSegmentBoundary();
175 
176  virtual void MakeSegmentObject();
177 
180  { return m_WorkingVD; }
181 
182 #if !defined( ITK_WRAPPING_PARSER ) // generates invalid iterator instantiation
183  // with msvc
187  void SetSeeds(int num, SeedsIterator begin)
188  {
189  m_NumberOfSeeds = num;
190  m_WorkingVD->SetSeeds(num, begin);
191  }
193 
194 #endif
195 
199  void SetSeeds(SeedsType & seeds)
200  {
201  m_NumberOfSeeds = seeds.size();
202  auto it = seeds.begin();
203  m_WorkingVD->SetSeeds(m_NumberOfSeeds, it);
204  }
206 
208  PointType GetSeed(int SeedID)
209  { return m_WorkingVD->GetSeed(SeedID); }
210 
212  void DrawDiagram(VDImagePointer result, unsigned char incolor,
213  unsigned char outcolor, unsigned char boundcolor);
214 
215  void BeforeNextStep();
216 
219  void GenerateInputRequestedRegion() override;
220 
223  void EnlargeOutputRequestedRegion(DataObject *output) override;
224 
225 protected:
227  ~VoronoiSegmentationImageFilterBase() override = default;
228  void PrintSelf(std::ostream & os, Indent indent) const override;
229 
230  void GenerateData() override; //general pipeline function.
231 
233  int m_NumberOfSeeds{200};
234  SizeValueType m_MinRegion{20};
235  int m_Steps{0};
236  int m_LastStepSeeds{0};
237  int m_NumberOfSeedsToAdded{0};
238  int m_NumberOfBoundary{0};
239 
240  std::vector< SizeValueType > m_NumberOfPixels;
241  std::vector< unsigned char > m_Label;
242 
243  double m_MeanDeviation{0.8};
244  bool m_UseBackgroundInAPrior{false};
245  bool m_OutputBoundary{false}; //if =1 then output the boundaries, if = 0 then
246  // output the object.
247  bool m_InteractiveSegmentation{false};
248 
250 
252 
253  std::vector< PointType > m_SeedsToAdded;
254 
255  // private methods:
256  // Classify all the voronoi cells as interior , exterior or boundary.
257  virtual void ClassifyDiagram();
258 
259  // Generate the seeds to be added by dividing the boundary cells.
260  virtual void GenerateAddingSeeds();
261 
262  // Compute the statistics of the pixels inside the cell.
263  void GetPixelIndexFromPolygon(PointTypeDeque VertList, IndexList *PixelPool);
264 
265  virtual bool TestHomogeneity(IndexList &)
266  { return true; }
267 
268  void FillPolygon(PointTypeDeque vertlist, OutputPixelType color = 1);
269 
270  // Draw a straight line to the output image.
271  void drawLine(PointType p1, PointType p2);
272 
273  // Draw the intermedia Voronoi Diagram structure.
274  void drawVDline(VDImagePointer result, PointType p1, PointType p2, unsigned char color);
275 };
276 } //end namespace
277 
278 #ifndef ITK_MANUAL_INSTANTIATION
279 #include "itkVoronoiSegmentationImageFilterBase.hxx"
280 #endif
281 
282 #endif
Base class for VoronoiSegmentationImageFilter.
unsigned long SizeValueType
Definition: itkIntTypes.h:83
typename Superclass::CellType CellType
std::vector< PointType > SeedsType
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
typename Superclass::CellAutoPointer CellAutoPointer
Implements the 2-Dimensional Voronoi Diagram.
Base class for all process objects that output image data.
typename InputImageType::Pointer InputImagePointer
TOutputImage OutputImageType
Represent a n-dimensional size (bounds) of a n-dimensional image.
Definition: itkSize.h:68
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
typename std::vector< VoronoiEdge >::iterator VoronoiEdgeIterator
typename INTvector::iterator NeighborIdIterator
typename InputImageType::ConstPointer InputImageConstPointer
Base class for all data objects in ITK.
typename SeedsType::iterator SeedsIterator
Templated n-dimensional image class.
Definition: itkImage.h:75