ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
itkVotingBinaryImageFilter.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 itkVotingBinaryImageFilter_h
19#define itkVotingBinaryImageFilter_h
20
22
23namespace itk
24{
39template <typename TInputImage, typename TOutputImage>
40class ITK_TEMPLATE_EXPORT VotingBinaryImageFilter : public ImageToImageFilter<TInputImage, TOutputImage>
41{
42public:
43 ITK_DISALLOW_COPY_AND_MOVE(VotingBinaryImageFilter);
44
46 static constexpr unsigned int InputImageDimension = TInputImage::ImageDimension;
47 static constexpr unsigned int OutputImageDimension = TOutputImage::ImageDimension;
48
50 using InputImageType = TInputImage;
51 using OutputImageType = TOutputImage;
52
58
60 itkNewMacro(Self);
61
63 itkOverrideGetNameOfClassMacro(VotingBinaryImageFilter);
64
66 using InputPixelType = typename InputImageType::PixelType;
67 using OutputPixelType = typename OutputImageType::PixelType;
68
69 using InputImageRegionType = typename InputImageType::RegionType;
70 using OutputImageRegionType = typename OutputImageType::RegionType;
71
72 using InputSizeType = typename InputImageType::SizeType;
73
75 itkSetMacro(Radius, InputSizeType);
76
78 itkGetConstReferenceMacro(Radius, InputSizeType);
79
82 itkSetMacro(BackgroundValue, InputPixelType);
83 itkSetMacro(ForegroundValue, InputPixelType);
85
88 itkGetConstReferenceMacro(BackgroundValue, InputPixelType);
89 itkGetConstReferenceMacro(ForegroundValue, InputPixelType);
91
94 itkGetConstReferenceMacro(BirthThreshold, unsigned int);
95 itkSetMacro(BirthThreshold, unsigned int);
97
100 itkGetConstReferenceMacro(SurvivalThreshold, unsigned int);
101 itkSetMacro(SurvivalThreshold, unsigned int);
103
110 void
112
113 itkConceptMacro(InputEqualityComparableCheck, (Concept::EqualityComparable<InputPixelType>));
114 itkConceptMacro(IntConvertibleToInputCheck, (Concept::Convertible<int, InputPixelType>));
117 itkConceptMacro(InputOStreamWritableCheck, (Concept::OStreamWritable<InputPixelType>));
118
119protected:
121 ~VotingBinaryImageFilter() override = default;
122 void
123 PrintSelf(std::ostream & os, Indent indent) const override;
124
135 void
136 DynamicThreadedGenerateData(const OutputImageRegionType & outputRegionForThread) override;
137
138private:
140
143
144 unsigned int m_BirthThreshold{};
145 unsigned int m_SurvivalThreshold{};
146};
147} // end namespace itk
148
149#ifndef ITK_MANUAL_INSTANTIATION
150# include "itkVotingBinaryImageFilter.hxx"
151#endif
152
153#endif
Control indentation during Print() invocation.
Definition itkIndent.h:50
Implements transparent reference counting.
~VotingBinaryImageFilter() override=default
void DynamicThreadedGenerateData(const OutputImageRegionType &outputRegionForThread) override
void PrintSelf(std::ostream &os, Indent indent) const override
void GenerateInputRequestedRegion() override
ImageToImageFilter< InputImageType, OutputImageType > Superclass
#define itkConceptMacro(name, concept)
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....