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 
21 #include "itkInPlaceImageFilter.h"
22 #include "itkImage.h"
23 #include <vector>
24 #include <mutex>
25 
26 namespace itk
27 {
81 template <typename TInputImage, typename TOutputImage>
82 class ITK_TEMPLATE_EXPORT RelabelComponentImageFilter : public InPlaceImageFilter<TInputImage, TOutputImage>
83 {
84 public:
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;
115  using SizeType = typename TInputImage::SizeType;
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);
166  itkSetMacro(MinimumObjectSize, ObjectSizeType);
167 
173  itkGetConstMacro(MinimumObjectSize, ObjectSizeType);
174 
177  itkSetMacro(SortByObjectSize, bool);
178  itkGetConstMacro(SortByObjectSize, bool);
179  itkBooleanMacro(SortByObjectSize);
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.
205  return this->m_SizeOfObjectsInPhysicalUnits;
206  }
207 
211  ObjectSizeType
213  {
214  if (obj > 0 && static_cast<SizeValueType>(obj) <= m_NumberOfObjects)
215  {
216  return m_SizeOfObjectsInPixels[obj - 1];
217  }
218  else
219  {
220  return 0;
221  }
222  }
228  float
230  {
231  if (obj > 0 && static_cast<SizeValueType>(obj) <= m_NumberOfObjects)
232  {
233  return m_SizeOfObjectsInPhysicalUnits[obj - 1];
234  }
235  else
236  {
237  return 0;
238  }
239  }
242 #ifdef ITK_USE_CONCEPT_CHECKING
243  // Begin concept checking
244  itkConceptMacro(InputEqualityComparableCheck, (Concept::EqualityComparable<InputPixelType>));
245  itkConceptMacro(UnsignedLongConvertibleToInputCheck, (Concept::Convertible<LabelType, InputPixelType>));
246  itkConceptMacro(OutputLongConvertibleToUnsignedLongCheck, (Concept::Convertible<OutputPixelType, LabelType>));
249  // End concept checking
250 #endif
251 
252 protected:
254  ~RelabelComponentImageFilter() override = default;
255 
259  void
260  GenerateData() override;
261 
262  void
263  ParallelComputeLabels(const RegionType & inputRegionForThread);
264 
268  void
269  GenerateInputRequestedRegion() override;
270 
272  void
273  PrintSelf(std::ostream & os, Indent indent) const override;
274 
276  {
278 
281  {
282  this->m_SizeInPixels += other.m_SizeInPixels;
283  return *this;
284  }
285  };
286 
287 private:
288  SizeValueType m_NumberOfObjects{ 0 };
289  SizeValueType m_NumberOfObjectsToPrint{ 10 };
290  SizeValueType m_OriginalNumberOfObjects{ 0 };
291  ObjectSizeType m_MinimumObjectSize{ 0 };
292  bool m_SortByObjectSize{ true };
293 
294  std::mutex m_Mutex{};
295 
296  using MapType = std::map<LabelType, RelabelComponentObjectType>;
297  MapType m_SizeMap{};
298 
299  ObjectSizeInPixelsContainerType m_SizeOfObjectsInPixels{};
300  ObjectSizeInPhysicalUnitsContainerType m_SizeOfObjectsInPhysicalUnits{};
301 };
302 } // end namespace itk
303 
304 #ifndef ITK_MANUAL_INSTANTIATION
305 # include "itkRelabelComponentImageFilter.hxx"
306 #endif
307 
308 #endif
itk::RelabelComponentImageFilter::RelabelComponentObjectType
Definition: itkRelabelComponentImageFilter.h:275
itk::RelabelComponentImageFilter::MapType
std::map< LabelType, RelabelComponentObjectType > MapType
Definition: itkRelabelComponentImageFilter.h:296
itk::RelabelComponentImageFilter::RegionType
typename TOutputImage::RegionType RegionType
Definition: itkRelabelComponentImageFilter.h:116
itk::RelabelComponentImageFilter::ObjectSizeInPhysicalUnitsContainerType
std::vector< float > ObjectSizeInPhysicalUnitsContainerType
Definition: itkRelabelComponentImageFilter.h:145
itk::InPlaceImageFilter
Base class for filters that take an image as input and overwrite that image as the output.
Definition: itkInPlaceImageFilter.h:77
itk::RelabelComponentImageFilter::RelabelComponentObjectType::operator+=
RelabelComponentObjectType & operator+=(const RelabelComponentObjectType &other)
Definition: itkRelabelComponentImageFilter.h:280
itk::GTest::TypedefsAndConstructors::Dimension2::SizeType
ImageBaseType::SizeType SizeType
Definition: itkGTestTypedefsAndConstructors.h:49
itkImage.h
itk::SmartPointer< Self >
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::RelabelComponentImageFilter::OutputInternalPixelType
typename TOutputImage::InternalPixelType OutputInternalPixelType
Definition: itkRelabelComponentImageFilter.h:103
itk::RelabelComponentImageFilter::InputInternalPixelType
typename TInputImage::InternalPixelType InputInternalPixelType
Definition: itkRelabelComponentImageFilter.h:105
itk::RelabelComponentImageFilter::GetSizeOfObjectInPixels
ObjectSizeType GetSizeOfObjectInPixels(LabelType obj) const
Definition: itkRelabelComponentImageFilter.h:212
itk::RelabelComponentImageFilter::GetSizeOfObjectInPhysicalUnits
float GetSizeOfObjectInPhysicalUnits(LabelType obj) const
Definition: itkRelabelComponentImageFilter.h:229
itk::RelabelComponentImageFilter::GetSizeOfObjectsInPhysicalUnits
const ObjectSizeInPhysicalUnitsContainerType & GetSizeOfObjectsInPhysicalUnits() const
Definition: itkRelabelComponentImageFilter.h:201
itk::Concept::SameDimension
Definition: itkConceptChecking.h:696
itk::GTest::TypedefsAndConstructors::Dimension2::IndexType
ImageBaseType::IndexType IndexType
Definition: itkGTestTypedefsAndConstructors.h:50
itk::ImageSource
Base class for all process objects that output image data.
Definition: itkImageSource.h:67
itk::RelabelComponentImageFilter::GetSizeOfObjectsInPixels
const ObjectSizeInPixelsContainerType & GetSizeOfObjectsInPixels() const
Definition: itkRelabelComponentImageFilter.h:188
itk::GTest::TypedefsAndConstructors::Dimension2::RegionType
ImageBaseType::RegionType RegionType
Definition: itkGTestTypedefsAndConstructors.h:54
itk::RelabelComponentImageFilter::ObjectSizeInPixelsContainerType
std::vector< ObjectSizeType > ObjectSizeInPixelsContainerType
Definition: itkRelabelComponentImageFilter.h:144
itk::ImageToImageFilter::InputImageType
TInputImage InputImageType
Definition: itkImageToImageFilter.h:129
itk::RelabelComponentImageFilter::OutputPixelType
typename TOutputImage::PixelType OutputPixelType
Definition: itkRelabelComponentImageFilter.h:102
itk::RelabelComponentImageFilter::RelabelComponentObjectType::m_SizeInPixels
ObjectSizeType m_SizeInPixels
Definition: itkRelabelComponentImageFilter.h:277
itkConceptMacro
#define itkConceptMacro(name, concept)
Definition: itkConceptChecking.h:65
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnatomicalOrientation.h:29
itk::RelabelComponentImageFilter
Relabel the components in an image such that consecutive labels are used.
Definition: itkRelabelComponentImageFilter.h:82
itkInPlaceImageFilter.h
itk::ProcessObject
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
Definition: itkProcessObject.h:139
itk::RelabelComponentImageFilter::SizeType
typename TInputImage::SizeType SizeType
Definition: itkRelabelComponentImageFilter.h:115
itk::RelabelComponentImageFilter::IndexType
typename TInputImage::IndexType IndexType
Definition: itkRelabelComponentImageFilter.h:114
itk::Concept::Convertible
Definition: itkConceptChecking.h:216
itk::RelabelComponentImageFilter::InputPixelType
typename TInputImage::PixelType InputPixelType
Definition: itkRelabelComponentImageFilter.h:104
itk::RelabelComponentImageFilter::ObjectSizeType
SizeValueType ObjectSizeType
Definition: itkRelabelComponentImageFilter.h:138
itk::RelabelComponentImageFilter::LabelType
InputPixelType LabelType
Definition: itkRelabelComponentImageFilter.h:135
itk::Concept::EqualityComparable
Definition: itkConceptChecking.h:306
itk::SizeValueType
unsigned long SizeValueType
Definition: itkIntTypes.h:86
itk::ImageSource::OutputImageType
TOutputImage OutputImageType
Definition: itkImageSource.h:90