ITK  5.1.0
Insight Toolkit
itkMorphologicalWatershedFromMarkersImageFilter.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  * 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 itkMorphologicalWatershedFromMarkersImageFilter_h
19 #define itkMorphologicalWatershedFromMarkersImageFilter_h
20 
21 #include "itkImageToImageFilter.h"
22 
23 namespace itk
24 {
81 template <typename TInputImage, typename TLabelImage>
83  : public ImageToImageFilter<TInputImage, TLabelImage>
84 {
85 public:
86  ITK_DISALLOW_COPY_AND_ASSIGN(MorphologicalWatershedFromMarkersImageFilter);
87 
93 
95  using InputImageType = TInputImage;
96  using LabelImageType = TLabelImage;
97  using InputImagePointer = typename InputImageType::Pointer;
98  using InputImageConstPointer = typename InputImageType::ConstPointer;
100  using InputImagePixelType = typename InputImageType::PixelType;
101  using LabelImagePointer = typename LabelImageType::Pointer;
102  using LabelImageConstPointer = typename LabelImageType::ConstPointer;
104  using LabelImagePixelType = typename LabelImageType::PixelType;
105 
107 
109  static constexpr unsigned int ImageDimension = TInputImage::ImageDimension;
110 
112  itkNewMacro(Self);
113 
116 
118  void
119  SetMarkerImage(const TLabelImage * input)
120  {
121  // Process object is not const-correct so the const casting is required.
122  this->SetNthInput(1, const_cast<TLabelImage *>(input));
123  }
124 
126  const LabelImageType *
128  {
129  return itkDynamicCastInDebugMode<LabelImageType *>(const_cast<DataObject *>(this->ProcessObject::GetInput(1)));
130  }
131 
133  void
134  SetInput1(const TInputImage * input)
135  {
136  this->SetInput(input);
137  }
138 
140  void
141  SetInput2(const TLabelImage * input)
142  {
143  this->SetMarkerImage(input);
144  }
145 
152  itkSetMacro(FullyConnected, bool);
153  itkGetConstReferenceMacro(FullyConnected, bool);
154  itkBooleanMacro(FullyConnected);
156 
162  itkSetMacro(MarkWatershedLine, bool);
163  itkGetConstReferenceMacro(MarkWatershedLine, bool);
164  itkBooleanMacro(MarkWatershedLine);
166 
167 protected:
169  ~MorphologicalWatershedFromMarkersImageFilter() override = default;
170  void
171  PrintSelf(std::ostream & os, Indent indent) const override;
172 
176  void
177  GenerateInputRequestedRegion() override;
178 
182  void
183  EnlargeOutputRequestedRegion(DataObject * itkNotUsed(output)) override;
184 
186  void
187  GenerateData() override;
188 
189 private:
190  bool m_FullyConnected{ false };
191 
192  bool m_MarkWatershedLine{ true };
193 }; // end of class
194 } // end namespace itk
195 
196 #ifndef ITK_MANUAL_INSTANTIATION
197 # include "itkMorphologicalWatershedFromMarkersImageFilter.hxx"
198 #endif
199 
200 #endif
itk::MorphologicalWatershedFromMarkersImageFilter::LabelImageConstPointer
typename LabelImageType::ConstPointer LabelImageConstPointer
Definition: itkMorphologicalWatershedFromMarkersImageFilter.h:102
itk::MorphologicalWatershedFromMarkersImageFilter
Morphological watershed transform from markers.
Definition: itkMorphologicalWatershedFromMarkersImageFilter.h:82
itk::MorphologicalWatershedFromMarkersImageFilter::LabelImageType
TLabelImage LabelImageType
Definition: itkMorphologicalWatershedFromMarkersImageFilter.h:96
itk::MorphologicalWatershedFromMarkersImageFilter::SetMarkerImage
void SetMarkerImage(const TLabelImage *input)
Definition: itkMorphologicalWatershedFromMarkersImageFilter.h:119
itk::MorphologicalWatershedFromMarkersImageFilter::InputImageConstPointer
typename InputImageType::ConstPointer InputImageConstPointer
Definition: itkMorphologicalWatershedFromMarkersImageFilter.h:98
itk::MorphologicalWatershedFromMarkersImageFilter::InputImagePointer
typename InputImageType::Pointer InputImagePointer
Definition: itkMorphologicalWatershedFromMarkersImageFilter.h:97
itk::SmartPointer< Self >
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::MorphologicalWatershedFromMarkersImageFilter::SetInput2
void SetInput2(const TLabelImage *input)
Definition: itkMorphologicalWatershedFromMarkersImageFilter.h:141
itk::MorphologicalWatershedFromMarkersImageFilter::LabelImagePixelType
typename LabelImageType::PixelType LabelImagePixelType
Definition: itkMorphologicalWatershedFromMarkersImageFilter.h:104
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::LightObject
Light weight base class for most itk classes.
Definition: itkLightObject.h:59
itk::MorphologicalWatershedFromMarkersImageFilter::LabelImageRegionType
typename LabelImageType::RegionType LabelImageRegionType
Definition: itkMorphologicalWatershedFromMarkersImageFilter.h:103
itk::GTest::TypedefsAndConstructors::Dimension2::RegionType
ImageBaseType::RegionType RegionType
Definition: itkGTestTypedefsAndConstructors.h:54
itk::MorphologicalWatershedFromMarkersImageFilter::LabelImagePointer
typename LabelImageType::Pointer LabelImagePointer
Definition: itkMorphologicalWatershedFromMarkersImageFilter.h:101
itkImageToImageFilter.h
itk::MorphologicalWatershedFromMarkersImageFilter::InputImageRegionType
typename InputImageType::RegionType InputImageRegionType
Definition: itkMorphologicalWatershedFromMarkersImageFilter.h:99
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkArray.h:26
itk::MorphologicalWatershedFromMarkersImageFilter::InputImageType
TInputImage InputImageType
Definition: itkMorphologicalWatershedFromMarkersImageFilter.h:95
itk::MorphologicalWatershedFromMarkersImageFilter::InputImagePixelType
typename InputImageType::PixelType InputImagePixelType
Definition: itkMorphologicalWatershedFromMarkersImageFilter.h:100
itk::MorphologicalWatershedFromMarkersImageFilter::GetMarkerImage
const LabelImageType * GetMarkerImage() const
Definition: itkMorphologicalWatershedFromMarkersImageFilter.h:127
itk::MorphologicalWatershedFromMarkersImageFilter::IndexType
typename LabelImageType::IndexType IndexType
Definition: itkMorphologicalWatershedFromMarkersImageFilter.h:106
itk::ProcessObject::GetInput
DataObject * GetInput(const DataObjectIdentifierType &key)
Return an input.
itk::MorphologicalWatershedFromMarkersImageFilter::SetInput1
void SetInput1(const TInputImage *input)
Definition: itkMorphologicalWatershedFromMarkersImageFilter.h:134
itk::DataObject
Base class for all data objects in ITK.
Definition: itkDataObject.h:293