ITK  5.4.0
Insight Toolkit
itkVoronoiSegmentationImageFilterBase.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright NumFOCUS
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  * https://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 : public ImageToImageFilter<TInputImage, TOutputImage>
60 {
61 public:
62  ITK_DISALLOW_COPY_AND_MOVE(VoronoiSegmentationImageFilterBase);
63 
69 
71  itkNewMacro(Self);
72 
74  itkOverrideGetNameOfClassMacro(VoronoiSegmentationImageFilterBase);
75 
77  static constexpr unsigned int ImageDimension = TInputImage::ImageDimension;
78 
80  using InputImageType = TInputImage;
84  using SizeType = typename TInputImage::SizeType;
86  using PixelType = typename TInputImage::PixelType;
87 
88  using OutputImageType = TOutputImage;
89  using OutputPixelType = typename TOutputImage::PixelType;
90 
97  using PointIdIterator = typename CellType::PointIdIterator;
103  using PointTypeVector = std::vector<PointType>;
104  using PointTypeDeque = std::deque<PointType>;
105  using BinaryObjectImage = TBinaryPriorImage;
107  using IndexList = std::vector<IndexType>;
108 
112 
114  itkSetMacro(NumberOfSeeds, int);
115  itkGetConstMacro(NumberOfSeeds, int);
119  itkSetMacro(MinRegion, SizeValueType);
120  itkGetConstMacro(MinRegion, SizeValueType);
125  itkSetMacro(Steps, int);
126  itkGetConstMacro(Steps, int);
130  itkGetConstMacro(LastStepSeeds, int);
131 
133  itkGetConstMacro(NumberOfSeedsToAdded, int);
134 
136  itkSetMacro(UseBackgroundInAPrior, bool);
137  itkGetConstMacro(UseBackgroundInAPrior, bool);
141  itkSetMacro(OutputBoundary, bool);
142  itkGetConstMacro(OutputBoundary, bool);
147  itkSetMacro(InteractiveSegmentation, bool);
148  itkGetConstMacro(InteractiveSegmentation, bool);
149  itkBooleanMacro(InteractiveSegmentation);
153  itkSetMacro(MeanDeviation, double);
154  itkGetConstMacro(MeanDeviation, double);
158  itkSetMacro(Size, SizeType);
159  itkGetConstMacro(Size, SizeType);
164  virtual void
166  {}
167 
169  void
170  RunSegment();
171 
173  void
174  RunSegmentOneStep();
175 
177  virtual void
178  MakeSegmentBoundary();
179 
180  virtual void
181  MakeSegmentObject();
182 
184  VoronoiPointer
186  {
187  return m_WorkingVD;
188  }
189 
190 #if !defined(ITK_WRAPPING_PARSER) // generates invalid iterator instantiation
191  // with msvc
195  void
196  SetSeeds(int num, SeedsIterator begin)
197  {
198  m_NumberOfSeeds = num;
199  m_WorkingVD->SetSeeds(num, begin);
200  }
203 #endif
204 
208  void
210  {
211  m_NumberOfSeeds = seeds.size();
212  auto it = seeds.begin();
213  m_WorkingVD->SetSeeds(m_NumberOfSeeds, it);
214  }
218  PointType
219  GetSeed(int SeedID)
220  {
221  return m_WorkingVD->GetSeed(SeedID);
222  }
223 
225  void
226  DrawDiagram(VDImagePointer result, unsigned char incolor, unsigned char outcolor, unsigned char boundcolor);
227 
228  void
229  BeforeNextStep();
230 
233  void
234  GenerateInputRequestedRegion() override;
235 
238  void
239  EnlargeOutputRequestedRegion(DataObject * output) override;
240 
241 protected:
243  ~VoronoiSegmentationImageFilterBase() override = default;
244  void
245  PrintSelf(std::ostream & os, Indent indent) const override;
246 
247  void
248  GenerateData() override; // general pipeline function.
249 
250  SizeType m_Size{};
251  int m_NumberOfSeeds{ 200 };
252  SizeValueType m_MinRegion{ 20 };
253  int m_Steps{ 0 };
254  int m_LastStepSeeds{ 0 };
255  int m_NumberOfSeedsToAdded{ 0 };
256  int m_NumberOfBoundary{ 0 };
257 
258  std::vector<SizeValueType> m_NumberOfPixels{};
259  std::vector<unsigned char> m_Label{};
260 
261  double m_MeanDeviation{ 0.8 };
262  bool m_UseBackgroundInAPrior{ false };
263  bool m_OutputBoundary{ false }; // if =1 then output the boundaries, if = 0 then
264  // output the object.
265  bool m_InteractiveSegmentation{ false };
266 
267  typename VoronoiDiagram::Pointer m_WorkingVD{};
268 
269  typename VoronoiDiagramGenerator::Pointer m_VDGenerator{};
270 
271  std::vector<PointType> m_SeedsToAdded{};
272 
273  // private methods:
274  // Classify all the voronoi cells as interior , exterior or boundary.
275  virtual void
276  ClassifyDiagram();
277 
278  // Generate the seeds to be added by dividing the boundary cells.
279  virtual void
280  GenerateAddingSeeds();
281 
282  // Compute the statistics of the pixels inside the cell.
283  void
284  GetPixelIndexFromPolygon(PointTypeDeque vertlist, IndexList * PixelPool);
285 
286  virtual bool
288  {
289  return true;
290  }
291 
292  void
293  FillPolygon(PointTypeDeque vertlist, OutputPixelType color = 1);
294 
295  // Draw a straight line to the output image.
296  void
297  drawLine(PointType p1, PointType p2);
298 
299  // Draw the intermediate Voronoi Diagram structure.
300  void
301  drawVDline(VDImagePointer result, PointType p1, PointType p2, unsigned char color);
302 };
303 } // namespace itk
304 
305 #ifndef ITK_MANUAL_INSTANTIATION
306 # include "itkVoronoiSegmentationImageFilterBase.hxx"
307 #endif
308 
309 #endif
Pointer
SmartPointer< Self > Pointer
Definition: itkAddImageFilter.h:93
itk::VoronoiSegmentationImageFilterBase::SetSeeds
void SetSeeds(int num, SeedsIterator begin)
Definition: itkVoronoiSegmentationImageFilterBase.h:196
ConstPointer
SmartPointer< const Self > ConstPointer
Definition: itkAddImageFilter.h:94
itk::VoronoiSegmentationImageFilterBase< TInputImage, TOutputImage >::VDImagePointer
typename VDImage::Pointer VDImagePointer
Definition: itkVoronoiSegmentationImageFilterBase.h:111
itk::VoronoiDiagram2D
Implements the 2-Dimensional Voronoi Diagram.
Definition: itkVoronoiDiagram2D.h:51
itk::VoronoiSegmentationImageFilterBase::GetVoronoiDiagram
VoronoiPointer GetVoronoiDiagram()
Definition: itkVoronoiSegmentationImageFilterBase.h:185
itk::Size
Represent a n-dimensional size (bounds) of a n-dimensional image.
Definition: itkSize.h:71
itk::GTest::TypedefsAndConstructors::Dimension2::PointType
ImageBaseType::PointType PointType
Definition: itkGTestTypedefsAndConstructors.h:51
itk::GTest::TypedefsAndConstructors::Dimension2::SizeType
ImageBaseType::SizeType SizeType
Definition: itkGTestTypedefsAndConstructors.h:49
itkImage.h
itk::VoronoiSegmentationImageFilterBase< TInputImage, TOutputImage >::SeedsType
typename VoronoiDiagram::SeedsType SeedsType
Definition: itkVoronoiSegmentationImageFilterBase.h:98
itk::SmartPointer< Self >
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::VoronoiSegmentationImageFilterBase< TInputImage, TOutputImage >::EdgeInfo
typename VoronoiDiagram::VoronoiEdge EdgeInfo
Definition: itkVoronoiSegmentationImageFilterBase.h:102
itk::VoronoiSegmentationImageFilterBase< TInputImage, TOutputImage >::EdgeIterator
typename VoronoiDiagram::VoronoiEdgeIterator EdgeIterator
Definition: itkVoronoiSegmentationImageFilterBase.h:101
itk::VoronoiSegmentationImageFilterBase< TInputImage, TOutputImage >::PointIdIterator
typename CellType::PointIdIterator PointIdIterator
Definition: itkVoronoiSegmentationImageFilterBase.h:97
itk::VoronoiDiagram2D::SeedsType
std::vector< PointType > SeedsType
Definition: itkVoronoiDiagram2D.h:115
itk::VoronoiSegmentationImageFilterBase< TInputImage, TOutputImage >::PointTypeDeque
std::deque< PointType > PointTypeDeque
Definition: itkVoronoiSegmentationImageFilterBase.h:104
itk::VoronoiSegmentationImageFilterBase< TInputImage, TOutputImage >::CellAutoPointer
typename VoronoiDiagram::CellAutoPointer CellAutoPointer
Definition: itkVoronoiSegmentationImageFilterBase.h:95
itk::VoronoiSegmentationImageFilterBase< TInputImage, TOutputImage >::VoronoiPointer
typename VoronoiDiagram::Pointer VoronoiPointer
Definition: itkVoronoiSegmentationImageFilterBase.h:96
itk::VoronoiSegmentationImageFilterBase< TInputImage, TOutputImage >::SeedsIterator
typename VoronoiDiagram::SeedsIterator SeedsIterator
Definition: itkVoronoiSegmentationImageFilterBase.h:99
itk::VoronoiSegmentationImageFilterBase< TInputImage, TOutputImage >::CellType
typename VoronoiDiagram::CellType CellType
Definition: itkVoronoiSegmentationImageFilterBase.h:94
itk::Mesh< TCoordType, 2, DefaultDynamicMeshTraits< TCoordType, 2, 2, TCoordType > >::CellAutoPointer
typename CellType::CellAutoPointer CellAutoPointer
Definition: itkMesh.h:229
itk::VoronoiSegmentationImageFilterBase< TInputImage, TOutputImage >::IndexList
std::vector< IndexType > IndexList
Definition: itkVoronoiSegmentationImageFilterBase.h:107
itk::GTest::TypedefsAndConstructors::Dimension2::IndexType
ImageBaseType::IndexType IndexType
Definition: itkGTestTypedefsAndConstructors.h:50
itk::ImageToImageFilter
Base class for filters that take an image as input and produce an image as output.
Definition: itkImageToImageFilter.h:108
itk::ImageSource
Base class for all process objects that output image data.
Definition: itkImageSource.h:67
itk::VoronoiSegmentationImageFilterBase::TakeAPrior
virtual void TakeAPrior(const BinaryObjectImage *)
Definition: itkVoronoiSegmentationImageFilterBase.h:165
itk::VoronoiSegmentationImageFilterBase
Base class for VoronoiSegmentationImageFilter.
Definition: itkVoronoiSegmentationImageFilterBase.h:59
itk::VoronoiSegmentationImageFilterBase< TInputImage, TOutputImage >::PointType
typename VoronoiDiagram::PointType PointType
Definition: itkVoronoiSegmentationImageFilterBase.h:93
itk::ImageToImageFilter::InputImagePointer
typename InputImageType::Pointer InputImagePointer
Definition: itkImageToImageFilter.h:130
itk::VoronoiSegmentationImageFilterBase< TInputImage, TOutputImage >::BinaryObjectImagePointer
typename BinaryObjectImage::Pointer BinaryObjectImagePointer
Definition: itkVoronoiSegmentationImageFilterBase.h:106
itk::GTest::TypedefsAndConstructors::Dimension2::RegionType
ImageBaseType::RegionType RegionType
Definition: itkGTestTypedefsAndConstructors.h:54
itk::VoronoiDiagram2D::VoronoiEdgeIterator
typename std::vector< VoronoiEdge >::iterator VoronoiEdgeIterator
Definition: itkVoronoiDiagram2D.h:177
itk::VoronoiSegmentationImageFilterBase::GetSeed
PointType GetSeed(int SeedID)
Definition: itkVoronoiSegmentationImageFilterBase.h:219
itk::ImageToImageFilter::InputImageType
TInputImage InputImageType
Definition: itkImageToImageFilter.h:129
itkImageToImageFilter.h
itk::VoronoiSegmentationImageFilterBase< TInputImage, TOutputImage >::NeighborIdIterator
typename VoronoiDiagram::NeighborIdIterator NeighborIdIterator
Definition: itkVoronoiSegmentationImageFilterBase.h:100
itk::VoronoiSegmentationImageFilterBase< TInputImage, TOutputImage >::PixelType
typename TInputImage::PixelType PixelType
Definition: itkVoronoiSegmentationImageFilterBase.h:86
itk::CellInterface
An abstract interface for cells.
Definition: itkCellInterface.h:96
itk::VoronoiSegmentationImageFilterBase::TestHomogeneity
virtual bool TestHomogeneity(IndexList &)
Definition: itkVoronoiSegmentationImageFilterBase.h:287
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnnulusOperator.h:24
itk::ProcessObject
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
Definition: itkProcessObject.h:139
itk::VoronoiSegmentationImageFilterBase< TInputImage, TOutputImage >::PointTypeVector
std::vector< PointType > PointTypeVector
Definition: itkVoronoiSegmentationImageFilterBase.h:103
itkVoronoiDiagram2DGenerator.h
itk::VoronoiSegmentationImageFilterBase< TInputImage, TOutputImage >::IndexType
typename TInputImage::IndexType IndexType
Definition: itkVoronoiSegmentationImageFilterBase.h:83
itk::VoronoiSegmentationImageFilterBase< TInputImage, TOutputImage >::OutputPixelType
typename TOutputImage::PixelType OutputPixelType
Definition: itkVoronoiSegmentationImageFilterBase.h:89
itk::VoronoiSegmentationImageFilterBase::SetSeeds
void SetSeeds(SeedsType &seeds)
Definition: itkVoronoiSegmentationImageFilterBase.h:209
itk::VoronoiDiagram2DGenerator
Implement the Sweep Line Algorithm for the construction of the 2D Voronoi Diagram.
Definition: itkVoronoiDiagram2DGenerator.h:50
itk::Image
Templated n-dimensional image class.
Definition: itkImage.h:88
itk::VoronoiSegmentationImageFilterBase< TInputImage, TOutputImage >::RegionType
typename TInputImage::RegionType RegionType
Definition: itkVoronoiSegmentationImageFilterBase.h:85
itk::VoronoiDiagram2D::NeighborIdIterator
typename INTvector::iterator NeighborIdIterator
Definition: itkVoronoiDiagram2D.h:121
itk::VoronoiSegmentationImageFilterBase< TInputImage, TOutputImage >::SizeType
typename TInputImage::SizeType SizeType
Definition: itkVoronoiSegmentationImageFilterBase.h:84
itk::VoronoiDiagram2D::SeedsIterator
typename SeedsType::iterator SeedsIterator
Definition: itkVoronoiDiagram2D.h:116
itk::ImageToImageFilter::InputImageConstPointer
typename InputImageType::ConstPointer InputImageConstPointer
Definition: itkImageToImageFilter.h:131
itk::SizeValueType
unsigned long SizeValueType
Definition: itkIntTypes.h:83
itk::ImageSource::OutputImageType
TOutputImage OutputImageType
Definition: itkImageSource.h:90
itk::DataObject
Base class for all data objects in ITK.
Definition: itkDataObject.h:293
itk::VoronoiDiagram2D::VoronoiEdge
Definition: itkVoronoiDiagram2D.h:164