ITK  5.0.0
Insight Segmentation and Registration Toolkit
itkMorphologicalWatershedImageFilter.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 itkMorphologicalWatershedImageFilter_h
19 #define itkMorphologicalWatershedImageFilter_h
20 
21 #include "itkImageToImageFilter.h"
22 
23 namespace itk
24 {
51 template< typename TInputImage, typename TOutputImage >
52 class ITK_TEMPLATE_EXPORT MorphologicalWatershedImageFilter:
53  public ImageToImageFilter< TInputImage, TOutputImage >
54 {
55 public:
56  ITK_DISALLOW_COPY_AND_ASSIGN(MorphologicalWatershedImageFilter);
57 
63 
65  using InputImageType = TInputImage;
66  using OutputImageType = TOutputImage;
67  using InputImagePointer = typename InputImageType::Pointer;
68  using InputImageConstPointer = typename InputImageType::ConstPointer;
70  using InputImagePixelType = typename InputImageType::PixelType;
71  using OutputImagePointer = typename OutputImageType::Pointer;
72  using OutputImageConstPointer = typename OutputImageType::ConstPointer;
74  using OutputImagePixelType = typename OutputImageType::PixelType;
75 
77  static constexpr unsigned int InputImageDimension = TInputImage::ImageDimension;
78  static constexpr unsigned int OutputImageDimension = TOutputImage::ImageDimension;
79 
81  itkNewMacro(Self);
82 
85 
92  itkSetMacro(FullyConnected, bool);
93  itkGetConstReferenceMacro(FullyConnected, bool);
94  itkBooleanMacro(FullyConnected);
96 
102  itkSetMacro(MarkWatershedLine, bool);
103  itkGetConstReferenceMacro(MarkWatershedLine, bool);
104  itkBooleanMacro(MarkWatershedLine);
106 
109  itkSetMacro(Level, InputImagePixelType);
110  itkGetConstMacro(Level, InputImagePixelType);
112 
113 protected:
115  ~MorphologicalWatershedImageFilter() override = default;
116  void PrintSelf(std::ostream & os, Indent indent) const override;
117 
121  void GenerateInputRequestedRegion() override;
122 
124  void EnlargeOutputRequestedRegion( DataObject *itkNotUsed(output) ) override;
125 
128  void GenerateData() override;
129 
130 private:
131  bool m_FullyConnected{ false };
132 
133  bool m_MarkWatershedLine{ true };
134 
136 }; // end of class
137 } // end namespace itk
138 
139 #ifndef ITK_MANUAL_INSTANTIATION
140 #include "itkMorphologicalWatershedImageFilter.hxx"
141 #endif
142 
143 #endif
typename OutputImageType::Pointer OutputImagePointer
typename OutputImageType::ConstPointer OutputImageConstPointer
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
Base class for all process objects that output image data.
typename OutputImageType::PixelType OutputImagePixelType
typename InputImageType::PixelType InputImagePixelType
typename InputImageType::Pointer InputImagePointer
typename OutputImageType::RegionType OutputImageRegionType
TOutputImage OutputImageType
typename InputImageType::RegionType InputImageRegionType
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
typename InputImageType::ConstPointer InputImageConstPointer
Base class for all data objects in ITK.
Watershed segmentation implementation with morphogical operators.