ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
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{
80
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;
114 using IndexType = typename TInputImage::IndexType;
115 using SizeType = typename TInputImage::SizeType;
116 using RegionType = typename TOutputImage::RegionType;
117
123
127 itkOverrideGetNameOfClassMacro(RelabelComponentImageFilter);
128
132 itkNewMacro(Self);
133
136
139
142 itkGetConstMacro(NumberOfObjects, SizeValueType);
143
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);
159
166 itkSetMacro(MinimumObjectSize, ObjectSizeType);
167
173 itkGetConstMacro(MinimumObjectSize, ObjectSizeType);
174
177 itkSetMacro(SortByObjectSize, bool);
178 itkGetConstMacro(SortByObjectSize, bool);
179 itkBooleanMacro(SortByObjectSize);
181
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
200 const ObjectSizeInPhysicalUnitsContainerType &
202 {
203 // The GetConstReferenceMacro can't be used here because this container
204 // doesn't have an ostream<< operator overloaded.
206 }
207
211 ObjectSizeType
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
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 itkConceptMacro(InputEqualityComparableCheck, (Concept::EqualityComparable<InputPixelType>));
237 itkConceptMacro(UnsignedLongConvertibleToInputCheck, (Concept::Convertible<LabelType, InputPixelType>));
238 itkConceptMacro(OutputLongConvertibleToUnsignedLongCheck, (Concept::Convertible<OutputPixelType, LabelType>));
241
242protected:
244 ~RelabelComponentImageFilter() override = default;
245
249 void
250 GenerateData() override;
251
252 void
253 ParallelComputeLabels(const RegionType & inputRegionForThread);
254
258 void
260
262 void
263 PrintSelf(std::ostream & os, Indent indent) const override;
264
266 {
268
271 {
272 this->m_SizeInPixels += other.m_SizeInPixels;
273 return *this;
274 }
275 };
276
277private:
282 bool m_SortByObjectSize{ true };
283
284 std::mutex m_Mutex{};
285
286 using MapType = std::map<LabelType, RelabelComponentObjectType>;
288
291};
292} // end namespace itk
293
294#ifndef ITK_MANUAL_INSTANTIATION
295# include "itkRelabelComponentImageFilter.hxx"
296#endif
297
298#endif
Control indentation during Print() invocation.
Definition itkIndent.h:50
ObjectSizeType GetSizeOfObjectInPixels(LabelType obj) const
void PrintSelf(std::ostream &os, Indent indent) const override
~RelabelComponentImageFilter() override=default
std::map< LabelType, RelabelComponentObjectType > MapType
void GenerateInputRequestedRegion() override
typename TOutputImage::RegionType RegionType
typename TOutputImage::InternalPixelType OutputInternalPixelType
typename TOutputImage::PixelType OutputPixelType
std::vector< ObjectSizeType > ObjectSizeInPixelsContainerType
InPlaceImageFilter< TInputImage, TOutputImage > Superclass
float GetSizeOfObjectInPhysicalUnits(LabelType obj) const
typename TInputImage::InternalPixelType InputInternalPixelType
typename TInputImage::PixelType InputPixelType
const ObjectSizeInPixelsContainerType & GetSizeOfObjectsInPixels() const
const ObjectSizeInPhysicalUnitsContainerType & GetSizeOfObjectsInPhysicalUnits() const
void ParallelComputeLabels(const RegionType &inputRegionForThread)
Implements transparent reference counting.
#define itkConceptMacro(name, concept)
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
unsigned long SizeValueType
Definition itkIntTypes.h:86
RelabelComponentObjectType & operator+=(const RelabelComponentObjectType &other)