OTB  10.0.0
Orfeo Toolbox
otbStreamingHistogramVectorImageFilter.hxx
Go to the documentation of this file.
1 /*
2  * Copyright (C) 1999-2011 Insight Software Consortium
3  * Copyright (C) 2005-2024 Centre National d'Etudes Spatiales (CNES)
4  *
5  * This file is part of Orfeo Toolbox
6  *
7  * https://www.orfeo-toolbox.org/
8  *
9  * Licensed under the Apache License, Version 2.0 (the "License");
10  * you may not use this file except in compliance with the License.
11  * You may obtain a copy of the License at
12  *
13  * http://www.apache.org/licenses/LICENSE-2.0
14  *
15  * Unless required by applicable law or agreed to in writing, software
16  * distributed under the License is distributed on an "AS IS" BASIS,
17  * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
18  * See the License for the specific language governing permissions and
19  * limitations under the License.
20  */
21 
22 #ifndef otbStreamingHistogramVectorImageFilter_hxx
23 #define otbStreamingHistogramVectorImageFilter_hxx
25 
26 #include "itkImageRegionIterator.h"
27 #include "itkImageRegionConstIteratorWithIndex.h"
28 #include "itkProgressReporter.h"
29 #include "otbMacro.h"
30 
31 namespace otb
32 {
33 
34 template <class TInputImage>
36  : m_ThreadHistogramList(),
37  m_Size(),
38  m_HistogramMin(),
39  m_HistogramMax(),
40  m_NoDataFlag(false),
41  m_NoDataValue(itk::NumericTraits<InternalPixelType>::Zero),
42  m_SubSamplingRate(1)
43 {
44  // first output is a copy of the image, DataObject created by
45  // superclass
46  //
47  // allocate the data objects for the outputs which are
48  // just decorators around pixel types and histogram list
49  this->DynamicMultiThreadingOff();
50  m_Size.Fill(255);
51  HistogramListPointerType output = static_cast<HistogramListType*>(this->MakeOutput(1).GetPointer());
52  this->itk::ProcessObject::SetNthOutput(1, output.GetPointer());
53 }
54 
55 template <class TInputImage>
57 {
58  itk::DataObject::Pointer ret;
59  switch (output)
60  {
61  case 0:
62  ret = static_cast<itk::DataObject*>(TInputImage::New().GetPointer());
63  break;
64  case 1:
65  ret = static_cast<itk::DataObject*>(HistogramListType::New().GetPointer());
66  break;
67  }
68  return ret;
69 }
70 
71 template <class TInputImage>
73 {
74  return static_cast<HistogramListType*>(this->itk::ProcessObject::GetOutput(1));
75 }
76 
77 template <class TInputImage>
80 {
81  return static_cast<const HistogramListType*>(this->itk::ProcessObject::GetOutput(1));
82 }
83 
84 template <class TInputImage>
86 {
87  Superclass::GenerateOutputInformation();
88  if (this->GetInput())
89  {
90  this->GetOutput()->CopyInformation(this->GetInput());
91  this->GetOutput()->SetLargestPossibleRegion(this->GetInput()->GetLargestPossibleRegion());
92 
93  if (this->GetOutput()->GetRequestedRegion().GetNumberOfPixels() == 0)
94  {
95  this->GetOutput()->SetRequestedRegion(this->GetOutput()->GetLargestPossibleRegion());
96  }
97  }
98 }
99 
100 template <class TInputImage>
102 {
103  // This is commented to prevent the streaming of the whole image for the first stream strip
104  // It shall not cause any problem because the output image of this filter is not intended to be used.
105  // InputImagePointer image = const_cast< TInputImage * >( this->GetInput() );
106  // this->GraftOutput( image );
107  // Nothing that needs to be allocated for the remaining outputs
108 }
109 
110 template <class TInputImage>
112 {
113  TInputImage* inputPtr = const_cast<TInputImage*>(this->GetInput());
114  inputPtr->UpdateOutputInformation();
115 
116  unsigned int numberOfThreads = this->GetNumberOfWorkUnits();
117  unsigned int numberOfComponent = inputPtr->GetNumberOfComponentsPerPixel();
118 
119  // TODO which is the good value ? (false in MVD2)
120  bool clipBins = false;
121 
122  // if histogram Min and Max have the wrong size : set to default [0, 255]
123  if (m_HistogramMin.Size() != numberOfComponent || m_HistogramMax.Size() != numberOfComponent)
124  {
125  m_HistogramMin.SetSize(numberOfComponent);
126  m_HistogramMax.SetSize(numberOfComponent);
127 
128  m_HistogramMin.Fill(itk::NumericTraits<InternalPixelType>::Zero);
129  m_HistogramMax.Fill(255);
130  }
131 
132  // Setup output histogram
133  HistogramListType* outputHisto = this->GetHistogramListOutput();
134  outputHisto->Clear();
135  for (unsigned int k = 0; k < numberOfComponent; ++k)
136  {
137  typename HistogramType::MeasurementVectorType bandMin, bandMax;
138  bandMin.SetSize(1);
139  bandMax.SetSize(1);
140  bandMin.Fill(m_HistogramMin[k]);
141  bandMax.Fill(m_HistogramMax[k]);
142 
143  typename HistogramType::Pointer histogram = HistogramType::New();
144  histogram->SetClipBinsAtEnds(clipBins);
145 
146  typename HistogramType::SizeType size;
147  size.SetSize(1);
148  size.Fill(m_Size[k]);
149  histogram->SetMeasurementVectorSize(1);
150  histogram->Initialize(size, bandMin, bandMax);
151 
152  outputHisto->PushBack(histogram);
153  }
154 
155 
156  // Setup HistogramLists for each thread
157  m_ThreadHistogramList.clear();
158  for (unsigned int i = 0; i < numberOfThreads; ++i)
159  {
160  HistogramListPointerType histoList = HistogramListType::New();
161  histoList->Clear();
162  for (unsigned int k = 0; k < numberOfComponent; ++k)
163  {
164  typename HistogramType::MeasurementVectorType bandMin, bandMax;
165  bandMin.SetSize(1);
166  bandMax.SetSize(1);
167  bandMin.Fill(m_HistogramMin[k]);
168  bandMax.Fill(m_HistogramMax[k]);
169 
170  typename HistogramType::Pointer histogram = HistogramType::New();
171  histogram->SetClipBinsAtEnds(clipBins);
172 
173  typename HistogramType::SizeType size;
174  size.SetSize(1);
175  size.Fill(m_Size[k]);
176  histogram->SetMeasurementVectorSize(1);
177  histogram->Initialize(size, bandMin, bandMax);
178 
179  histoList->PushBack(histogram);
180  }
181  m_ThreadHistogramList.push_back(histoList);
182  }
183 }
184 
185 template <class TInputImage>
187 {
188  HistogramListType* outputHisto = this->GetHistogramListOutput();
189 
190  int numberOfThreads = this->GetNumberOfWorkUnits();
191  unsigned int numberOfComponent = this->GetInput()->GetNumberOfComponentsPerPixel();
192 
193  // copy histograms to output
194  for (int i = 0; i < numberOfThreads; ++i)
195  {
196  for (unsigned int j = 0; j < numberOfComponent; ++j)
197  {
198  HistogramType* outHisto = outputHisto->GetNthElement(j);
199  HistogramType* threadHisto = m_ThreadHistogramList[i]->GetNthElement(j);
200 
201  typename HistogramType::Iterator iterOutput = outHisto->Begin();
202  typename HistogramType::Iterator iterThread = threadHisto->Begin();
203 
204  while (iterOutput != outHisto->End() && iterThread != threadHisto->End())
205  {
206  iterOutput.SetFrequency(iterOutput.GetFrequency() + iterThread.GetFrequency());
207 
208  ++iterOutput;
209  ++iterThread;
210  }
211  }
212  }
213 }
214 
215 template <class TInputImage>
216 void PersistentHistogramVectorImageFilter<TInputImage>::ThreadedGenerateData(const RegionType& outputRegionForThread, itk::ThreadIdType threadId)
217 {
221  InputImagePointer inputPtr = const_cast<TInputImage*>(this->GetInput());
222  // support progress methods/callbacks
223  itk::ProgressReporter progress(this, threadId, outputRegionForThread.GetNumberOfPixels());
225 
226  typename HistogramType::IndexType index;
227 
228  itk::ImageRegionConstIteratorWithIndex<TInputImage> it(inputPtr, outputRegionForThread);
229  it.GoToBegin();
230 
231  bool skipSample = false;
232 
233  // do the work
234  while (!it.IsAtEnd())
235  {
236  if (m_SubSamplingRate > 1)
237  {
238  skipSample = false;
239  for (unsigned int i = 0; i < InputImageDimension; ++i)
240  {
241  if (it.GetIndex()[i] % m_SubSamplingRate != 0)
242  {
243  skipSample = true;
244  break;
245  }
246  }
247  if (skipSample)
248  {
249  ++it;
250  progress.CompletedPixel();
251  continue;
252  }
253  }
254 
255  PixelType vectorValue = it.Get();
256 
257  bool skipSampleNoData = false;
258  if (m_NoDataFlag)
259  {
260  unsigned int itComp = 0;
261  while (itComp < vectorValue.GetSize())
262  {
263  if (vectorValue[itComp] == m_NoDataValue)
264  {
265  skipSampleNoData = true;
266  itComp++;
267  }
268  else
269  {
270  skipSampleNoData = false;
271  break;
272  }
273  }
274  }
275 
276  if (!skipSampleNoData)
277  {
278  for (unsigned int j = 0; j < vectorValue.GetSize(); ++j)
279  {
280  typename HistogramType::MeasurementVectorType value;
281  value.SetSize(1);
282  value.Fill(vectorValue[j]);
283 
284  m_ThreadHistogramList[threadId]->GetNthElement(j)->GetIndex(value, index);
285  if (!m_ThreadHistogramList[threadId]->GetNthElement(j)->IsIndexOutOfBounds(index))
286  {
287  // if the measurement vector is out of bound then
288  // the GetIndex method has returned an index set to the max size of
289  // the invalid dimension - even if the hvector is less than the minimum
290  // bin value.
291  // If the index isn't valid, we don't increase the frequency.
292  // See the comments in Histogram->GetIndex() for more info.
293  m_ThreadHistogramList[threadId]->GetNthElement(j)->IncreaseFrequencyOfIndex(index, 1);
294  }
295  }
296  }
297 
298  ++it;
299  progress.CompletedPixel();
300  }
301 }
302 
303 template <class TImage>
304 void PersistentHistogramVectorImageFilter<TImage>::PrintSelf(std::ostream& os, itk::Indent indent) const
305 {
306  Superclass::PrintSelf(os, indent);
307 
308  os << indent << "Histogram minimum: " << this->GetHistogramMin() << std::endl;
309  os << indent << "Histogram maximum: " << this->GetHistogramMax() << std::endl;
310  os << indent << "Number of bins: " << m_Size[0] << std::endl;
311  if (m_NoDataFlag)
312  {
313  os << indent << "Use NoData: true" << std::endl;
314  }
315  else
316  {
317  os << indent << "Use NoData: false" << std::endl;
318  }
319  os << indent << "NoData value: " << this->GetNoDataValue() << std::endl;
320 }
321 
322 } // end namespace otb
323 #endif
This class is a generic all-purpose wrapping around an std::vector<itk::SmartPointer<ObjectType> >.
Definition: otbObjectList.h:41
ObjectPointerType GetNthElement(unsigned int index) const
void PushBack(ObjectType *element)
DataObjectPointer MakeOutput(DataObjectPointerArraySizeType idx) override
itk::ProcessObject::DataObjectPointerArraySizeType DataObjectPointerArraySizeType
itk::Statistics::Histogram< HistogramMeasurementRealType, DFContainerType > HistogramType
void PrintSelf(std::ostream &os, itk::Indent indent) const override
void ThreadedGenerateData(const RegionType &outputRegionForThread, itk::ThreadIdType threadId) override
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.