ITK  4.8.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 {
149 template< typename TInputImage >
151  public ImageToImageFilter< TInputImage, Image< IdentifierType,
152  TInputImage::ImageDimension > >
153 {
154 public:
157 
159  typedef TInputImage InputImageType;
160 
162  itkStaticConstMacro (ImageDimension, unsigned int,
163  TInputImage::ImageDimension);
164 
167 
169  typedef typename InputImageType::RegionType RegionType;
170  typedef typename InputImageType::SizeType SizeType;
171  typedef typename InputImageType::IndexType IndexType;
172 
175 
177  typedef typename InputImageType::PixelType ScalarType;
178 
181 
184 
186  itkNewMacro(Self);
187 
189  void GenerateData() ITK_OVERRIDE;
190 
193  using Superclass::SetInput;
194  void SetInput(const InputImageType *input) ITK_OVERRIDE
195  {
196  // if the input is changed, we'll need to clear the cached tree
197  // when we execute
198  if ( input != this->GetInput(0) )
199  {
200  m_InputChanged = true;
201  }
203 
204  // processObject is not const-correct so a const_cast is needed here
205  this->ProcessObject::SetNthInput( 0, const_cast< InputImageType * >( input ) );
206  m_Segmenter->SetInputImage( const_cast< InputImageType * >( input ) );
207  }
208 
209  virtual void SetInput(unsigned int i, const TInputImage *image) ITK_OVERRIDE
210  {
211  if ( i != 0 )
212  { itkExceptionMacro(<< "Filter has only one input."); }
213  else
214  { this->SetInput(image); }
215  }
216 
219  void SetThreshold(double);
220 
221  itkGetConstMacro(Threshold, double);
222 
225  void SetLevel(double);
226 
227  itkGetConstMacro(Level, double);
228 
232  {
233  m_Segmenter->Update();
234  return m_Segmenter->GetOutputImage();
235  }
237 
241  {
243  }
244 
245  // Override since the filter produces all of its output
246  void EnlargeOutputRequestedRegion(DataObject *data) ITK_OVERRIDE;
247 
248 #ifdef ITK_USE_CONCEPT_CHECKING
249  // Begin concept checking
250  itkConceptMacro( InputEqualityComparableCheck,
252  itkConceptMacro( InputAdditiveOperatorsCheck,
254  itkConceptMacro( DoubleInputMultiplyOperatorCheck,
256  itkConceptMacro( InputLessThanComparableCheck,
258  // End concept checking
259 #endif
260 
261 protected:
264  void PrintSelf(std::ostream & os, Indent indent) const ITK_OVERRIDE;
265 
268  virtual void PrepareOutputs() ITK_OVERRIDE;
269 
270 private:
271  WatershedImageFilter(const Self &); //purposely not implemented
272  void operator=(const Self &); //purposely not implemented
273 
277  double m_Threshold;
278 
283  double m_Level;
284 
289  typename watershed::Segmenter< InputImageType >::Pointer m_Segmenter;
290 
291  typename watershed::SegmentTreeGenerator< ScalarType >::Pointer m_TreeGenerator;
292 
293  typename watershed::Relabeler< ScalarType, itkGetStaticConstMacro(ImageDimension) >::Pointer m_Relabeler;
294 
295  unsigned long m_ObserverTag;
296 
300 
302 };
303 } // end namespace itk
304 
305 #ifndef ITK_MANUAL_INSTANTIATION
306 #include "itkWatershedImageFilter.hxx"
307 #endif
308 
309 #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.
OutputImageType * GetOutputImage(void)
watershed::Relabeler< ScalarType, itkGetStaticConstMacro(ImageDimension) >::Pointer m_Relabeler
Light weight base class for most itk classes.
void EnlargeOutputRequestedRegion(DataObject *data) override
void PrintSelf(std::ostream &os, Indent indent) const override
ImageToImageFilter< InputImageType, OutputImageType > Superclass
static const unsigned int ImageDimension
void SetInput(const InputImageType *input) override
InputImageType::SizeType SizeType
watershed::SegmentTreeGenerator< ScalarType >::SegmentTreeType * GetSegmentTree()
watershed::SegmentTreeGenerator< ScalarType >::Pointer m_TreeGenerator
InputImageType::PixelType ScalarType
InputImageType::RegionType RegionType
Generate a unique, increasing time value.
Definition: itkTimeStamp.h:58
virtual void PrepareOutputs() override
virtual void SetInput(unsigned int i, const TInputImage *image) override
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
InputImageType::IndexType IndexType
watershed::Segmenter< InputImageType >::Pointer m_Segmenter
void GenerateData() override
virtual void SetNthInput(DataObjectPointerArraySizeType num, DataObject *input)
watershed::Segmenter< InputImageType >::OutputImageType * GetBasicSegmentation()
#define itkConceptMacro(name, concept)
Image< IdentifierType, itkGetStaticConstMacro(ImageDimension) > OutputImageType
Base class for all data objects in ITK.
Templated n-dimensional image class.
Definition: itkImage.h:75
void SetInputImage(InputImageType *img)