ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
itkWatershedImageFilter.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 itkWatershedImageFilter_h
19#define itkWatershedImageFilter_h
20
21
26
27namespace itk
28{
148template <typename TInputImage>
149class ITK_TEMPLATE_EXPORT WatershedImageFilter
150 : public ImageToImageFilter<TInputImage, Image<IdentifierType, TInputImage::ImageDimension>>
151{
152public:
153 ITK_DISALLOW_COPY_AND_MOVE(WatershedImageFilter);
154
157
159 using InputImageType = TInputImage;
160
162 static constexpr unsigned int ImageDimension = TInputImage::ImageDimension;
163
166
168 using RegionType = typename InputImageType::RegionType;
169 using SizeType = typename InputImageType::SizeType;
170 using IndexType = typename InputImageType::IndexType;
171
174
176 using ScalarType = typename InputImageType::PixelType;
177
180
182 itkOverrideGetNameOfClassMacro(WatershedImageFilter);
183
185 itkNewMacro(Self);
186
188 void
189 GenerateData() override;
190
194 void
195 SetInput(const InputImageType * input) override
196 {
197 // if the input is changed, we'll need to clear the cached tree
198 // when we execute
199 if (input != this->GetInput(0))
200 {
201 m_InputChanged = true;
202 }
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 void
210 SetInput(unsigned int i, const TInputImage * image) override
211 {
212 if (i != 0)
213 {
214 itkExceptionMacro("Filter has only one input.");
215 }
216 else
217 {
218 this->SetInput(image);
219 }
220 }
221
224 void
226
227 itkGetConstMacro(Threshold, double);
228
231 void
232 SetLevel(double);
233
234 itkGetConstMacro(Level, double);
235
239 {
240 m_Segmenter->Update();
241 return m_Segmenter->GetOutputImage();
242 }
243
247 {
248 return m_TreeGenerator->GetOutputSegmentTree();
249 }
250
251 // Override since the filter produces all of its output
252 void
254
255 itkConceptMacro(InputEqualityComparableCheck, (Concept::EqualityComparable<ScalarType>));
256 itkConceptMacro(InputAdditiveOperatorsCheck, (Concept::AdditiveOperators<ScalarType>));
258 itkConceptMacro(InputLessThanComparableCheck, (Concept::LessThanComparable<ScalarType>));
259
260protected:
262 ~WatershedImageFilter() override = default;
263 void
264 PrintSelf(std::ostream & os, Indent indent) const override;
265
268 void
269 PrepareOutputs() override;
270
271private:
275 double m_Threshold{ 0.0 };
276
281 double m_Level{ 0.0 };
282
288
290
292
293 unsigned long m_ObserverTag{};
294
298
300};
301} // end namespace itk
302
303#ifndef ITK_MANUAL_INSTANTIATION
304# include "itkWatershedImageFilter.hxx"
305#endif
306
307#endif
Base class for all data objects in ITK.
virtual void SetInput(const InputImageType *input)
const InputImageType * GetInput() const
Templated n-dimensional image class.
Definition itkImage.h:89
Control indentation during Print() invocation.
Definition itkIndent.h:50
virtual void SetNthInput(DataObjectPointerArraySizeType idx, DataObject *input)
Implements transparent reference counting.
Generate a unique, increasing time value.
typename InputImageType::IndexType IndexType
void SetInput(const InputImageType *input) override
watershed::SegmentTreeGenerator< ScalarType >::SegmentTreeType * GetSegmentTree()
watershed::SegmentTreeGenerator< ScalarType >::Pointer m_TreeGenerator
void GenerateData() override
void PrintSelf(std::ostream &os, Indent indent) const override
ImageToImageFilter< InputImageType, OutputImageType > Superclass
typename InputImageType::PixelType ScalarType
watershed::Relabeler< ScalarType, Self::ImageDimension >::Pointer m_Relabeler
void SetInput(unsigned int i, const TInputImage *image) override
void PrepareOutputs() override
watershed::Segmenter< InputImageType >::Pointer m_Segmenter
void EnlargeOutputRequestedRegion(DataObject *data) override
~WatershedImageFilter() override=default
Image< IdentifierType, Self::ImageDimension > OutputImageType
watershed::Segmenter< InputImageType >::OutputImageType * GetBasicSegmentation()
typename InputImageType::RegionType RegionType
Image< IdentifierType, Self::ImageDimension > OutputImageType
#define itkConceptMacro(name, concept)
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....