ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
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{
53template <typename TInputImage>
54class ITK_TEMPLATE_EXPORT StatisticsImageFilter : public ImageSink<TInputImage>
55{
56public:
57 ITK_DISALLOW_COPY_AND_MOVE(StatisticsImageFilter);
58
64
66 itkNewMacro(Self);
67
69 itkOverrideGetNameOfClassMacro(StatisticsImageFilter);
70
72 using InputImagePointer = typename TInputImage::Pointer;
73
74 using RegionType = typename TInputImage::RegionType;
75 using SizeType = typename TInputImage::SizeType;
76 using IndexType = typename TInputImage::IndexType;
77 using PixelType = typename TInputImage::PixelType;
78
80 static constexpr unsigned int ImageDimension = TInputImage::ImageDimension;
81
84
87
91
93 itkGetDecoratedOutputMacro(Minimum, PixelType);
94
96 itkGetDecoratedOutputMacro(Maximum, PixelType);
97
99 itkGetDecoratedOutputMacro(Mean, RealType);
100
102 itkGetDecoratedOutputMacro(Sigma, RealType);
103
105 itkGetDecoratedOutputMacro(Variance, RealType);
106
108 itkGetDecoratedOutputMacro(Sum, RealType);
109
111 itkGetDecoratedOutputMacro(SumOfSquares, RealType);
112
113 // Change the access from protected to public to expose streaming option, a using statement can not be used due to
114 // limitations of wrapping.
115 void
116 SetNumberOfStreamDivisions(const unsigned int n) override
117 {
119 }
120 [[nodiscard]] unsigned int
122 {
124 }
125
131 MakeOutput(const DataObjectIdentifierType & name) override;
132
133 itkConceptMacro(InputHasNumericTraitsCheck, (Concept::HasNumericTraits<PixelType>));
134
135protected:
137 ~StatisticsImageFilter() override = default;
138 void
139 PrintSelf(std::ostream & os, Indent indent) const override;
140
142 void
144
147 void
149
150 void
152
153 itkSetDecoratedOutputMacro(Minimum, PixelType);
154 itkSetDecoratedOutputMacro(Maximum, PixelType);
155 itkSetDecoratedOutputMacro(Mean, RealType);
156 itkSetDecoratedOutputMacro(Sigma, RealType);
157 itkSetDecoratedOutputMacro(Variance, RealType);
158 itkSetDecoratedOutputMacro(Sum, RealType);
159 itkSetDecoratedOutputMacro(SumOfSquares, RealType);
160
161private:
164
168
169 std::mutex m_Mutex{};
170}; // end of class
171} // end namespace itk
172
173#ifndef ITK_MANUAL_INSTANTIATION
174# include "itkStatisticsImageFilter.hxx"
175#endif
176
177#endif
Perform more precise accumulation of floating point numbers.
SmartPointer< Self > Pointer
virtual void SetNumberOfStreamDivisions(unsigned int _arg)
virtual unsigned int GetNumberOfStreamDivisions() const
Control indentation during Print() invocation.
Definition itkIndent.h:50
DataObject::DataObjectIdentifierType DataObjectIdentifierType
virtual DataObjectPointer MakeOutput(DataObjectPointerArraySizeType idx)
Make a DataObject of the correct type to used as the specified output.
Decorates any "simple" data type (data types without smart pointers) with a DataObject API.
Implements transparent reference counting.
SimpleDataObjectDecorator< RealType > RealObjectType
static constexpr unsigned int ImageDimension
SimpleDataObjectDecorator< PixelType > PixelObjectType
CompensatedSummation< RealType > m_ThreadSum
typename NumericTraits< PixelType >::RealType RealType
ProcessObject::DataObjectIdentifierType DataObjectIdentifierType
typename TInputImage::RegionType RegionType
CompensatedSummation< RealType > m_SumOfSquares
ImageSink< TInputImage > Superclass
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::Pointer InputImagePointer
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
SmartPointer< const Self > ConstPointer
void SetNumberOfStreamDivisions(const unsigned int n) override
unsigned int GetNumberOfStreamDivisions() const override
void BeforeStreamedGenerateData() override
#define itkConceptMacro(name, concept)
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
unsigned long SizeValueType
Definition itkIntTypes.h:86