ITK  6.0.0
Insight Toolkit
itkRelabelComponentImageFilter.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 itkRelabelComponentImageFilter_h
19#define itkRelabelComponentImageFilter_h
20
22#include "itkImage.h"
23#include <vector>
24#include <mutex>
25
26namespace itk
27{
81template <typename TInputImage, typename TOutputImage>
82class ITK_TEMPLATE_EXPORT RelabelComponentImageFilter : public InPlaceImageFilter<TInputImage, TOutputImage>
83{
84public:
85 ITK_DISALLOW_COPY_AND_MOVE(RelabelComponentImageFilter);
86
92
96 using typename Superclass::InputImagePointer;
97
102 using OutputPixelType = typename TOutputImage::PixelType;
103 using OutputInternalPixelType = typename TOutputImage::InternalPixelType;
104 using InputPixelType = typename TInputImage::PixelType;
105 using InputInternalPixelType = typename TInputImage::InternalPixelType;
106 static constexpr unsigned int ImageDimension = TOutputImage::ImageDimension;
107 static constexpr unsigned int InputImageDimension = TInputImage::ImageDimension;
108
112 using InputImageType = TInputImage;
113 using OutputImageType = TOutputImage;
116 using RegionType = typename TOutputImage::RegionType;
117
123
127 itkOverrideGetNameOfClassMacro(RelabelComponentImageFilter);
128
132 itkNewMacro(Self);
133
136
139
142 itkGetConstMacro(NumberOfObjects, SizeValueType);
144 using ObjectSizeInPixelsContainerType = std::vector<ObjectSizeType>;
145 using ObjectSizeInPhysicalUnitsContainerType = std::vector<float>;
146
152 itkGetConstMacro(OriginalNumberOfObjects, SizeValueType);
153
156 itkSetMacro(NumberOfObjectsToPrint, SizeValueType);
157 itkGetConstReferenceMacro(NumberOfObjectsToPrint, SizeValueType);
166 itkSetMacro(MinimumObjectSize, ObjectSizeType);
167
173 itkGetConstMacro(MinimumObjectSize, ObjectSizeType);
174
177 itkSetMacro(SortByObjectSize, bool);
178 itkGetConstMacro(SortByObjectSize, bool);
179 itkBooleanMacro(SortByObjectSize);
188 GetSizeOfObjectsInPixels() const
189 {
190 // The GetConstReferenceMacro can't be used here because this container
191 // doesn't have an ostream<< operator overloaded.
192 return this->m_SizeOfObjectsInPixels;
193 }
194
201 GetSizeOfObjectsInPhysicalUnits() const
202 {
203 // The GetConstReferenceMacro can't be used here because this container
204 // doesn't have an ostream<< operator overloaded.
205 return this->m_SizeOfObjectsInPhysicalUnits;
206 }
207
212 GetSizeOfObjectInPixels(LabelType obj) const
213 {
214 if (obj > 0 && static_cast<SizeValueType>(obj) <= m_NumberOfObjects)
215 {
216 return m_SizeOfObjectsInPixels[obj - 1];
217 }
218
219 return 0;
220 }
221
225 float
226 GetSizeOfObjectInPhysicalUnits(LabelType obj) const
227 {
228 if (obj > 0 && static_cast<SizeValueType>(obj) <= m_NumberOfObjects)
229 {
230 return m_SizeOfObjectsInPhysicalUnits[obj - 1];
231 }
232
233 return 0;
234 }
235
236#ifdef ITK_USE_CONCEPT_CHECKING
237 // Begin concept checking
238 itkConceptMacro(InputEqualityComparableCheck, (Concept::EqualityComparable<InputPixelType>));
239 itkConceptMacro(UnsignedLongConvertibleToInputCheck, (Concept::Convertible<LabelType, InputPixelType>));
240 itkConceptMacro(OutputLongConvertibleToUnsignedLongCheck, (Concept::Convertible<OutputPixelType, LabelType>));
243 // End concept checking
244#endif
245
246protected:
248 ~RelabelComponentImageFilter() override = default;
249
253 void
254 GenerateData() override;
255
256 void
257 ParallelComputeLabels(const RegionType & inputRegionForThread);
258
262 void
263 GenerateInputRequestedRegion() override;
264
266 void
267 PrintSelf(std::ostream & os, Indent indent) const override;
271 ObjectSizeType m_SizeInPixels;
272
274 operator+=(const RelabelComponentObjectType & other)
275 {
276 this->m_SizeInPixels += other.m_SizeInPixels;
277 return *this;
278 }
279 };
280
281private:
282 SizeValueType m_NumberOfObjects{ 0 };
283 SizeValueType m_NumberOfObjectsToPrint{ 10 };
284 SizeValueType m_OriginalNumberOfObjects{ 0 };
285 ObjectSizeType m_MinimumObjectSize{ 0 };
286 bool m_SortByObjectSize{ true };
288 std::mutex m_Mutex{};
290 using MapType = std::map<LabelType, RelabelComponentObjectType>;
291 MapType m_SizeMap{};
293 ObjectSizeInPixelsContainerType m_SizeOfObjectsInPixels{};
294 ObjectSizeInPhysicalUnitsContainerType m_SizeOfObjectsInPhysicalUnits{};
295};
296} // end namespace itk
297
298#ifndef ITK_MANUAL_INSTANTIATION
299# include "itkRelabelComponentImageFilter.hxx"
300#endif
301
302#endif
Base class for all process objects that output image data.
TOutputImage OutputImageType
Base class for filters that take an image as input and overwrite that image as the output.
Control indentation during Print() invocation.
Definition: itkIndent.h:50
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
Relabel the components in an image such that consecutive labels are used.
std::map< LabelType, RelabelComponentObjectType > MapType
typename TOutputImage::RegionType RegionType
typename TOutputImage::InternalPixelType OutputInternalPixelType
typename TOutputImage::PixelType OutputPixelType
std::vector< ObjectSizeType > ObjectSizeInPixelsContainerType
typename TInputImage::IndexType IndexType
typename TInputImage::InternalPixelType InputInternalPixelType
typename TInputImage::PixelType InputPixelType
#define itkConceptMacro(name, concept)
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
unsigned long SizeValueType
Definition: itkIntTypes.h:86