OTB  9.0.0
Orfeo Toolbox
otbComputeGainLutFilter.hxx
Go to the documentation of this file.
1 /*
2  * Copyright (C) 2005-2022 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 otbComputeGainLutFilter_hxx
22 #define otbComputeGainLutFilter_hxx
23 
25 #include "itkImageRegionIterator.h"
26 
27 #include <limits>
28 #include <numeric>
29 
30 namespace otb
31 {
32 
33 template <class TInputImage, class TOutputImage>
35 {
36  m_NbBin = 256;
37  m_NbPixel = 0;
38  m_Min = std::numeric_limits<double>::quiet_NaN();
39  m_Max = std::numeric_limits<double>::quiet_NaN();
40  m_Step = -1;
41 }
42 
43 template <class TInputImage, class TOutputImage>
45 {
46  m_NbBin = this->GetInput()->GetNumberOfComponentsPerPixel();
47  m_Step = static_cast<double>(m_Max - m_Min) / static_cast<double>(m_NbBin - 1);
48 }
49 
50 template <class TInputImage, class TOutputImage>
52  itk::ThreadIdType itkNotUsed(threadId))
53 {
54  assert(m_Step > 0);
55  assert(m_NbBin > 0);
56  // TODO error
57  // itk::ProgressReporter progress(this , threadId ,
58  // outputRegionForThread.GetNumberOfPixels() );
59 
60  typename InputImageType::ConstPointer input(this->GetInput());
61  typename OutputImageType::Pointer output(this->GetOutput());
62 
63  itk::ImageRegionConstIterator<InputImageType> it(input, outputRegionForThread);
64 
65  itk::ImageRegionIterator<OutputImageType> oit(output, outputRegionForThread);
66 
67  HistoType target;
68  target.SetSize(m_NbBin);
69 
70  HistoType currentHisto;
71 
72  LutType lut;
73  lut.SetSize(m_NbBin);
74 
75  for (it.GoToBegin(), oit.GoToBegin(); !oit.IsAtEnd() || !it.IsAtEnd(); ++oit, ++it)
76  {
77  currentHisto = it.Get();
78  target.Fill(0);
79  lut.Fill(-1);
80  if (IsValid(currentHisto))
81  {
82  CreateTarget(currentHisto, target);
83  Equalized(currentHisto, target, lut);
84  }
85  oit.Set(lut);
86  }
87  assert(oit.IsAtEnd() && it.IsAtEnd());
88 }
89 
90 template <class TInputImage, class TOutputImage>
91 typename TOutputImage::InternalPixelType ComputeGainLutFilter<TInputImage, TOutputImage>::PostProcess(unsigned int countValue, unsigned int countMapValue)
92 {
93  double denum(countValue * m_Step + m_Min);
94  if (denum == 0)
95  return 0;
96  return static_cast<OutputPixelType>((countMapValue * m_Step + m_Min) / denum);
97 }
98 
99 template <class TInputImage, class TOutputImage>
101 {
102  unsigned int countValue(0), countMapValue(0);
103  lut[countValue] = 1; // Black stays black
104  ++countValue;
105  unsigned int countInput(inputHisto[0] + inputHisto[countValue]);
106  lut[m_NbBin - 1] = 1; // White stays white
107  unsigned int countTarget(targetHisto[countMapValue]);
108 
109  while ((countMapValue < m_NbBin) && countValue < (m_NbBin - 1))
110  {
111  if (countInput > countTarget)
112  {
113  ++countMapValue;
114  countTarget += targetHisto[countMapValue];
115  }
116  else
117  {
118  lut[countValue] = PostProcess(countValue, countMapValue);
119  ++countValue;
120  countInput += inputHisto[countValue];
121  }
122  }
123  for (unsigned int i = 0; i < m_NbBin; i++)
124  {
125  if (lut[i] == -1)
126  {
127  lut[i] = 1;
128  }
129  }
130 }
131 
132 template <class TInputImage, class TOutputImage>
134 {
135  unsigned int nbPixel(0);
136  for (unsigned int i = 0; i < m_NbBin; i++)
137  {
138  nbPixel += inputHisto[i];
139  }
140  unsigned int rest(nbPixel % m_NbBin), height(nbPixel / m_NbBin);
141  targetHisto.Fill(height);
142  for (unsigned int i = 0; i < rest; i++)
143  {
144  ++targetHisto[(m_NbBin - rest) / 2 + i];
145  }
146 }
147 
148 template <class TInputImage, class TOutputImage>
150 {
151  long acc = std::accumulate(&inputHisto[0], &inputHisto[m_NbBin - 1], 0);
152  return acc >= (0.5 * m_NbPixel);
153 }
154 
158 template <class TInputImage, class TOutputImage>
159 void ComputeGainLutFilter<TInputImage, TOutputImage>::PrintSelf(std::ostream& os, itk::Indent indent) const
160 {
161  Superclass::PrintSelf(os, indent);
162  os << indent << "Minimum: " << m_Min << std::endl;
163  os << indent << "Maximum: " << m_Max << std::endl;
164  os << indent << "Step: " << m_Step << std::endl;
165  os << indent << "Number of bin: " << m_NbBin << std::endl;
166  os << indent << "Number of pixel by histogram: " << m_NbPixel << std::endl;
167 }
169 
170 } // End namespace otb
171 
172 #endif
otb::ComputeGainLutFilter::IsValid
bool IsValid(const HistoType &inputHisto)
Definition: otbComputeGainLutFilter.hxx:149
otb::ComputeGainLutFilter::ComputeGainLutFilter
ComputeGainLutFilter()
Definition: otbComputeGainLutFilter.hxx:34
otb::ComputeGainLutFilter::OutputPixelType
OutputImageType::InternalPixelType OutputPixelType
Definition: otbComputeGainLutFilter.h:58
otb
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.
Definition: otbJoinContainer.h:32
otb::ComputeGainLutFilter::CreateTarget
void CreateTarget(const HistoType &inputHisto, HistoType &targetHisto)
Definition: otbComputeGainLutFilter.hxx:133
otb::ComputeGainLutFilter::PrintSelf
void PrintSelf(std::ostream &os, itk::Indent indent) const override
Definition: otbComputeGainLutFilter.hxx:159
otb::ComputeGainLutFilter::PostProcess
OutputPixelType PostProcess(unsigned int countMapValue, unsigned int countValue)
Definition: otbComputeGainLutFilter.hxx:91
otb::ComputeGainLutFilter::ThreadedGenerateData
void ThreadedGenerateData(const OutputImageRegionType &outputRegionForThread, itk::ThreadIdType threadId) override
Definition: otbComputeGainLutFilter.hxx:51
otb::ComputeGainLutFilter::BeforeThreadedGenerateData
void BeforeThreadedGenerateData() override
Definition: otbComputeGainLutFilter.hxx:44
otb::ComputeGainLutFilter::Equalized
void Equalized(const HistoType &inputHisto, HistoType &targetHisto, LutType &lut)
Definition: otbComputeGainLutFilter.hxx:100
otb::ComputeGainLutFilter::OutputImageRegionType
OutputImageType::RegionType OutputImageRegionType
Definition: otbComputeGainLutFilter.h:59
otb::ComputeGainLutFilter::LutType
OutputImageType::PixelType LutType
Definition: otbComputeGainLutFilter.h:57
otbComputeGainLutFilter.h
otb::ComputeGainLutFilter::HistoType
InputImageType::PixelType HistoType
Definition: otbComputeGainLutFilter.h:56