ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
itkBinaryPruningImageFilter.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 itkBinaryPruningImageFilter_h
19#define itkBinaryPruningImageFilter_h
20
24
25namespace itk
26{
53
54template <typename TInputImage, typename TOutputImage>
55class ITK_TEMPLATE_EXPORT BinaryPruningImageFilter : public ImageToImageFilter<TInputImage, TOutputImage>
56{
57public:
58 ITK_DISALLOW_COPY_AND_MOVE(BinaryPruningImageFilter);
59
65
67 itkNewMacro(Self);
68
70 itkOverrideGetNameOfClassMacro(BinaryPruningImageFilter);
71
73 using InputImageType = TInputImage;
74
76 using OutputImageType = TOutputImage;
77
79 using RegionType = typename InputImageType::RegionType;
80
82 using IndexType = typename RegionType::IndexType;
83
85 using PixelType = typename InputImageType::PixelType;
86
88 using SizeType = typename RegionType::SizeType;
89
91 using InputImagePointer = typename InputImageType::ConstPointer;
92
94 using OutputImagePointer = typename OutputImageType::Pointer;
95
98
102
104 itkSetMacro(Iteration, unsigned int);
105 itkGetConstMacro(Iteration, unsigned int);
107
109 static constexpr unsigned int InputImageDimension = TInputImage::ImageDimension;
110 static constexpr unsigned int OutputImageDimension = TOutputImage::ImageDimension;
111
115 itkConceptMacro(IntConvertibleToPixelTypeCheck, (Concept::Convertible<int, PixelType>));
117
118protected:
120 ~BinaryPruningImageFilter() override = default;
121 void
122 PrintSelf(std::ostream & os, Indent indent) const override;
123
125 void
126 GenerateData() override;
127
129 void
131
133 void
135
136private:
137 unsigned int m_Iteration{};
138}; // end of BinaryThinningImageFilter class
139} // end namespace itk
140
141#ifndef ITK_MANUAL_INSTANTIATION
142# include "itkBinaryPruningImageFilter.hxx"
143#endif
144
145#endif
OutputImageType * GetPruning()
SmartPointer< const Self > ConstPointer
void PrintSelf(std::ostream &os, Indent indent) const override
static constexpr unsigned int OutputImageDimension
typename InputImageType::ConstPointer InputImagePointer
typename RegionType::IndexType IndexType
static constexpr unsigned int InputImageDimension
typename OutputImageType::Pointer OutputImagePointer
~BinaryPruningImageFilter() override=default
typename InputImageType::PixelType PixelType
typename InputImageType::RegionType RegionType
ImageToImageFilter< TInputImage, TOutputImage > Superclass
NeighborhoodIterator< TInputImage > NeighborhoodIteratorType
typename RegionType::SizeType SizeType
Control indentation during Print() invocation.
Definition itkIndent.h:50
Defines iteration of a local N-dimensional neighborhood of pixels across an itk::Image.
Implements transparent reference counting.
#define itkConceptMacro(name, concept)
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....