ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
itkAdaptiveHistogramEqualizationImageFilter.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 itkAdaptiveHistogramEqualizationImageFilter_h
19#define itkAdaptiveHistogramEqualizationImageFilter_h
20
23#include "itkImage.h"
24
25namespace itk
26{
69template <typename TImageType, typename TKernel = Neighborhood<bool, TImageType::ImageDimension>>
72 TImageType,
73 TImageType,
74 TKernel,
75 typename Function::AdaptiveEqualizationHistogram<typename TImageType::PixelType, typename TImageType::PixelType>>
76
77{
78public:
79 ITK_DISALLOW_COPY_AND_MOVE(AdaptiveHistogramEqualizationImageFilter);
80
86 TImageType,
87 TImageType,
88 TKernel,
92
93 static constexpr unsigned int ImageDimension = TImageType::ImageDimension;
94
96 itkNewMacro(Self);
97
99 itkOverrideGetNameOfClassMacro(AdaptiveHistogramEqualizationImageFilter);
100
102 using ImageType = TImageType;
103 using InputPixelType = typename ImageType::PixelType;
104 using ImageSizeType = typename ImageType::SizeType;
105
109 itkSetMacro(Alpha, float);
110 itkGetConstMacro(Alpha, float);
112
117 itkSetMacro(Beta, float);
118 itkGetConstMacro(Beta, float);
120
121#if !defined(ITK_FUTURE_LEGACY_REMOVE)
126 virtual void
127 SetUseLookupTable(const bool _arg)
128 {
129 itkDebugMacro("setting UseLookupTable to " << _arg);
130 itkGenericLegacyReplaceBodyMacro("UseLookupTable", "", "nothing");
131 if (this->m_UseLookupTable != _arg)
132 {
133 this->m_UseLookupTable = _arg;
134 this->Modified();
135 }
136 }
138
139 itkGetConstMacro(UseLookupTable, bool);
140 itkBooleanMacro(UseLookupTable);
141#endif
142
143 void
145 {
146 h.SetAlpha(this->m_Alpha);
147 h.SetBeta(this->m_Beta);
148 h.SetMinimum(this->m_InputMinimum);
149 h.SetMaximum(this->m_InputMaximum);
150
151 typename Superclass::HistogramType::RealType kernelSize = 1;
152 for (unsigned int i = 0; i < ImageDimension; ++i)
153 {
154 kernelSize *= (2 * this->GetRadius()[i] + 1);
155 }
156 h.SetKernelSize(kernelSize);
157 }
158
159protected:
172
174 void
175 PrintSelf(std::ostream & os, Indent indent) const override;
176
180 void
182
183private:
184 float m_Alpha{};
185 float m_Beta{};
186
189
191};
192} // end namespace itk
193
194#ifndef ITK_MANUAL_INSTANTIATION
195# include "itkAdaptiveHistogramEqualizationImageFilter.hxx"
196#endif
197
198#endif
void ConfigureHistogram(typename Superclass::HistogramType &h) override
~AdaptiveHistogramEqualizationImageFilter() override=default
void PrintSelf(std::ostream &os, Indent indent) const override
MovingHistogramImageFilter< TImageType, TImageType, TKernel, typename Function::AdaptiveEqualizationHistogram< typename TImageType::PixelType, typename TImageType::PixelType > > Superclass
virtual const RadiusType & GetRadius() const
Control indentation during Print() invocation.
Definition itkIndent.h:50
void SetRadius(const RadiusType &radius) override
static constexpr T max(const T &)
static constexpr T min(const T &)
virtual void Modified() const
Implements transparent reference counting.
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....