ITK  5.0.0
Insight Segmentation and Registration Toolkit
itkWatershedImageFilter.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 itkWatershedImageFilter_h
19 #define itkWatershedImageFilter_h
20 
21 
22 #include "itkImageToImageFilter.h"
24 #include "itkWatershedRelabeler.h"
26 
27 namespace itk
28 {
144 template< typename TInputImage >
145 class ITK_TEMPLATE_EXPORT WatershedImageFilter:
146  public ImageToImageFilter< TInputImage, Image< IdentifierType,
147  TInputImage::ImageDimension > >
148 {
149 public:
150  ITK_DISALLOW_COPY_AND_ASSIGN(WatershedImageFilter);
151 
154 
156  using InputImageType = TInputImage;
157 
159  static constexpr unsigned int ImageDimension = TInputImage::ImageDimension;
160 
163 
168 
171 
173  using ScalarType = typename InputImageType::PixelType;
174 
177 
180 
182  itkNewMacro(Self);
183 
185  void GenerateData() override;
186 
189  using Superclass::SetInput;
190  void SetInput(const InputImageType *input) override
191  {
192  // if the input is changed, we'll need to clear the cached tree
193  // when we execute
194  if ( input != this->GetInput(0) )
195  {
196  m_InputChanged = true;
197  }
199 
200  // processObject is not const-correct so a const_cast is needed here
201  this->ProcessObject::SetNthInput( 0, const_cast< InputImageType * >( input ) );
202  m_Segmenter->SetInputImage( const_cast< InputImageType * >( input ) );
203  }
204 
205  void SetInput(unsigned int i, const TInputImage *image) override
206  {
207  if ( i != 0 )
208  { itkExceptionMacro(<< "Filter has only one input."); }
209  else
210  { this->SetInput(image); }
211  }
212 
215  void SetThreshold(double);
216 
217  itkGetConstMacro(Threshold, double);
218 
221  void SetLevel(double);
222 
223  itkGetConstMacro(Level, double);
224 
228  {
229  m_Segmenter->Update();
230  return m_Segmenter->GetOutputImage();
231  }
233 
237  {
238  return m_TreeGenerator->GetOutputSegmentTree();
239  }
240 
241  // Override since the filter produces all of its output
242  void EnlargeOutputRequestedRegion(DataObject *data) override;
243 
244 #ifdef ITK_USE_CONCEPT_CHECKING
245  // Begin concept checking
246  itkConceptMacro( InputEqualityComparableCheck,
248  itkConceptMacro( InputAdditiveOperatorsCheck,
250  itkConceptMacro( DoubleInputMultiplyOperatorCheck,
252  itkConceptMacro( InputLessThanComparableCheck,
254  // End concept checking
255 #endif
256 
257 protected:
259  ~WatershedImageFilter() override = default;
260  void PrintSelf(std::ostream & os, Indent indent) const override;
261 
264  void PrepareOutputs() override;
265 
266 private:
270  double m_Threshold{0.0};
271 
276  double m_Level{0.0};
277 
283 
285 
287 
288  unsigned long m_ObserverTag;
289 
293 
295 };
296 } // end namespace itk
297 
298 #ifndef ITK_MANUAL_INSTANTIATION
299 #include "itkWatershedImageFilter.hxx"
300 #endif
301 
302 #endif
virtual void Update()
Bring this filter up-to-date.
A low-level image analysis algorithm that automatically produces a hierarchy of segmented, labeled images from a scalar-valued image input.
watershed::Relabeler< ScalarType, Self::ImageDimension >::Pointer m_Relabeler
Light weight base class for most itk classes.
void SetInput(unsigned int i, const TInputImage *image) override
typename InputImageType::PixelType ScalarType
void SetInput(const InputImageType *input) override
watershed::SegmentTreeGenerator< ScalarType >::SegmentTreeType * GetSegmentTree()
watershed::SegmentTreeGenerator< ScalarType >::Pointer m_TreeGenerator
Generate a unique, increasing time value.
Definition: itkTimeStamp.h:60
typename InputImageType::RegionType RegionType
typename InputImageType::IndexType IndexType
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
watershed::Segmenter< InputImageType >::Pointer m_Segmenter
virtual void SetNthInput(DataObjectPointerArraySizeType num, DataObject *input)
watershed::Segmenter< InputImageType >::OutputImageType * GetBasicSegmentation()
#define itkConceptMacro(name, concept)
typename InputImageType::SizeType SizeType
Base class for all data objects in ITK.
Templated n-dimensional image class.
Definition: itkImage.h:75