ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
itkLabelStatisticsImageFilter.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 itkLabelStatisticsImageFilter_h
19#define itkLabelStatisticsImageFilter_h
20
21#include "itkImageSink.h"
22#include "itkNumericTraits.h"
24#include "itkHistogram.h"
25#include "itkPrintHelper.h"
26#include <mutex>
27#include <unordered_map>
28#include <vector>
29
30namespace itk
31{
61template <typename TInputImage, typename TLabelImage>
62class ITK_TEMPLATE_EXPORT LabelStatisticsImageFilter : public ImageSink<TInputImage>
63{
64public:
65 ITK_DISALLOW_COPY_AND_MOVE(LabelStatisticsImageFilter);
66
72
74 itkNewMacro(Self);
75
77 itkOverrideGetNameOfClassMacro(LabelStatisticsImageFilter);
78
80 using InputImagePointer = typename TInputImage::Pointer;
81 using RegionType = typename TInputImage::RegionType;
82 using SizeType = typename TInputImage::SizeType;
83 using IndexType = typename TInputImage::IndexType;
84 using PixelType = typename TInputImage::PixelType;
85
87 using LabelImageType = TLabelImage;
88 using LabelImagePointer = typename TLabelImage::Pointer;
89 using LabelRegionType = typename TLabelImage::RegionType;
90 using LabelSizeType = typename TLabelImage::SizeType;
91 using LabelIndexType = typename TLabelImage::IndexType;
92 using LabelPixelType = typename TLabelImage::PixelType;
93
95 static constexpr unsigned int ImageDimension = TInputImage::ImageDimension;
96
99
102
105
107 using BoundingBoxType = std::vector<IndexValueType>;
108
112
118 {
119 public:
120 // default constructor
122 // initialized to the default values
125 , m_Maximum(NumericTraits<RealType>::NonpositiveMin())
126 , m_Mean(RealType{})
127 , m_Sum(RealType{})
129 , m_Sigma(RealType{})
131 {
132 // Set such that the first pixel encountered can be compared
133
134 const unsigned int imageDimension = Self::ImageDimension;
135 m_BoundingBox.resize(imageDimension * 2);
136 for (unsigned int i = 0; i < imageDimension * 2; i += 2)
137 {
140 }
141 m_Histogram = nullptr;
142 }
143
144 // constructor with histogram enabled
145 LabelStatistics(int size, RealType lowerBound, RealType upperBound)
148 , m_Maximum(NumericTraits<RealType>::NonpositiveMin())
149 , m_Mean(RealType{})
150 , m_Sum(RealType{})
152 , m_Sigma(RealType{})
154 {
155 // Set such that the first pixel encountered can be compared
156 const unsigned int imageDimension = Self::ImageDimension;
157 m_BoundingBox.resize(imageDimension * 2);
158 for (unsigned int i = 0; i < imageDimension * 2; i += 2)
159 {
162 }
163
164 // Histogram
166 typename HistogramType::SizeType hsize;
169 hsize.SetSize(1);
170 lb.SetSize(1);
171 ub.SetSize(1);
172 m_Histogram->SetMeasurementVectorSize(1);
173 hsize[0] = size;
174 lb[0] = lowerBound;
175 ub[0] = upperBound;
176 m_Histogram->Initialize(hsize, lb, ub);
177 }
178
179 // need copy constructor because of smart pointer to histogram
192
194
195 // added for completeness
198 {
199 if (this != &l)
200 {
201 m_Count = l.m_Count;
204 m_Mean = l.m_Mean;
205 m_Sum = l.m_Sum;
207 m_Sigma = l.m_Sigma;
211 }
212 return *this;
213 }
214
215 friend std::ostream &
216 operator<<(std::ostream & os, const LabelStatistics & labelStatistics)
217 {
218 using namespace print_helper;
219
220 os << "Count: " << static_cast<NumericTraits<IdentifierType>::PrintType>(labelStatistics.m_Count) << std::endl;
221 os << "Minimum: " << static_cast<typename NumericTraits<RealType>::PrintType>(labelStatistics.m_Minimum)
222 << std::endl;
223 os << "Maximum: " << static_cast<typename NumericTraits<RealType>::PrintType>(labelStatistics.m_Maximum)
224 << std::endl;
225 os << "Mean: " << static_cast<typename NumericTraits<RealType>::PrintType>(labelStatistics.m_Mean) << std::endl;
226 os << "Sum: " << static_cast<typename NumericTraits<RealType>::PrintType>(labelStatistics.m_Sum) << std::endl;
227 os << "SumOfSquares: " << static_cast<typename NumericTraits<RealType>::PrintType>(labelStatistics.m_SumOfSquares)
228 << std::endl;
229 os << "Sigma: " << static_cast<typename NumericTraits<RealType>::PrintType>(labelStatistics.m_Sigma) << std::endl;
230 os << "Variance: " << static_cast<typename NumericTraits<RealType>::PrintType>(labelStatistics.m_Variance)
231 << std::endl;
232 os << "BoundingBox: " << labelStatistics.m_BoundingBox << std::endl;
233
234 os << "Histogram: ";
235 if (labelStatistics.m_Histogram)
236 {
237 labelStatistics.m_Histogram->Print(os);
238 }
239 else
240 {
241 os << "nullptr" << std::endl;
242 }
243
244 return os;
245 }
246
257 };
258
260 using MapType = std::unordered_map<LabelPixelType, LabelStatistics>;
261 using MapIterator = typename MapType::iterator;
262 using MapConstIterator = typename MapType::const_iterator;
264
266 using ValidLabelValuesContainerType = std::vector<LabelPixelType>;
267
268 // macros for Histogram enables
269 itkSetMacro(UseHistograms, bool);
270 itkGetConstMacro(UseHistograms, bool);
271 itkBooleanMacro(UseHistograms);
272
273
274 virtual const ValidLabelValuesContainerType &
276 {
277 return m_ValidLabelValues;
278 }
279
282 itkSetInputMacro(LabelInput, TLabelImage);
283 itkGetInputMacro(LabelInput, TLabelImage);
287 bool
289 {
290 return m_LabelStatistics.find(label) != m_LabelStatistics.end();
291 }
292
294 [[nodiscard]] MapSizeType
296 {
297 return static_cast<MapSizeType>(m_LabelStatistics.size());
298 }
299
300 [[nodiscard]] MapSizeType
302 {
303 return static_cast<MapSizeType>(this->GetNumberOfObjects());
304 }
305
309
313
317
322
326
330
336
340
343 GetSum(LabelPixelType label) const;
344
348
352
354 void
355 SetHistogramParameters(const int numBins, RealType lowerBound, RealType upperBound);
356
357 // Change the access from protected to public to expose streaming option, a using statement can not be used due to
358 // limitations of wrapping.
359 void
360 SetNumberOfStreamDivisions(const unsigned int n) override
361 {
363 }
364 [[nodiscard]] unsigned int
366 {
368 }
369
370 itkConceptMacro(InputHasNumericTraitsCheck, (Concept::HasNumericTraits<PixelType>));
371
372protected:
374 ~LabelStatisticsImageFilter() override = default;
375 void
376 PrintSelf(std::ostream & os, Indent indent) const override;
377
378 void
380 {
381 this->AllocateOutputs();
382 m_LabelStatistics.clear();
383 }
384
387 void
389
390 void
392
393private:
394 void
396
399
401
403
406
407 std::mutex m_Mutex{};
408
409}; // end of class
410} // end namespace itk
411
412#ifndef ITK_MANUAL_INSTANTIATION
413# include "itkLabelStatisticsImageFilter.hxx"
414#endif
415
416#endif
void SetSize(SizeValueType sz)
SmartPointer< Self > Pointer
virtual void AllocateOutputs()
virtual void SetNumberOfStreamDivisions(unsigned int _arg)
virtual unsigned int GetNumberOfStreamDivisions() const
Control indentation during Print() invocation.
Definition itkIndent.h:50
LabelStatistics(int size, RealType lowerBound, RealType upperBound)
LabelStatistics & operator=(const LabelStatistics &l)
void SetNumberOfStreamDivisions(const unsigned int n) override
~LabelStatisticsImageFilter() override=default
HistogramPointer GetHistogram(LabelPixelType label) const
RealType GetMean(LabelPixelType label) const
unsigned int GetNumberOfStreamDivisions() const override
typename TLabelImage::RegionType LabelRegionType
itk::Statistics::Histogram< RealType > HistogramType
typename HistogramType::Pointer HistogramPointer
typename TInputImage::IndexType IndexType
SimpleDataObjectDecorator< RealType > RealObjectType
RealType GetSigma(LabelPixelType label) const
BoundingBoxType GetBoundingBox(LabelPixelType label) const
void MergeMap(MapType &, MapType &) const
MapSizeType GetCount(LabelPixelType label) const
typename TLabelImage::Pointer LabelImagePointer
RealType GetMinimum(LabelPixelType label) const
typename TInputImage::RegionType RegionType
RegionType GetRegion(LabelPixelType label) const
bool HasLabel(LabelPixelType label) const
typename TLabelImage::PixelType LabelPixelType
itkSetInputMacro(LabelInput, TLabelImage)
std::vector< IndexValueType > BoundingBoxType
ValidLabelValuesContainerType m_ValidLabelValues
void AfterStreamedGenerateData() override
std::vector< LabelPixelType > ValidLabelValuesContainerType
typename TLabelImage::IndexType LabelIndexType
itkGetInputMacro(LabelInput, TLabelImage)
static constexpr unsigned int ImageDimension
typename TLabelImage::SizeType LabelSizeType
RealType GetSum(LabelPixelType label) const
typename TInputImage::PixelType PixelType
RealType GetMaximum(LabelPixelType label) const
typename TInputImage::Pointer InputImagePointer
std::unordered_map< LabelPixelType, LabelStatistics > MapType
virtual const ValidLabelValuesContainerType & GetValidLabelValues() const
typename TInputImage::SizeType SizeType
void PrintSelf(std::ostream &os, Indent indent) const override
RealType GetMedian(LabelPixelType label) const
typename NumericTraits< PixelType >::RealType RealType
void SetHistogramParameters(const int numBins, RealType lowerBound, RealType upperBound)
typename MapType::const_iterator MapConstIterator
void ThreadedStreamedGenerateData(const RegionType &) override
RealType GetVariance(LabelPixelType label) const
Define additional traits for native types such as int or float.
static constexpr T NonpositiveMin()
static constexpr T max(const T &)
Decorates any "simple" data type (data types without smart pointers) with a DataObject API.
Implements transparent reference counting.
ObjectType * Print(std::ostream &os) const
This class stores measurement vectors in the context of n-dimensional histogram.
Array< itk::SizeValueType > SizeType
#define itkConceptMacro(name, concept)
std::ostream & operator<<(std::ostream &os, const itk::DOMNode &object)
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
SizeValueType IdentifierType
Definition itkIntTypes.h:90