ITK  4.13.0
Insight Segmentation and Registration Toolkit
itkLabelStatisticsImageFilter.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 itkLabelStatisticsImageFilter_h
19 #define itkLabelStatisticsImageFilter_h
20 
21 #include "itkImageToImageFilter.h"
22 #include "itkNumericTraits.h"
24 #include "itksys/hash_map.hxx"
25 #include "itkHistogram.h"
26 #include "itkFastMutexLock.h"
27 #include <vector>
28 
29 namespace itk
30 {
59 template< typename TInputImage, typename TLabelImage >
60 class ITK_TEMPLATE_EXPORT LabelStatisticsImageFilter:
61  public ImageToImageFilter< TInputImage, TInputImage >
62 {
63 public:
69 
71  itkNewMacro(Self);
72 
75 
77  typedef typename TInputImage::Pointer InputImagePointer;
78  typedef typename TInputImage::RegionType RegionType;
79  typedef typename TInputImage::SizeType SizeType;
81  typedef typename TInputImage::PixelType PixelType;
82 
84  typedef TLabelImage LabelImageType;
85  typedef typename TLabelImage::Pointer LabelImagePointer;
86  typedef typename TLabelImage::RegionType LabelRegionType;
89  typedef typename TLabelImage::PixelType LabelPixelType;
90 
92  itkStaticConstMacro(ImageDimension, unsigned int,
93  TInputImage::ImageDimension);
94 
97 
100 
103 
105  typedef std::vector< IndexValueType > BoundingBoxType;
106 
110 
116  {
117 public:
118 
119  // default constructor
121  {
122  // initialized to the default values
125  m_SumOfSquares = NumericTraits< RealType >::ZeroValue();
126 
127  // Set such that the first pixel encountered can be compared
128  m_Minimum = NumericTraits< RealType >::max();
130 
131  // Default these to zero
135 
136  const unsigned int imageDimension = itkGetStaticConstMacro(ImageDimension);
137  m_BoundingBox.resize(imageDimension * 2);
138  for ( unsigned int i = 0; i < imageDimension * 2; i += 2 )
139  {
140  m_BoundingBox[i] = NumericTraits< IndexValueType >::max();
141  m_BoundingBox[i + 1] = NumericTraits< IndexValueType >::NonpositiveMin();
142  }
143  m_Histogram = ITK_NULLPTR;
144  }
145 
146  // constructor with histogram enabled
147  LabelStatistics(int size, RealType lowerBound, RealType upperBound)
148  {
149  // initialized to the default values
152  m_SumOfSquares = NumericTraits< RealType >::ZeroValue();
153 
154  // Set such that the first pixel encountered can be compared
155  m_Minimum = NumericTraits< RealType >::max();
157 
158  // Default these to zero
162 
163  const unsigned int imageDimension = itkGetStaticConstMacro(ImageDimension);
164  m_BoundingBox.resize(imageDimension * 2);
165  for ( unsigned int i = 0; i < imageDimension * 2; i += 2 )
166  {
167  m_BoundingBox[i] = NumericTraits< IndexValueType >::max();
168  m_BoundingBox[i + 1] = NumericTraits< IndexValueType >::NonpositiveMin();
169  }
170 
171  // Histogram
172  m_Histogram = HistogramType::New();
173  typename HistogramType::SizeType hsize;
176  hsize.SetSize(1);
177  lb.SetSize(1);
178  ub.SetSize(1);
179  m_Histogram->SetMeasurementVectorSize(1);
180  hsize[0] = size;
181  lb[0] = lowerBound;
182  ub[0] = upperBound;
183  m_Histogram->Initialize(hsize, lb, ub);
184  }
185 
186  // need copy constructor because of smart pointer to histogram
188  {
189  m_Count = l.m_Count;
190  m_Minimum = l.m_Minimum;
191  m_Maximum = l.m_Maximum;
192  m_Mean = l.m_Mean;
193  m_Sum = l.m_Sum;
194  m_SumOfSquares = l.m_SumOfSquares;
195  m_Sigma = l.m_Sigma;
196  m_Variance = l.m_Variance;
197  m_BoundingBox = l.m_BoundingBox;
198  m_Histogram = l.m_Histogram;
199  }
200 
201  // added for completeness
202  LabelStatistics &operator= (const LabelStatistics& l)
203  {
204  if(this != &l)
205  {
206  m_Count = l.m_Count;
207  m_Minimum = l.m_Minimum;
208  m_Maximum = l.m_Maximum;
209  m_Mean = l.m_Mean;
210  m_Sum = l.m_Sum;
211  m_SumOfSquares = l.m_SumOfSquares;
212  m_Sigma = l.m_Sigma;
213  m_Variance = l.m_Variance;
214  m_BoundingBox = l.m_BoundingBox;
215  m_Histogram = l.m_Histogram;
216  }
217  return *this;
218  }
219 
230  };
231 
233  typedef itksys::hash_map< LabelPixelType, LabelStatistics > MapType;
234  typedef typename itksys::hash_map< LabelPixelType, LabelStatistics >::iterator MapIterator;
235  typedef typename itksys::hash_map< LabelPixelType, LabelStatistics >::const_iterator MapConstIterator;
237 
239  typedef std::vector<LabelPixelType> ValidLabelValuesContainerType;
240 
241  // macros for Histogram enables
242  itkSetMacro(UseHistograms, bool);
243  itkGetConstMacro(UseHistograms, bool);
244  itkBooleanMacro(UseHistograms);
245 
246 
248  {
249  return m_ValidLabelValues;
250  }
251 
253  void SetLabelInput(const TLabelImage *input)
254  {
255  // Process object is not const-correct so the const casting is required.
256  this->SetNthInput( 1, const_cast< TLabelImage * >( input ) );
257  }
258 
261  {
262  return itkDynamicCastInDebugMode< LabelImageType * >( const_cast< DataObject * >( this->ProcessObject::GetInput(1) ) );
263  }
264 
267  bool HasLabel(LabelPixelType label) const
268  {
269  return m_LabelStatistics.find(label) != m_LabelStatistics.end();
270  }
271 
274  {
275  return static_cast<MapSizeType>( m_LabelStatistics.size() );
276  }
277 
279  {
280  return static_cast<MapSizeType>( this->GetNumberOfObjects() );
281  }
282 
284  RealType GetMinimum(LabelPixelType label) const;
285 
287  RealType GetMaximum(LabelPixelType label) const;
288 
290  RealType GetMean(LabelPixelType label) const;
291 
294  RealType GetMedian(LabelPixelType label) const;
295 
297  RealType GetSigma(LabelPixelType label) const;
298 
300  RealType GetVariance(LabelPixelType label) const;
301 
305  BoundingBoxType GetBoundingBox(LabelPixelType label) const;
306 
308  RegionType GetRegion(LabelPixelType label) const;
309 
311  RealType GetSum(LabelPixelType label) const;
312 
314  MapSizeType GetCount(LabelPixelType label) const;
315 
317  HistogramPointer GetHistogram(LabelPixelType label) const;
318 
320  void SetHistogramParameters(const int numBins, RealType lowerBound,
321  RealType upperBound);
322 
323 #ifdef ITK_USE_CONCEPT_CHECKING
324  // Begin concept checking
325  itkConceptMacro( InputHasNumericTraitsCheck,
327  // End concept checking
328 #endif
329 
330 protected:
332  ~LabelStatisticsImageFilter() ITK_OVERRIDE {}
333  void PrintSelf(std::ostream & os, Indent indent) const ITK_OVERRIDE;
334 
337  void AllocateOutputs() ITK_OVERRIDE;
338 
340  void BeforeThreadedGenerateData() ITK_OVERRIDE;
341 
344  void AfterThreadedGenerateData() ITK_OVERRIDE;
345 
347  void ThreadedGenerateData(const RegionType &
348  outputRegionForThread,
349  ThreadIdType threadId) ITK_OVERRIDE;
350 
351  // Override since the filter produces all of its output
352  void EnlargeOutputRequestedRegion(DataObject *data) ITK_OVERRIDE;
353 
354 private:
355  ITK_DISALLOW_COPY_AND_ASSIGN(LabelStatisticsImageFilter);
356 
357  std::vector< MapType > m_LabelStatisticsPerThread;
358  MapType m_LabelStatistics;
359  ValidLabelValuesContainerType m_ValidLabelValues;
360 
361  bool m_UseHistograms;
362 
363  typename HistogramType::SizeType m_NumBins;
364 
365  RealType m_LowerBound;
366  RealType m_UpperBound;
368 }; // end of class
369 } // end namespace itk
370 
371 #ifndef ITK_MANUAL_INSTANTIATION
372 #include "itkLabelStatisticsImageFilter.hxx"
373 #endif
374 
375 #endif
Critical section locking class that can be allocated on the stack.
itksys::hash_map< LabelPixelType, LabelStatistics >::iterator MapIterator
bool HasLabel(LabelPixelType label) const
void SetLabelInput(const TLabelImage *input)
Light weight base class for most itk classes.
itksys::hash_map< LabelPixelType, LabelStatistics > MapType
virtual const ValidLabelValuesContainerType & GetValidLabelValues() const
std::vector< IndexValueType > BoundingBoxType
NumericTraits< PixelType >::RealType RealType
SimpleDataObjectDecorator< RealType > RealObjectType
This class stores measurement vectors in the context of n-dimensional histogram.
Definition: itkHistogram.h:77
Given an intensity image and a label map, compute min, max, variance and mean of the pixels associate...
static ITK_CONSTEXPR_FUNC T max(const T &)
SizeValueType IdentifierType
Definition: itkIntTypes.h:147
void SetSize(SizeValueType sz)
Decorates any &quot;simple&quot; data type (data types without smart pointers) with a DataObject API...
Superclass::MeasurementVectorType MeasurementVectorType
Definition: itkHistogram.h:101
itk::Statistics::Histogram< RealType > HistogramType
std::vector< LabelPixelType > ValidLabelValuesContainerType
unsigned int ThreadIdType
Definition: itkIntTypes.h:159
DataObject * GetInput(const DataObjectIdentifierType &key)
Return an input.
itksys::hash_map< LabelPixelType, LabelStatistics >::const_iterator MapConstIterator
static ITK_CONSTEXPR_FUNC T NonpositiveMin()
Base class for filters that take an image as input and produce an image as output.
LabelStatistics(int size, RealType lowerBound, RealType upperBound)
Control indentation during Print() invocation.
Definition: itkIndent.h:49
const LabelImageType * GetLabelInput() const
Define additional traits for native types such as int or float.
#define itkConceptMacro(name, concept)
Base class for all data objects in ITK.
ImageToImageFilter< TInputImage, TInputImage > Superclass