ITK  5.0.0
Insight Segmentation and Registration Toolkit
itkAdaptiveHistogramEqualizationImageFilter.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright Insight Software Consortium
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  * http://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 
25 namespace itk
26 {
70 template< typename TImageType , typename TKernel = Neighborhood<bool, TImageType::ImageDimension> >
71 class ITK_TEMPLATE_EXPORT AdaptiveHistogramEqualizationImageFilter:
72  public MovingHistogramImageFilter< TImageType,
73  TImageType,
74  TKernel,
75  typename Function::AdaptiveEqualizationHistogram< typename TImageType::PixelType,
76  typename TImageType::PixelType > >
77 
78 {
79 public:
80  ITK_DISALLOW_COPY_AND_ASSIGN(AdaptiveHistogramEqualizationImageFilter);
81 
86  using Superclass = MovingHistogramImageFilter< TImageType,
87  TImageType,
88  TKernel,
89  typename Function::AdaptiveEqualizationHistogram< typename TImageType::PixelType,
90  typename TImageType::PixelType > >;
93 
94  static constexpr unsigned int ImageDimension = TImageType::ImageDimension;
95 
97  itkNewMacro(Self);
98 
101 
103  using ImageType = TImageType;
104  using InputPixelType = typename ImageType::PixelType;
106 
110  itkSetMacro(Alpha, float);
111  itkGetConstMacro(Alpha, float);
113 
118  itkSetMacro(Beta, float);
119  itkGetConstMacro(Beta, float);
121 
122 #if ! defined ( ITK_FUTURE_LEGACY_REMOVE )
123 
127  virtual void 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  }
137  itkGetConstMacro(UseLookupTable, bool);
138  itkBooleanMacro(UseLookupTable);
139 #endif
140 
141 
142  void ConfigureHistogram( typename Superclass::HistogramType &h) override
143  {
144  h.SetAlpha( this->m_Alpha );
145  h.SetBeta( this->m_Beta );
146  h.SetMinimum( this->m_InputMinimum );
147  h.SetMaximum( this->m_InputMaximum );
148 
149  typename Superclass::HistogramType::RealType kernelSize = 1;
150  for ( unsigned int i = 0; i < ImageDimension; i++ )
151  {
152  kernelSize *= ( 2 * this->GetRadius()[i] + 1 );
153  }
154  h.SetKernelSize(kernelSize);
155  }
156 
157 protected:
159  {
160  m_Alpha = .3;
161  m_Beta = .3;
162 
163  this->SetRadius(5);
164 
165  m_InputMinimum = NumericTraits< InputPixelType >::min();
166  m_InputMaximum = NumericTraits< InputPixelType >::max();
167 
168  m_UseLookupTable = false;
169  }
170 
171  ~AdaptiveHistogramEqualizationImageFilter() override = default;
172  void PrintSelf(std::ostream & os, Indent indent) const override;
173 
177  void BeforeThreadedGenerateData() override;
178 
179 private:
180  float m_Alpha;
181  float m_Beta;
182 
185 
187 
188 };
189 } // end namespace itk
190 
191 #ifndef ITK_MANUAL_INSTANTIATION
192 #include "itkAdaptiveHistogramEqualizationImageFilter.hxx"
193 #endif
194 
195 #endif
Define numeric traits for std::vector.
Implements a generic moving histogram algorithm.
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
Base class for all process objects that output image data.
void ConfigureHistogram(typename Superclass::HistogramType &h) override
typename TInputImage::PixelType InputPixelType
Base class for filters that take an image as input and produce an image as output.
Control indentation during Print() invocation.
Definition: itkIndent.h:49