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{
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
73 using InputImagePointer = typename TInputImage::Pointer;
74
75 using RegionType = typename TInputImage::RegionType;
76 using SizeType = typename TInputImage::SizeType;
77 using IndexType = typename TInputImage::IndexType;
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 {
120 }
121 unsigned int
123 {
125 }
126
132 MakeOutput(const DataObjectIdentifierType & name) override;
133
134 itkConceptMacro(InputHasNumericTraitsCheck, (Concept::HasNumericTraits<PixelType>));
135
136protected:
138 ~StatisticsImageFilter() override = default;
139 void
140 PrintSelf(std::ostream & os, Indent indent) const override;
141
143 void
145
148 void
150
151 void
153
154 itkSetDecoratedOutputMacro(Minimum, PixelType);
155 itkSetDecoratedOutputMacro(Maximum, PixelType);
156 itkSetDecoratedOutputMacro(Mean, RealType);
157 itkSetDecoratedOutputMacro(Sigma, RealType);
158 itkSetDecoratedOutputMacro(Variance, RealType);
159 itkSetDecoratedOutputMacro(Sum, RealType);
160 itkSetDecoratedOutputMacro(SumOfSquares, RealType);
161
162private:
165
169
170 std::mutex m_Mutex{};
171}; // end of class
172} // end namespace itk
173
174#ifndef ITK_MANUAL_INSTANTIATION
175# include "itkStatisticsImageFilter.hxx"
176#endif
177
178#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
typename DataObject::Pointer DataObjectPointer
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