ITK  4.13.0
Insight Segmentation and Registration Toolkit
itkStatisticsImageFilter.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 itkStatisticsImageFilter_h
19 #define itkStatisticsImageFilter_h
20 
21 #include "itkImageToImageFilter.h"
22 #include "itkNumericTraits.h"
23 #include "itkArray.h"
25 
26 namespace itk
27 {
48 template< typename TInputImage >
49 class ITK_TEMPLATE_EXPORT StatisticsImageFilter:
50  public ImageToImageFilter< TInputImage, TInputImage >
51 {
52 public:
58 
60  itkNewMacro(Self);
61 
64 
66  typedef typename TInputImage::Pointer InputImagePointer;
67 
68  typedef typename TInputImage::RegionType RegionType;
69  typedef typename TInputImage::SizeType SizeType;
71  typedef typename TInputImage::PixelType PixelType;
72 
74  itkStaticConstMacro(ImageDimension, unsigned int,
75  TInputImage::ImageDimension);
76 
79 
82 
86 
89  { return this->GetMinimumOutput()->Get(); }
90  PixelObjectType * GetMinimumOutput();
92 
93  const PixelObjectType * GetMinimumOutput() const;
94 
97  { return this->GetMaximumOutput()->Get(); }
98  PixelObjectType * GetMaximumOutput();
100 
101  const PixelObjectType * GetMaximumOutput() const;
102 
105  { return this->GetMeanOutput()->Get(); }
106  RealObjectType * GetMeanOutput();
108 
109  const RealObjectType * GetMeanOutput() const;
110 
113  { return this->GetSigmaOutput()->Get(); }
114  RealObjectType * GetSigmaOutput();
116 
117  const RealObjectType * GetSigmaOutput() const;
118 
121  { return this->GetVarianceOutput()->Get(); }
122  RealObjectType * GetVarianceOutput();
124 
125  const RealObjectType * GetVarianceOutput() const;
126 
128  RealType GetSum() const
129  { return this->GetSumOutput()->Get(); }
130  RealObjectType * GetSumOutput();
132 
133  const RealObjectType * GetSumOutput() const;
134 
138  using Superclass::MakeOutput;
139  virtual DataObjectPointer MakeOutput(DataObjectPointerArraySizeType idx) ITK_OVERRIDE;
140 
141 #ifdef ITK_USE_CONCEPT_CHECKING
142  // Begin concept checking
143  itkConceptMacro( InputHasNumericTraitsCheck,
145  // End concept checking
146 #endif
147 
148 protected:
150  ~StatisticsImageFilter() ITK_OVERRIDE {}
151  void PrintSelf(std::ostream & os, Indent indent) const ITK_OVERRIDE;
152 
156  void AllocateOutputs() ITK_OVERRIDE;
157 
159  void BeforeThreadedGenerateData() ITK_OVERRIDE;
160 
163  void AfterThreadedGenerateData() ITK_OVERRIDE;
164 
166  void ThreadedGenerateData(const RegionType &
167  outputRegionForThread,
168  ThreadIdType threadId) ITK_OVERRIDE;
169 
170  // Override since the filter needs all the data for the algorithm
171  void GenerateInputRequestedRegion() ITK_OVERRIDE;
172 
173  // Override since the filter produces all of its output
174  void EnlargeOutputRequestedRegion(DataObject *data) ITK_OVERRIDE;
175 
176 private:
177  ITK_DISALLOW_COPY_AND_ASSIGN(StatisticsImageFilter);
178 
179  Array< RealType > m_ThreadSum;
180  Array< RealType > m_SumOfSquares;
181  Array< SizeValueType > m_Count;
182  Array< PixelType > m_ThreadMin;
183  Array< PixelType > m_ThreadMax;
184 }; // end of class
185 } // end namespace itk
186 
187 #ifndef ITK_MANUAL_INSTANTIATION
188 #include "itkStatisticsImageFilter.hxx"
189 #endif
190 
191 #endif
Array class with size defined at construction time.
Definition: itkArray.h:50
TInputImage::Pointer InputImagePointer
Light weight base class for most itk classes.
NumericTraits< PixelType >::RealType RealType
unsigned long SizeValueType
Definition: itkIntTypes.h:143
TInputImage::PixelType PixelType
SimpleDataObjectDecorator< RealType > RealObjectType
Decorates any &quot;simple&quot; data type (data types without smart pointers) with a DataObject API...
DataObjectPointerArray::size_type DataObjectPointerArraySizeType
Compute min. max, variance and mean of an Image.
TInputImage::IndexType IndexType
unsigned int ThreadIdType
Definition: itkIntTypes.h:159
ProcessObject::DataObjectPointerArraySizeType DataObjectPointerArraySizeType
SmartPointer< const Self > ConstPointer
TInputImage::RegionType RegionType
SimpleDataObjectDecorator< PixelType > PixelObjectType
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
Define additional traits for native types such as int or float.
ImageToImageFilter< TInputImage, TInputImage > Superclass
#define itkConceptMacro(name, concept)
Base class for all data objects in ITK.