OTB  10.0.0
Orfeo Toolbox
otbCLHistogramEqualizationFilter.hxx
Go to the documentation of this file.
1 /*
2  * Copyright (C) 2005-2024 Centre National d'Etudes Spatiales (CNES)
3  *
4  * This file is part of Orfeo Toolbox
5  *
6  * https://www.orfeo-toolbox.org/
7  *
8  * Licensed under the Apache License, Version 2.0 (the "License");
9  * you may not use this file except in compliance with the License.
10  * You may obtain a copy of the License at
11  *
12  * http://www.apache.org/licenses/LICENSE-2.0
13  *
14  * Unless required by applicable law or agreed to in writing, software
15  * distributed under the License is distributed on an "AS IS" BASIS,
16  * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
17  * See the License for the specific language governing permissions and
18  * limitations under the License.
19  */
20 
21 #ifndef otbCLHistogramEqualizationFilter_hxx
22 #define otbCLHistogramEqualizationFilter_hxx
23 
25 #include "itkImageRegionIterator.h"
26 
27 #include <limits>
28 
29 namespace otb
30 {
31 template <class TInputImage, class TOutputImage>
33  : m_HistoFilter(HistoFilter::New()),
34  m_GainLutFilter(GainLutFilter::New()),
35  m_ApplyGainFilter(ApplyGainFilter::New()),
36  m_StreamingImageFilter(StreamingImageFilter::New()),
37  m_BufferFilter(BufferFilter::New()),
38  m_Min(std::numeric_limits<InputPixelType>::quiet_NaN()),
39  m_Max(std::numeric_limits<InputPixelType>::quiet_NaN()),
40  m_NoData(std::numeric_limits<InputPixelType>::quiet_NaN()),
41  m_NbBin(256),
42  m_Threshold(-1),
43  m_NoDataFlag(false)
44 {
45  m_ThumbSize.Fill(0);
46  m_GainLutFilter->SetInput(m_HistoFilter->GetHistoOutput());
47  m_StreamingImageFilter->SetInput(m_GainLutFilter->GetOutput());
48  m_ApplyGainFilter->SetInputLut(m_StreamingImageFilter->GetOutput());
49  m_ApplyGainFilter->SetInputImage(m_BufferFilter->GetOutput());
50 }
51 
52 template <class TInputImage, class TOutputImage>
54 {
55  const InputImageType* input = this->GetInput();
56  m_HistoFilter->SetInput(input);
57  m_BufferFilter->SetInput(input);
58  m_ApplyGainFilter->GetOutput()->UpdateOutputInformation();
59  this->GetOutput()->CopyInformation(m_ApplyGainFilter->GetOutput());
60 }
61 
62 template <class TInputImage, class TOutputImage>
64 {
65  m_ApplyGainFilter->GetOutput()->SetRequestedRegion(static_cast<OutputImageType*>(output)->GetRequestedRegion());
66  m_ApplyGainFilter->GetOutput()->PropagateRequestedRegion();
67 }
68 
69 template <class TInputImage, class TOutputImage>
71 {
72  m_ApplyGainFilter->GraftOutput(this->GetOutput());
73  m_ApplyGainFilter->Update();
74  this->GraftOutput(m_ApplyGainFilter->GetOutput());
75 }
76 
80 template <class TInputImage, class TOutputImage>
81 void CLHistogramEqualizationFilter<TInputImage, TOutputImage>::PrintSelf(std::ostream& os, itk::Indent indent) const
82 {
83  Superclass::PrintSelf(os, indent);
84  os << indent << "Minimum : " << m_Min << std::endl;
85  os << indent << "Maximum : " << m_Max << std::endl;
86  os << indent << "Bin Number : " << m_NbBin << std::endl;
87  os << indent << "Thumbnail size : " << m_ThumbSize << std::endl;
88  os << indent << "Threshold value : " << m_Threshold << std::endl;
89  os << indent << "Is no data activated : " << m_NoDataFlag << std::endl;
90  os << indent << "No Data : " << m_NoData << std::endl;
91 }
93 
94 
95 } // End namespace otb
96 
97 #endif
Apply gain on the input image with a bilineare interpolation.
itk::StreamingImageFilter< LutType, LutType > StreamingImageFilter
void PrintSelf(std::ostream &os, itk::Indent indent) const override
StreamingImageFilter::Pointer m_StreamingImageFilter
void PropagateRequestedRegion(itk::DataObject *output) override
Compute the gain for each pixel value from a histogram.
Compute local histogram with several parameters.
This filter has the only purpose to recall regions.
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.