ITK  6.0.0
Insight Toolkit
itkStatisticsImageFilter.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 itkStatisticsImageFilter_h
19#define itkStatisticsImageFilter_h
20
21#include "itkImageSink.h"
22#include "itkNumericTraits.h"
23#include "itkArray.h"
25#include <mutex>
27
28namespace itk
29{
54template <typename TInputImage>
55class ITK_TEMPLATE_EXPORT StatisticsImageFilter : public ImageSink<TInputImage>
56{
57public:
58 ITK_DISALLOW_COPY_AND_MOVE(StatisticsImageFilter);
59
65
67 itkNewMacro(Self);
68
70 itkOverrideGetNameOfClassMacro(StatisticsImageFilter);
71
74
78 using PixelType = typename TInputImage::PixelType;
79
81 static constexpr unsigned int ImageDimension = TInputImage::ImageDimension;
82
85
88
92
94 itkGetDecoratedOutputMacro(Minimum, PixelType);
95
97 itkGetDecoratedOutputMacro(Maximum, PixelType);
98
100 itkGetDecoratedOutputMacro(Mean, RealType);
101
103 itkGetDecoratedOutputMacro(Sigma, RealType);
104
106 itkGetDecoratedOutputMacro(Variance, RealType);
107
109 itkGetDecoratedOutputMacro(Sum, RealType);
110
112 itkGetDecoratedOutputMacro(SumOfSquares, RealType);
113
114 // Change the access from protected to public to expose streaming option, a using statement can not be used due to
115 // limitations of wrapping.
116 void
117 SetNumberOfStreamDivisions(const unsigned int n) override
118 {
119 Superclass::SetNumberOfStreamDivisions(n);
120 }
121 unsigned int
123 {
124 return Superclass::GetNumberOfStreamDivisions();
125 }
126
130 using Superclass::MakeOutput;
132 MakeOutput(const DataObjectIdentifierType & name) override;
133
134#ifdef ITK_USE_CONCEPT_CHECKING
135 // Begin concept checking
136 itkConceptMacro(InputHasNumericTraitsCheck, (Concept::HasNumericTraits<PixelType>));
137 // End concept checking
138#endif
139
140protected:
142 ~StatisticsImageFilter() override = default;
143 void
144 PrintSelf(std::ostream & os, Indent indent) const override;
145
147 void
149
152 void
154
155 void
157
158 itkSetDecoratedOutputMacro(Minimum, PixelType);
159 itkSetDecoratedOutputMacro(Maximum, PixelType);
160 itkSetDecoratedOutputMacro(Mean, RealType);
161 itkSetDecoratedOutputMacro(Sigma, RealType);
162 itkSetDecoratedOutputMacro(Variance, RealType);
163 itkSetDecoratedOutputMacro(Sum, RealType);
164 itkSetDecoratedOutputMacro(SumOfSquares, RealType);
165
166private:
169
170 SizeValueType m_Count{ 1 };
171 PixelType m_ThreadMin{ 1 };
172 PixelType m_ThreadMax{ 1 };
173
174 std::mutex m_Mutex{};
175}; // end of class
176} // end namespace itk
177
178#ifndef ITK_MANUAL_INSTANTIATION
179# include "itkStatisticsImageFilter.hxx"
180#endif
181
182#endif
SmartPointer< Self > Pointer
typename InputImageType::Pointer InputImagePointer
Definition: itkImageSink.h:74
Control indentation during Print() invocation.
Definition: itkIndent.h:50
Define additional traits for native types such as int or float.
DataObject::DataObjectIdentifierType DataObjectIdentifierType
Decorates any "simple" data type (data types without smart pointers) with a DataObject API.
Compute min, max, variance and mean of an Image.
typename NumericTraits< PixelType >::RealType RealType
typename TInputImage::RegionType RegionType
void AfterStreamedGenerateData() override
void PrintSelf(std::ostream &os, Indent indent) const override
typename TInputImage::PixelType PixelType
void ThreadedStreamedGenerateData(const RegionType &) override
~StatisticsImageFilter() override=default
typename TInputImage::SizeType SizeType
DataObjectPointer MakeOutput(const DataObjectIdentifierType &name) override
Make a DataObject of the correct type to used as the specified output.
typename TInputImage::IndexType IndexType
void SetNumberOfStreamDivisions(const unsigned int n) override
unsigned int GetNumberOfStreamDivisions() const override
void BeforeStreamedGenerateData() override
Base class interface to process data on multiple requested input chunks.
SmartPointer< Self > Pointer
#define itkConceptMacro(name, concept)
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
unsigned long SizeValueType
Definition: itkIntTypes.h:86