ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
itkBinaryThinningImageFilter.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 itkBinaryThinningImageFilter_h
19#define itkBinaryThinningImageFilter_h
20
23
24namespace itk
25{
55
56template <typename TInputImage, typename TOutputImage>
57class ITK_TEMPLATE_EXPORT BinaryThinningImageFilter : public ImageToImageFilter<TInputImage, TOutputImage>
58{
59public:
60 ITK_DISALLOW_COPY_AND_MOVE(BinaryThinningImageFilter);
61
67
69 itkNewMacro(Self);
70
72 itkOverrideGetNameOfClassMacro(BinaryThinningImageFilter);
73
75 using InputImageType = TInputImage;
76
78 using OutputImageType = TOutputImage;
79
81 using RegionType = typename InputImageType::RegionType;
82
84 using IndexType = typename RegionType::IndexType;
85
87 using PixelType = typename InputImageType::PixelType;
88
90 using SizeType = typename RegionType::SizeType;
91
93 using InputImagePointer = typename InputImageType::ConstPointer;
94
96 using OutputImagePointer = typename OutputImageType::Pointer;
97
100
103
107
109 static constexpr unsigned int InputImageDimension = TInputImage::ImageDimension;
110 static constexpr unsigned int OutputImageDimension = TOutputImage::ImageDimension;
111
113 itkConceptMacro(InputAdditiveOperatorsCheck, (Concept::AdditiveOperators<PixelType>));
114 itkConceptMacro(InputConvertibleToIntCheck, (Concept::Convertible<PixelType, int>));
115 itkConceptMacro(IntConvertibleToInputCheck, (Concept::Convertible<int, PixelType>));
117
118protected:
120 ~BinaryThinningImageFilter() override = default;
121
123 void
124 GenerateData() override;
125
127 void
129
131 void
133}; // end of BinaryThinningImageFilter
134 // class
135} // end namespace itk
136
137#ifndef ITK_MANUAL_INSTANTIATION
138# include "itkBinaryThinningImageFilter.hxx"
139#endif
140
141#endif
OutputImageType * GetThinning()
typename OutputImageType::Pointer OutputImagePointer
typename InputImageType::ConstPointer InputImagePointer
typename RegionType::SizeType SizeType
static constexpr unsigned int OutputImageDimension
typename InputImageType::PixelType PixelType
static constexpr unsigned int InputImageDimension
NeighborhoodIterator< TInputImage > NeighborhoodIteratorType
~BinaryThinningImageFilter() override=default
ImageToImageFilter< TInputImage, TOutputImage > Superclass
typename InputImageType::RegionType RegionType
typename RegionType::IndexType IndexType
typename OutputImageType::PixelType OutputImagePixelType
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....