ITK  5.0.0
Insight Segmentation and Registration Toolkit
itkRelabelComponentImageFilter.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright Insight Software Consortium
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  * http://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 
25 namespace itk
26 {
80 template< typename TInputImage, typename TOutputImage >
81 class ITK_TEMPLATE_EXPORT RelabelComponentImageFilter:
82  public InPlaceImageFilter< TInputImage, TOutputImage >
83 {
84 public:
85  ITK_DISALLOW_COPY_AND_ASSIGN(RelabelComponentImageFilter);
86 
92 
96  using InputImagePointer = 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 
128 
132  itkNewMacro(Self);
133 
136 
139 
142  itkGetConstMacro(NumberOfObjects, LabelType);
143 
144  using ObjectSizeInPixelsContainerType = std::vector< ObjectSizeType >;
145  using ObjectSizeInPhysicalUnitsContainerType = std::vector< float >;
146 
152  itkGetConstMacro(OriginalNumberOfObjects, LabelType);
153 
156  itkSetMacro(NumberOfObjectsToPrint, LabelType);
157  itkGetConstReferenceMacro(NumberOfObjectsToPrint, LabelType);
159 
166  itkSetMacro(MinimumObjectSize, ObjectSizeType);
167 
173  itkGetConstMacro(MinimumObjectSize, ObjectSizeType);
174 
177  itkSetMacro(SortByObjectSize, bool);
178  itkGetConstMacro(SortByObjectSize, bool);
179  itkBooleanMacro(SortByObjectSize);
181 
188  {
189  // The GetConstReferenceMacro can't be used here becase this container
190  // doesn't have an ostream<< operator overloaded.
191  return this->m_SizeOfObjectsInPixels;
192  }
193 
200  {
201  // The GetConstReferenceMacro can't be used here becase this container
202  // doesn't have an ostream<< operator overloaded.
203  return this->m_SizeOfObjectsInPhysicalUnits;
204  }
205 
210  {
211  if ( obj > 0 && obj <= m_NumberOfObjects )
212  {
213  return m_SizeOfObjectsInPixels[obj - 1];
214  }
215  else
216  {
217  return 0;
218  }
219  }
221 
226  {
227  if ( obj > 0 && obj <= m_NumberOfObjects )
228  {
229  return m_SizeOfObjectsInPhysicalUnits[obj - 1];
230  }
231  else
232  {
233  return 0;
234  }
235  }
237 
238 #ifdef ITK_USE_CONCEPT_CHECKING
239  // Begin concept checking
240  itkConceptMacro( InputEqualityComparableCheck,
242  itkConceptMacro( UnsignedLongConvertibleToInputCheck,
244  itkConceptMacro( OutputLongConvertibleToUnsignedLongCheck,
246  itkConceptMacro( InputConvertibleToOutputCheck,
248  itkConceptMacro( SameDimensionCheck,
250  // End concept checking
251 #endif
252 
253 protected:
254 
256 
257  { this->InPlaceOff(); }
258  ~RelabelComponentImageFilter() override = default;
259 
263  void GenerateData() override;
264 
268  void GenerateInputRequestedRegion() override;
269 
271  void PrintSelf(std::ostream & os, Indent indent) const override;
272 
277  };
278 
279  // put the function objects here for sorting in descending order
281  {
282  public:
284  const RelabelComponentObjectType & b)
285  {
286  if ( a.m_SizeInPixels > b.m_SizeInPixels )
287  {
288  return true;
289  }
290  else if ( a.m_SizeInPixels < b.m_SizeInPixels )
291  {
292  return false;
293  }
294  // size in pixels and physical units are the same, sort based on
295  // original object number
296  else if ( a.m_ObjectNumber < b.m_ObjectNumber )
297  {
298  return true;
299  }
300  else
301  {
302  return false;
303  }
304  }
305  };
306 
307 private:
308 
309  LabelType m_NumberOfObjects{0};
310  LabelType m_NumberOfObjectsToPrint{10};
311  LabelType m_OriginalNumberOfObjects{0};
312  ObjectSizeType m_MinimumObjectSize{0};
313  bool m_SortByObjectSize{true};
314 
317 };
318 } // end namespace itk
319 
320 #ifndef ITK_MANUAL_INSTANTIATION
321 #include "itkRelabelComponentImageFilter.hxx"
322 #endif
323 
324 #endif
const ObjectSizeInPixelsContainerType & GetSizeOfObjectsInPixels() const
ObjectSizeInPhysicalUnitsContainerType m_SizeOfObjectsInPhysicalUnits
unsigned long SizeValueType
Definition: itkIntTypes.h:83
typename TInputImage::InternalPixelType InputInternalPixelType
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
Base class for all process objects that output image data.
typename TOutputImage::PixelType OutputPixelType
typename TInputImage::IndexType IndexType
typename InputImageType::Pointer InputImagePointer
typename TInputImage::PixelType InputPixelType
SizeValueType IdentifierType
Definition: itkIntTypes.h:87
typename TOutputImage::InternalPixelType OutputInternalPixelType
TOutputImage OutputImageType
bool operator()(const RelabelComponentObjectType &a, const RelabelComponentObjectType &b)
ObjectSizeType GetSizeOfObjectInPixels(LabelType obj) const
float GetSizeOfObjectInPhysicalUnits(LabelType obj) const
Relabel the components in an image such that consecutive labels are used.
std::vector< ObjectSizeType > ObjectSizeInPixelsContainerType
const ObjectSizeInPhysicalUnitsContainerType & GetSizeOfObjectsInPhysicalUnits() const
Base class for filters that take an image as input and produce an image as output.
Control indentation during Print() invocation.
Definition: itkIndent.h:49
Base class for filters that take an image as input and overwrite that image as the output...
#define itkConceptMacro(name, concept)
ObjectSizeInPixelsContainerType m_SizeOfObjectsInPixels
typename TOutputImage::RegionType RegionType