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{
62template <typename TInputImage, typename TLabelImage>
63class ITK_TEMPLATE_EXPORT LabelStatisticsImageFilter : public ImageSink<TInputImage>
64{
65public:
66 ITK_DISALLOW_COPY_AND_MOVE(LabelStatisticsImageFilter);
67
73
75 itkNewMacro(Self);
76
78 itkOverrideGetNameOfClassMacro(LabelStatisticsImageFilter);
79
81 using InputImagePointer = typename TInputImage::Pointer;
82 using RegionType = typename TInputImage::RegionType;
83 using SizeType = typename TInputImage::SizeType;
84 using IndexType = typename TInputImage::IndexType;
85 using PixelType = typename TInputImage::PixelType;
86
88 using LabelImageType = TLabelImage;
89 using LabelImagePointer = typename TLabelImage::Pointer;
90 using LabelRegionType = typename TLabelImage::RegionType;
91 using LabelSizeType = typename TLabelImage::SizeType;
92 using LabelIndexType = typename TLabelImage::IndexType;
93 using LabelPixelType = typename TLabelImage::PixelType;
94
96 static constexpr unsigned int ImageDimension = TInputImage::ImageDimension;
97
100
103
106
108 using BoundingBoxType = std::vector<IndexValueType>;
109
113
119 {
120 public:
121 // default constructor
123 {
124 // initialized to the default values
126 m_Sum = RealType{};
128
129 // Set such that the first pixel encountered can be compared
132
133 // Default these to zero
134 m_Mean = RealType{};
135 m_Sigma = RealType{};
137
138 const unsigned int imageDimension = Self::ImageDimension;
139 m_BoundingBox.resize(imageDimension * 2);
140 for (unsigned int i = 0; i < imageDimension * 2; i += 2)
141 {
144 }
145 m_Histogram = nullptr;
146 }
147
148 // constructor with histogram enabled
149 LabelStatistics(int size, RealType lowerBound, RealType upperBound)
150 {
151 // initialized to the default values
153 m_Sum = RealType{};
155
156 // Set such that the first pixel encountered can be compared
159
160 // Default these to zero
161 m_Mean = RealType{};
162 m_Sigma = RealType{};
164
165 const unsigned int imageDimension = Self::ImageDimension;
166 m_BoundingBox.resize(imageDimension * 2);
167 for (unsigned int i = 0; i < imageDimension * 2; i += 2)
168 {
171 }
172
173 // Histogram
175 typename HistogramType::SizeType hsize;
178 hsize.SetSize(1);
179 lb.SetSize(1);
180 ub.SetSize(1);
181 m_Histogram->SetMeasurementVectorSize(1);
182 hsize[0] = size;
183 lb[0] = lowerBound;
184 ub[0] = upperBound;
185 m_Histogram->Initialize(hsize, lb, ub);
186 }
187
188 // need copy constructor because of smart pointer to histogram
202
204
205 // added for completeness
208 {
209 if (this != &l)
210 {
211 m_Count = l.m_Count;
214 m_Mean = l.m_Mean;
215 m_Sum = l.m_Sum;
217 m_Sigma = l.m_Sigma;
221 }
222 return *this;
223 }
224
225 friend std::ostream &
226 operator<<(std::ostream & os, const LabelStatistics & labelStatistics)
227 {
228 using namespace print_helper;
229
230 os << "Count: " << static_cast<typename NumericTraits<IdentifierType>::PrintType>(labelStatistics.m_Count)
231 << std::endl;
232 os << "Minimum: " << static_cast<typename NumericTraits<RealType>::PrintType>(labelStatistics.m_Minimum)
233 << std::endl;
234 os << "Maximum: " << static_cast<typename NumericTraits<RealType>::PrintType>(labelStatistics.m_Maximum)
235 << std::endl;
236 os << "Mean: " << static_cast<typename NumericTraits<RealType>::PrintType>(labelStatistics.m_Mean) << std::endl;
237 os << "Sum: " << static_cast<typename NumericTraits<RealType>::PrintType>(labelStatistics.m_Sum) << std::endl;
238 os << "SumOfSquares: " << static_cast<typename NumericTraits<RealType>::PrintType>(labelStatistics.m_SumOfSquares)
239 << std::endl;
240 os << "Sigma: " << static_cast<typename NumericTraits<RealType>::PrintType>(labelStatistics.m_Sigma) << std::endl;
241 os << "Variance: " << static_cast<typename NumericTraits<RealType>::PrintType>(labelStatistics.m_Variance)
242 << std::endl;
243 os << "BoundingBox: " << labelStatistics.m_BoundingBox << std::endl;
244
245 os << "Histogram: ";
246 if (labelStatistics.m_Histogram)
247 {
248 labelStatistics.m_Histogram->Print(os);
249 }
250 else
251 {
252 os << "nullptr" << std::endl;
253 }
254
255 return os;
256 }
257
268 };
269
271 using MapType = std::unordered_map<LabelPixelType, LabelStatistics>;
272 using MapIterator = typename MapType::iterator;
273 using MapConstIterator = typename MapType::const_iterator;
275
277 using ValidLabelValuesContainerType = std::vector<LabelPixelType>;
278
279 // macros for Histogram enables
280 itkSetMacro(UseHistograms, bool);
281 itkGetConstMacro(UseHistograms, bool);
282 itkBooleanMacro(UseHistograms);
283
284
285 virtual const ValidLabelValuesContainerType &
287 {
288 return m_ValidLabelValues;
289 }
290
292 itkSetInputMacro(LabelInput, TLabelImage);
293 itkGetInputMacro(LabelInput, TLabelImage);
295
298 bool
300 {
301 return m_LabelStatistics.find(label) != m_LabelStatistics.end();
302 }
303
305 MapSizeType
307 {
308 return static_cast<MapSizeType>(m_LabelStatistics.size());
309 }
310
311 MapSizeType
313 {
314 return static_cast<MapSizeType>(this->GetNumberOfObjects());
315 }
316
320
324
328
333
337
341
347
351
354 GetSum(LabelPixelType label) const;
355
359
363
365 void
366 SetHistogramParameters(const int numBins, RealType lowerBound, RealType upperBound);
367
368 // Change the access from protected to public to expose streaming option, a using statement can not be used due to
369 // limitations of wrapping.
370 void
371 SetNumberOfStreamDivisions(const unsigned int n) override
372 {
374 }
375 unsigned int
377 {
379 }
380
381 itkConceptMacro(InputHasNumericTraitsCheck, (Concept::HasNumericTraits<PixelType>));
382
383protected:
385 ~LabelStatisticsImageFilter() override = default;
386 void
387 PrintSelf(std::ostream & os, Indent indent) const override;
388
389 void
391 {
392 this->AllocateOutputs();
393 m_LabelStatistics.clear();
394 }
395
398 void
400
401 void
403
404private:
405 void
407
410
412
414
417
418 std::mutex m_Mutex{};
419
420}; // end of class
421} // end namespace itk
422
423#ifndef ITK_MANUAL_INSTANTIATION
424# include "itkLabelStatisticsImageFilter.hxx"
425#endif
426
427#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)
typename DataObject::Pointer DataObjectPointer
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
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.
This class stores measurement vectors in the context of n-dimensional histogram.
Array< itk::SizeValueType > SizeType
#define itkConceptMacro(name, concept)
ITKTransform_EXPORT std::ostream & operator<<(std::ostream &out, const TransformBaseTemplateEnums::TransformCategory value)
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
SizeValueType IdentifierType
Definition itkIntTypes.h:90