ITK  6.0.0
Insight Toolkit
itkMinimumMaximumImageFilter.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 itkMinimumMaximumImageFilter_h
19#define itkMinimumMaximumImageFilter_h
20
21#include "itkImageSink.h"
23#include <mutex>
24
25#include <vector>
26
27#include "itkNumericTraits.h"
28
29namespace itk
30{
47template <typename TInputImage>
48class ITK_TEMPLATE_EXPORT MinimumMaximumImageFilter : public ImageSink<TInputImage>
49{
50public:
51 ITK_DISALLOW_COPY_AND_MOVE(MinimumMaximumImageFilter);
52
54 static constexpr unsigned int InputImageDimension = TInputImage::ImageDimension;
55
61
64
68 using PixelType = typename TInputImage::PixelType;
69
72
74 itkNewMacro(Self);
75
77 itkOverrideGetNameOfClassMacro(MinimumMaximumImageFilter);
78
80 using InputImageType = TInputImage;
81
84
86 itkGetDecoratedOutputMacro(Minimum, PixelType);
87
89 itkGetDecoratedOutputMacro(Maximum, PixelType);
90
94 using Superclass::MakeOutput;
96 MakeOutput(const DataObjectIdentifierType & name) override;
97
98
99 // Change the access from protected to public to expose streaming option, a using statement can not be used due to
100 // limitations of wrapping.
101 void
102 SetNumberOfStreamDivisions(const unsigned int n) override
103 {
104 Superclass::SetNumberOfStreamDivisions(n);
105 }
106 unsigned int
108 {
109 return Superclass::GetNumberOfStreamDivisions();
110 }
111
112#ifdef ITK_USE_CONCEPT_CHECKING
113 // Begin concept checking
115 itkConceptMacro(GreaterThanComparableCheck, (Concept::GreaterThanComparable<PixelType>));
117 // End concept checking
118#endif
119
120protected:
122 ~MinimumMaximumImageFilter() override = default;
123 void
124 PrintSelf(std::ostream & os, Indent indent) const override;
125
127 void
129
132 void
134
135 void
137
138
139 itkSetDecoratedOutputMacro(Minimum, PixelType);
140 itkSetDecoratedOutputMacro(Maximum, PixelType);
141
142private:
143 PixelType m_ThreadMin{};
144 PixelType m_ThreadMax{};
145
146 std::mutex m_Mutex{};
147};
148} // end namespace itk
149
150#ifndef ITK_MANUAL_INSTANTIATION
151# include "itkMinimumMaximumImageFilter.hxx"
152#endif
153
154#endif
SmartPointer< Self > Pointer
typename InputImageType::Pointer InputImagePointer
Definition: itkImageSink.h:74
TInputImage InputImageType
Definition: itkImageSink.h:73
Control indentation during Print() invocation.
Definition: itkIndent.h:50
Computes the minimum and the maximum intensity values of an image.
typename TInputImage::IndexType IndexType
DataObjectPointer MakeOutput(const DataObjectIdentifierType &name) override
Make a DataObject of the correct type to used as the specified output.
typename TInputImage::PixelType PixelType
~MinimumMaximumImageFilter() override=default
unsigned int GetNumberOfStreamDivisions() const override
typename TInputImage::SizeType SizeType
void ThreadedStreamedGenerateData(const RegionType &) override
void SetNumberOfStreamDivisions(const unsigned int n) override
void PrintSelf(std::ostream &os, Indent indent) const override
typename TInputImage::RegionType RegionType
void AfterStreamedGenerateData() override
void BeforeStreamedGenerateData() override
DataObject::DataObjectIdentifierType DataObjectIdentifierType
Decorates any "simple" data type (data types without smart pointers) with a DataObject API.
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....