OTB  10.0.0
Orfeo Toolbox
otbComputeHistoFilter.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 otbComputeHistoFilter_hxx
22 #define otbComputeHistoFilter_hxx
23 
24 #include "otbComputeHistoFilter.h"
25 #include "otbMacro.h"
26 #include <limits>
27 
28 namespace otb
29 {
30 
31 template <class TInputImage, class TOutputImage>
33 {
34  this->DynamicMultiThreadingOff();
35  this->SetNumberOfRequiredOutputs(2);
36  this->SetNthOutput(0, this->MakeOutput(0));
37  this->SetNthOutput(1, this->MakeOutput(1));
38  m_Min = std::numeric_limits<InputPixelType>::quiet_NaN();
39  m_Max = std::numeric_limits<InputPixelType>::quiet_NaN();
40  m_NoDataFlag = false;
41  m_NoData = std::numeric_limits<InputPixelType>::quiet_NaN();
42  m_NbBin = 256;
43  m_Threshold = -1;
44  m_ThumbSize.Fill(0);
45  m_ValidThreads = 1;
46  m_Step = -1;
47 }
48 
49 template <class TInputImage, class TOutputImage>
50 itk::ProcessObject::DataObjectPointer ComputeHistoFilter<TInputImage, TOutputImage>::MakeOutput(itk::ProcessObject::DataObjectPointerArraySizeType idx)
51 {
52  itk::DataObject::Pointer output;
53 
54  switch (idx)
55  {
56  case 0:
57  output = (OutputImageType::New()).GetPointer();
58  break;
59  case 1:
60  output = (OutputImageType::New()).GetPointer();
61  break;
62  default:
63  std::cerr << "No output " << idx << std::endl;
64  output = NULL;
65  break;
66  }
67  return output;
68 }
69 
70 template <class TInputImage, class TOutputImage>
71 itk::ProcessObject::DataObjectPointer ComputeHistoFilter<TInputImage, TOutputImage>::MakeOutput(const itk::ProcessObject::DataObjectIdentifierType& name)
72 {
73  return Superclass::MakeOutput(name);
74 }
75 
76 template <class TInputImage, class TOutputImage>
78 {
79  typename Superclass::InputImagePointer inputPtr(const_cast<InputImageType*>(this->GetInput()));
80  inputPtr->SetRequestedRegion(this->GetOutput()->GetRequestedRegion());
81  if (inputPtr->GetRequestedRegion().GetNumberOfPixels() == 0)
82  {
83  inputPtr->SetRequestedRegionToLargestPossibleRegion();
84  }
85 }
86 
87 template <class TInputImage, class TOutputImage>
89 {
90  Superclass::GenerateOutputInformation();
91  typename InputImageType::ConstPointer input(this->GetInput());
92  typename OutputImageType::Pointer output(this->GetHistoOutput());
93  typename OutputImageType::Pointer outImage(this->GetOutput());
94 
95  typename InputImageType::RegionType inputLargestRegion(input->GetLargestPossibleRegion());
96  outImage->SetLargestPossibleRegion(inputLargestRegion);
97 
98  typename OutputImageType::IndexType start;
99  typename OutputImageType::SizeType size;
100 
101  start.Fill(0);
102 
103  assert(m_ThumbSize[0] != 0);
104  assert(m_ThumbSize[1] != 0);
105 
106  // TODO throw error if 0
107 
108  size[0] = std::ceil(inputLargestRegion.GetSize()[0] / static_cast<double>(m_ThumbSize[0]));
109  size[1] = std::ceil(inputLargestRegion.GetSize()[1] / static_cast<double>(m_ThumbSize[1]));
110 
111  typename OutputImageType::RegionType region;
112  region.SetSize(size);
113  region.SetIndex(start);
114  output->SetNumberOfComponentsPerPixel(m_NbBin);
115  output->SetLargestPossibleRegion(region);
116  typename InputImageType::SpacingType inputSpacing(input->GetSignedSpacing());
117  typename InputImageType::PointType inputOrigin(input->GetOrigin());
118 
119  typename OutputImageType::SpacingType histoSpacing;
120  histoSpacing[0] = inputSpacing[0] * m_ThumbSize[0];
121  histoSpacing[1] = inputSpacing[1] * m_ThumbSize[1];
122  output->SetSignedSpacing(histoSpacing);
123 
124  typename OutputImageType::PointType histoOrigin;
125  histoOrigin[0] = histoSpacing[0] / 2 + inputOrigin[0] - inputSpacing[0] / 2;
126  histoOrigin[1] = histoSpacing[1] / 2 + inputOrigin[1] - inputSpacing[1] / 2;
127  output->SetOrigin(histoOrigin);
128 }
129 
130 template <class TInputImage, class TOutputImage>
132 {
133  if (GetHistoOutput()->GetRequestedRegion().GetNumberOfPixels() == 0)
134  {
135  GetHistoOutput()->SetRequestedRegionToLargestPossibleRegion();
136  }
137  typename OutputImageType::Pointer outImage(this->GetOutput());
138  SetRequestedRegion(outImage);
139 }
140 
141 template <class TInputImage, class TOutputImage>
143 {
144  this->AllocateOutputs();
145 
146  // Set up the multithreaded processing
147  typename itk::ImageSource<OutputImageType>::ThreadStruct str;
148  str.Filter = this;
149 
150  // Get the output pointer
151  const OutputImageType* outputPtr(this->GetOutput());
152  const itk::ImageRegionSplitterBase* splitter(this->GetImageRegionSplitter());
153  m_ValidThreads = splitter->GetNumberOfSplits(outputPtr->GetRequestedRegion(), this->GetNumberOfWorkUnits());
154 
155  this->BeforeThreadedGenerateData();
156 
157  this->GetMultiThreader()->SetNumberOfWorkUnits(m_ValidThreads);
158  this->GetMultiThreader()->SetSingleMethod(this->ThreaderCallback, &str);
159  // multithread the execution
160  this->GetMultiThreader()->SingleMethodExecute();
161 
162  this->AfterThreadedGenerateData();
163 }
164 
165 template <class TInputImage, class TOutputImage>
167 {
168  // Initializing output
169  typename OutputImageType::Pointer output(this->GetHistoOutput());
170  typename OutputImageType::PixelType zeroPixel(m_NbBin);
171  zeroPixel.Fill(0);
172  output->FillBuffer(zeroPixel);
173 
174  // Initializing shared variable with thread number parameter
175  SizeType outSize(output->GetRequestedRegion().GetSize());
176  m_HistoThread.resize(m_ValidThreads * outSize[0] * outSize[1]);
177  m_HistoThread.shrink_to_fit();
178  std::fill(m_HistoThread.begin(), m_HistoThread.end(), zeroPixel);
179 
180  m_Step = static_cast<double>(m_Max - m_Min) / static_cast<double>(m_NbBin - 1);
181 }
182 
183 template <class TInputImage, class TOutputImage>
184 void ComputeHistoFilter<TInputImage, TOutputImage>::ThreadedGenerateData(const OutputImageRegionType& outputRegionForThread, itk::ThreadIdType threadId)
185 {
186  assert(m_Step > 0);
187  typename InputImageType::ConstPointer input(this->GetInput());
188  typename OutputImageType::Pointer output(this->GetHistoOutput());
189 
190  OutputImageRegionType histoRegion(GetHistoOutput()->GetRequestedRegion());
191  SizeType outSize(histoRegion.GetSize());
192  IndexType outIndex(histoRegion.GetIndex());
193 
194  typename InputImageType::RegionType region;
195 
196  unsigned int threadIndex(threadId * outSize[0] * outSize[1]), pixel(0);
197 
198  for (unsigned int nthHisto = 0; nthHisto < outSize[0] * outSize[1]; nthHisto++)
199  {
200  IndexType start;
201  start[0] = (outIndex[0] + nthHisto % outSize[0]) * m_ThumbSize[0];
202  start[1] = (outIndex[1] + nthHisto / outSize[0]) * m_ThumbSize[1];
203  region.SetSize(m_ThumbSize);
204  region.SetIndex(start);
205 
206  if (!region.Crop(outputRegionForThread))
207  continue;
208 
209  typename itk::ImageRegionConstIterator<InputImageType> it(input, region);
210  InputPixelType currentPixel(0);
211  for (it.GoToBegin(); !it.IsAtEnd(); ++it)
212  {
213  currentPixel = it.Get();
214  if ((currentPixel == m_NoData && m_NoDataFlag) || currentPixel > m_Max || currentPixel < m_Min)
215  continue;
216 
217  pixel = static_cast<unsigned int>(std::round((currentPixel - m_Min) / m_Step));
218  ++m_HistoThread[threadIndex + nthHisto][pixel];
219  }
220  }
221 }
222 
223 template <class TInputImage, class TOutputImage>
225 {
226  typename OutputImageType::Pointer output(this->GetHistoOutput());
227  typename itk::ImageRegionIterator<OutputImageType> oit(output, output->GetRequestedRegion());
228  OutputImageRegionType histoRegion(GetHistoOutput()->GetRequestedRegion());
229  SizeType outSize(histoRegion.GetSize());
230  IndexType outIndex(histoRegion.GetIndex());
231 
232  unsigned int agreg(0), total(0);
233 
234  for (oit.GoToBegin(); !oit.IsAtEnd(); ++oit)
235  {
236  total = 0;
237 
238  for (unsigned int i = 0; i < m_NbBin; i++)
239  {
240  agreg = 0;
241 
242  for (unsigned int threadId = 0; threadId < m_ValidThreads; threadId++)
243  {
244  agreg += m_HistoThread[threadId * outSize[0] * outSize[1] + (oit.GetIndex()[1] - outIndex[1]) * outSize[0] + ((oit.GetIndex()[0] - outIndex[0]))][i];
245  }
246  oit.Get()[i] = agreg;
247  total += agreg;
248  }
249  if (m_Threshold > 0)
250  ApplyThreshold(oit, total);
251  }
252 }
253 
254 template <class TInputImage, class TOutputImage>
255 void ComputeHistoFilter<TInputImage, TOutputImage>::ApplyThreshold(typename itk::ImageRegionIterator<OutputImageType> oit, unsigned int total)
256 {
257  unsigned int rest(0);
258  unsigned int height(static_cast<unsigned int>(m_Threshold * (total / m_NbBin)));
259  // Do i need to static_cast in a an assignation?
260  // Warning!!!! Need to handle out of bound int!!!
261 
262  for (unsigned int i = 0; i < m_NbBin; i++)
263  {
264  if (static_cast<unsigned int>(oit.Get()[i]) > height)
265  {
266  rest += oit.Get()[i] - height;
267  oit.Get()[i] = height;
268  }
269  }
270  height = rest / m_NbBin;
271  rest = rest % m_NbBin;
272  for (unsigned int i = 0; i < m_NbBin; i++)
273  {
274  oit.Get()[i] += height;
275  if (i > (m_NbBin - rest) / 2 && i <= (m_NbBin - rest) / 2 + rest)
276  {
277  ++oit.Get()[i];
278  }
279  }
280 }
281 
282 template <class TInputImage, class TOutputImage>
284 {
285  assert(this->itk::ProcessObject::GetOutput(1));
286 
287  return dynamic_cast<TOutputImage*>(this->itk::ProcessObject::GetOutput(1));
288 }
289 
290 template <class TInputImage, class TOutputImage>
292 {
293  OutputImageRegionType histoRegion(GetHistoOutput()->GetRequestedRegion());
294 
295  IndexType start;
296  start[0] = histoRegion.GetIndex()[0] * m_ThumbSize[0];
297  start[1] = histoRegion.GetIndex()[1] * m_ThumbSize[1];
298 
299  SizeType size;
300  size[0] = histoRegion.GetSize()[0] * m_ThumbSize[0];
301  size[1] = histoRegion.GetSize()[1] * m_ThumbSize[1];
302 
303  typename OutputImageType::RegionType outputRequestedRegion;
304  outputRequestedRegion.SetIndex(start);
305  outputRequestedRegion.SetSize(size);
306 
307  outputRequestedRegion.Crop(image->GetLargestPossibleRegion());
308  image->SetRequestedRegion(outputRequestedRegion);
309 }
310 
311 template <class TInputImage, class TOutputImage>
312 void ComputeHistoFilter<TInputImage, TOutputImage>::PrintSelf(std::ostream& os, itk::Indent indent) const
313 {
314  Superclass::PrintSelf(os, indent);
315  os << indent << "Is no data activated: " << m_NoDataFlag << std::endl;
316  os << indent << "No Data: " << m_NoData << std::endl;
317  os << indent << "Minimum: " << m_Min << std::endl;
318  os << indent << "Maximum: " << m_Max << std::endl;
319  os << indent << "Step: " << m_Step << std::endl;
320  os << indent << "Number of bin: " << m_NbBin << std::endl;
321  os << indent << "Thumbnail size: " << m_ThumbSize << std::endl;
322  os << indent << "Threshold value: " << m_Threshold << std::endl;
323 }
324 
325 } // End namespace otb
326 
327 #endif
void AfterThreadedGenerateData() override
void GenerateOutputInformation() override
virtual void ThreadedGenerateData(const OutputImageRegionType &outputRegionForThread, itk::ThreadIdType threadId) override
void BeforeThreadedGenerateData() override
virtual itk::ProcessObject::DataObjectPointer MakeOutput(itk::ProcessObject::DataObjectPointerArraySizeType idx) override
OutputImageType::Pointer GetHistoOutput()
InputImageType::InternalPixelType InputPixelType
void PrintSelf(std::ostream &os, itk::Indent indent) const override
void SetRequestedRegion(itk::ImageBase< 2 > *image)
void GenerateInputRequestedRegion() override
InputImageType::IndexType IndexType
InputImageType::SizeType SizeType
void ApplyThreshold(typename itk::ImageRegionIterator< OutputImageType > oit, unsigned int total)
void GenerateOutputRequestedRegion(itk::DataObject *output) override
OutputImageType::RegionType OutputImageRegionType
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.