ITK  4.13.0
Insight Segmentation and Registration Toolkit
itkHistogram.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 itkHistogram_h
19 #define itkHistogram_h
20 
21 #include <vector>
22 
23 #include "itkArray.h"
24 #include "itkSample.h"
27 
28 namespace itk
29 {
30 namespace Statistics
31 {
32 
75 template< typename TMeasurement = float,
76  typename TFrequencyContainer = DenseFrequencyContainer2 >
77 class ITK_TEMPLATE_EXPORT Histogram:
78  public Sample< Array< TMeasurement > >
79 {
80 public:
81 
82  // This type serves as the indirect definition of MeasurementVectorType
84 
86  typedef Histogram Self;
90 
92  itkTypeMacro(Histogram, Sample);
93 
95  itkNewMacro(Self);
96 
98  typedef TMeasurement MeasurementType;
99 
102  typedef typename Superclass::InstanceIdentifier InstanceIdentifier;
103  typedef typename Superclass::MeasurementVectorSizeType MeasurementVectorSizeType;
104 
106 
108  typedef TFrequencyContainer FrequencyContainerType;
109  typedef typename FrequencyContainerType::Pointer FrequencyContainerPointer;
110 
112  typedef typename FrequencyContainerType::AbsoluteFrequencyType AbsoluteFrequencyType;
113  typedef typename FrequencyContainerType::TotalAbsoluteFrequencyType TotalAbsoluteFrequencyType;
114  typedef typename FrequencyContainerType::RelativeFrequencyType RelativeFrequencyType;
115  typedef typename FrequencyContainerType::TotalRelativeFrequencyType TotalRelativeFrequencyType;
116 
120 
124 
126  typedef std::vector< MeasurementType > BinMinVectorType;
127  typedef std::vector< MeasurementType > BinMaxVectorType;
128  typedef std::vector< BinMinVectorType > BinMinContainerType;
129  typedef std::vector< BinMaxVectorType > BinMaxContainerType;
130 
134  void Initialize(const SizeType & size);
135 
139  void Initialize(const SizeType & size, MeasurementVectorType & lowerBound,
140  MeasurementVectorType & upperBound);
141 
143  void SetToZero();
144 
148  itkLegacyMacro(const IndexType & GetIndex(const MeasurementVectorType & measurement) const);
149 
153  bool GetIndex(const MeasurementVectorType & measurement,
154  IndexType & index) const;
155 
159  const IndexType & GetIndex(InstanceIdentifier id) const;
160 
163  itkGetConstMacro(ClipBinsAtEnds, bool);
164 
167  itkSetMacro(ClipBinsAtEnds, bool);
168 
171  bool IsIndexOutOfBounds(const IndexType & index) const;
172 
176  InstanceIdentifier GetInstanceIdentifier(const IndexType & index) const;
177 
179  InstanceIdentifier Size() const ITK_OVERRIDE;
180 
182  const SizeType & GetSize() const;
183 
185  SizeValueType GetSize(unsigned int dimension) const;
186 
188  const MeasurementType & GetBinMin(unsigned int dimension,
189  InstanceIdentifier nbin) const;
190 
192  const MeasurementType & GetBinMax(unsigned int dimension,
193  InstanceIdentifier nbin) const;
194 
196  void SetBinMin(unsigned int dimension, InstanceIdentifier nbin,
197  MeasurementType min);
198 
200  void SetBinMax(unsigned int dimension,
202 
205  const MeasurementType & GetBinMinFromValue(unsigned int dimension,
206  float value) const;
207 
210  const MeasurementType & GetBinMaxFromValue(unsigned int dimension,
211  float value) const;
212 
214  const BinMinVectorType & GetDimensionMins(unsigned int dimension) const;
215 
217  const BinMaxVectorType & GetDimensionMaxs(unsigned int dimension) const;
218 
220  const BinMinContainerType & GetMins() const;
221 
223  const BinMaxContainerType & GetMaxs() const;
224 
226  const MeasurementVectorType & GetHistogramMinFromIndex(const IndexType & index) const;
227 
229  const MeasurementVectorType & GetHistogramMaxFromIndex(const IndexType & index) const;
230 
232  AbsoluteFrequencyType GetFrequency(InstanceIdentifier id) const ITK_OVERRIDE;
233 
235  AbsoluteFrequencyType GetFrequency(const IndexType & index) const;
236 
238  void SetFrequency(AbsoluteFrequencyType value);
239 
242  bool SetFrequency(InstanceIdentifier id, AbsoluteFrequencyType value);
243 
246  bool SetFrequencyOfIndex(const IndexType & index,
247  AbsoluteFrequencyType value);
248 
251  bool SetFrequencyOfMeasurement(const MeasurementVectorType & measurement,
252  AbsoluteFrequencyType value);
253 
257  bool IncreaseFrequency(InstanceIdentifier id, AbsoluteFrequencyType value);
258 
262  bool IncreaseFrequencyOfIndex(const IndexType & index,
263  AbsoluteFrequencyType value);
264 
272  bool IncreaseFrequencyOfMeasurement(
273  const MeasurementVectorType & measurement,
274  AbsoluteFrequencyType value);
275 #ifdef ITKV3_COMPATIBILITY
276  //In ITKv4 the member functions are given unique names to dis-ambiguate
277  //the intended behavior. Make aliases of the overloaded calls
278  //for ITKv3 backwards compatibility.
279  bool IncreaseFrequency(const IndexType & index,
280  AbsoluteFrequencyType value)
281  {
282  return IncreaseFrequencyOfIndex(index,value);
283  }
285 
286  bool IncreaseFrequency(
287  const MeasurementVectorType & measurement,
288  AbsoluteFrequencyType value)
289  {
290  return IncreaseFrequencyOfMeasurement(measurement,value);
291  }
292 
293 #endif
294 
298  const MeasurementVectorType & GetMeasurementVector(InstanceIdentifier id) const ITK_OVERRIDE;
299 
301  const MeasurementVectorType & GetMeasurementVector(const IndexType & index) const;
302 
305  MeasurementType GetMeasurement(InstanceIdentifier n,
306  unsigned int dimension) const;
307 
309  TotalAbsoluteFrequencyType GetTotalFrequency() const ITK_OVERRIDE;
310 
312  AbsoluteFrequencyType GetFrequency(InstanceIdentifier n,
313  unsigned int dimension) const;
314 
330  double Quantile(unsigned int dimension, double p) const;
332 
334  double Mean(unsigned int dimension) const;
335 
337  virtual void Graft(const DataObject *) ITK_OVERRIDE;
338 
339 protected:
340  void PrintSelf(std::ostream & os, Indent indent) const ITK_OVERRIDE;
341 
342 public:
343 
349  {
350 public:
351 
352  friend class Histogram;
353 
354  ConstIterator(const Self *histogram)
355  {
356  m_Id = 0;
357  m_Histogram = histogram;
358  }
359 
360  ConstIterator(const ConstIterator & it)
361  {
362  m_Id = it.m_Id;
363  m_Histogram = it.m_Histogram;
364  }
365 
366  ConstIterator & operator=(const ConstIterator & it)
367  {
368  m_Id = it.m_Id;
369  m_Histogram = it.m_Histogram;
370  return *this;
371  }
372 
374  {
375  return m_Histogram->GetFrequency(m_Id);
376  }
377 
379  {
380  return m_Id;
381  }
382 
384  {
385  return m_Histogram->GetMeasurementVector(m_Id);
386  }
387 
388  const IndexType & GetIndex() const
389  {
390  return m_Histogram->GetIndex(m_Id);
391  }
392 
393  ConstIterator & operator++()
394  {
395  ++m_Id;
396  return *this;
397  }
398 
399  bool operator!=(const ConstIterator & it)
400  {
401  return ( m_Id != it.m_Id );
402  }
403 
404  bool operator==(const ConstIterator & it)
405  {
406  return ( m_Id == it.m_Id );
407  }
408 
409 protected:
410 
411  ConstIterator(InstanceIdentifier id, const Self *histogram):
412  m_Id(id), m_Histogram(histogram)
413  {}
414 
415  // ConstIterator pointing DenseFrequencyContainer
417 
418  // Pointer of DenseFrequencyContainer
420 private:
421  ConstIterator() ITK_DELETED_FUNCTION;
422  }; // end of iterator class
423 
428  class Iterator:public ConstIterator
429  {
430 public:
431 
432  Iterator(Self *histogram):ConstIterator(histogram)
433  {}
434 
436  ConstIterator(id, histogram)
437  {}
438 
439  Iterator(const Iterator & it):ConstIterator(it)
440  {}
441 
443  {
444  this->ConstIterator::operator=(it);
445  return *this;
446  }
447 
449  {
450  Self *histogram = const_cast< Self * >( this->m_Histogram );
451 
452  return histogram->SetFrequency(this->m_Id, value);
453  }
454 
455 private:
456  // To ensure const-correctness these method must not be in the public API.
457  // The are not implemented, since they should never be called.
458  Iterator() ITK_DELETED_FUNCTION;
459  Iterator(const Self *histogram) ITK_DELETED_FUNCTION;
460  Iterator(InstanceIdentifier id, const Self *histogram) ITK_DELETED_FUNCTION;
461  Iterator(const ConstIterator & it) ITK_DELETED_FUNCTION;
462  ConstIterator & operator=(const ConstIterator & it) ITK_DELETED_FUNCTION;
463  }; // end of iterator class
464 
465  Iterator Begin()
466  {
467  Iterator iter(0, this);
468 
469  return iter;
470  }
471 
473  {
474  return Iterator(m_OffsetTable[this->GetMeasurementVectorSize()], this);
475  }
476 
477  ConstIterator Begin() const
478  {
479  ConstIterator iter(0, this);
480 
481  return iter;
482  }
483 
484  ConstIterator End() const
485  {
486  return ConstIterator(m_OffsetTable[this->GetMeasurementVectorSize()], this);
487  }
488 
489 protected:
490  Histogram();
491  virtual ~Histogram() ITK_OVERRIDE {}
492 
493  // The number of bins for each dimension
495 
496 private:
497  ITK_DISALLOW_COPY_AND_ASSIGN(Histogram);
498 
499  typedef std::vector< InstanceIdentifier > OffsetTableType;
500  OffsetTableType m_OffsetTable;
502  unsigned int m_NumberOfInstances;
503 
504  // This method is provided here just to avoid a "hidden" warning
505  // related to the virtual method available in DataObject.
506  virtual void Initialize() ITK_OVERRIDE {}
507 
508  // lower bound of each bin
509  std::vector< std::vector< MeasurementType > > m_Min;
510 
511  // upper bound of each bin
512  std::vector< std::vector< MeasurementType > > m_Max;
513 
516 
518 };
519 } // end of namespace Statistics
520 } // end of namespace itk
521 
522 #ifndef ITK_MANUAL_INSTANTIATION
523 #include "itkHistogram.hxx"
524 #endif
525 
526 #endif
Sample< ArrayType > Superclass
Definition: itkHistogram.h:87
Array class with size defined at construction time.
Definition: itkArray.h:50
IndexType::ValueType IndexValueType
Definition: itkHistogram.h:119
::itk::IndexValueType ValueType
Definition: itkArray.h:55
ConstIterator(InstanceIdentifier id, const Self *histogram)
Definition: itkHistogram.h:411
SmartPointer< Self > Pointer
Definition: itkHistogram.h:88
ConstIterator End() const
Definition: itkHistogram.h:484
ConstIterator(const ConstIterator &it)
Definition: itkHistogram.h:360
Represent the size (bounds) of a n-dimensional image.
Definition: itkSize.h:52
class that walks through the elements of the histogram.
Definition: itkHistogram.h:348
Superclass::InstanceIdentifier InstanceIdentifier
Definition: itkHistogram.h:102
virtual ~Histogram() override
Definition: itkHistogram.h:491
std::vector< BinMinVectorType > BinMinContainerType
Definition: itkHistogram.h:128
This class stores measurement vectors in the context of n-dimensional histogram.
Definition: itkHistogram.h:77
Superclass::MeasurementVectorSizeType MeasurementVectorSizeType
Definition: itkHistogram.h:103
FrequencyContainerType::AbsoluteFrequencyType AbsoluteFrequencyType
Definition: itkHistogram.h:112
TFrequencyContainer FrequencyContainerType
Definition: itkHistogram.h:108
SmartPointer< const Self > ConstPointer
Definition: itkHistogram.h:89
virtual void Initialize() override
Definition: itkHistogram.h:506
bool SetFrequency(const AbsoluteFrequencyType value)
Definition: itkHistogram.h:448
const MeasurementVectorType & GetMeasurementVector() const
Definition: itkHistogram.h:383
std::vector< MeasurementType > BinMaxVectorType
Definition: itkHistogram.h:127
Array< ::itk::IndexValueType > IndexType
Definition: itkHistogram.h:118
std::vector< std::vector< MeasurementType > > m_Max
Definition: itkHistogram.h:512
std::vector< std::vector< MeasurementType > > m_Min
Definition: itkHistogram.h:509
const IndexType & GetIndex() const
Definition: itkHistogram.h:388
OffsetTableType m_OffsetTable
Definition: itkHistogram.h:500
ConstIterator & operator=(const ConstIterator &it)
Definition: itkHistogram.h:366
FrequencyContainerType::RelativeFrequencyType RelativeFrequencyType
Definition: itkHistogram.h:114
Iterator(InstanceIdentifier id, Self *histogram)
Definition: itkHistogram.h:435
unsigned int m_NumberOfInstances
Definition: itkHistogram.h:502
MeasurementVectorType ValueType
Definition: itkHistogram.h:105
Superclass::MeasurementVectorType MeasurementVectorType
Definition: itkHistogram.h:101
std::vector< BinMaxVectorType > BinMaxContainerType
Definition: itkHistogram.h:129
ConstIterator Begin() const
Definition: itkHistogram.h:477
Iterator & operator=(const Iterator &it)
Definition: itkHistogram.h:442
bool operator==(const ConstIterator &it)
Definition: itkHistogram.h:404
TMeasurement MeasurementType
Definition: itkHistogram.h:95
MeasurementVectorType m_TempMeasurementVector
Definition: itkHistogram.h:514
FrequencyContainerPointer m_FrequencyContainer
Definition: itkHistogram.h:501
TMeasurementVector MeasurementVectorType
Definition: itkSample.h:71
class that walks through the elements of the histogram.
Definition: itkHistogram.h:428
FrequencyContainerType::TotalAbsoluteFrequencyType TotalAbsoluteFrequencyType
Definition: itkHistogram.h:113
A collection of measurements for statistical analysis.
Definition: itkSample.h:61
Control indentation during Print() invocation.
Definition: itkIndent.h:49
std::vector< MeasurementType > BinMinVectorType
Definition: itkHistogram.h:126
FrequencyContainerType::Pointer FrequencyContainerPointer
Definition: itkHistogram.h:109
Array< ::itk::SizeValueType > SizeType
Definition: itkHistogram.h:122
Base class for all data objects in ITK.
bool operator!=(const ConstIterator &it)
Definition: itkHistogram.h:399
AbsoluteFrequencyType GetFrequency() const
Definition: itkHistogram.h:373
Array< TMeasurement > ArrayType
Definition: itkHistogram.h:83
InstanceIdentifier GetInstanceIdentifier() const
Definition: itkHistogram.h:378
SizeType::ValueType SizeValueType
Definition: itkHistogram.h:123
std::vector< InstanceIdentifier > OffsetTableType
Definition: itkHistogram.h:497
FrequencyContainerType::TotalRelativeFrequencyType TotalRelativeFrequencyType
Definition: itkHistogram.h:115