OTB  10.0.0
Orfeo Toolbox
otbStreamingStatisticsMapFromLabelImageFilter.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 otbStreamingStatisticsMapFromLabelImageFilter_hxx
23 #define otbStreamingStatisticsMapFromLabelImageFilter_hxx
25 
26 #include "itkInputDataObjectIterator.h"
27 #include "itkImageRegionIterator.h"
28 #include "itkProgressReporter.h"
29 #include "otbMacro.h"
30 #include <cmath>
31 #include <utility>
32 
33 namespace otb
34 {
35 
36 template <class TInputVectorImage, class TLabelImage>
38  : m_UseNoDataValue()
39 {
40  // first output is a copy of the image, DataObject created by
41  // superclass
42  //
43  // allocate the data objects for the outputs which are
44  // just decorators around pixel types
45  this->DynamicMultiThreadingOff();
46  typename PixelValueMapObjectType::Pointer output
47  = static_cast<PixelValueMapObjectType*>(this->MakeOutput(1).GetPointer());
48  this->itk::ProcessObject::SetNthOutput(1, output.GetPointer());
49 
50  this->Reset();
51 }
52 
53 template <class TInputVectorImage, class TLabelImage>
54 typename itk::DataObject::Pointer
56 {
57  return static_cast<itk::DataObject*>(PixelValueMapObjectType::New().GetPointer());
58 }
59 
60 template <class TInputVectorImage, class TLabelImage>
62 {
63  // Process object is not const-correct so the const_cast is required here
64  this->itk::ProcessObject::SetNthInput(1, const_cast<LabelImageType*>(input));
65 }
66 
67 template <class TInputVectorImage, class TLabelImage>
70 {
71  return static_cast<const TLabelImage*>(this->itk::ProcessObject::GetInput(1));
72 }
73 
74 template <class TInputVectorImage, class TLabelImage>
77 {
78  return m_MeanRadiometricValue;
79 }
80 
81 template <class TInputVectorImage, class TLabelImage>
84 {
85  return m_StDevRadiometricValue;
86 }
87 
88 template <class TInputVectorImage, class TLabelImage>
91 {
92  return m_MinRadiometricValue;
93 }
94 
95 template <class TInputVectorImage, class TLabelImage>
98 {
99  return m_MaxRadiometricValue;
100 }
101 
102 template <class TInputVectorImage, class TLabelImage>
105 {
106  return m_LabelPopulation;
107 }
108 
109 template <class TInputVectorImage, class TLabelImage>
111 {
112  Superclass::GenerateOutputInformation();
113 
114  if (this->GetInput())
115  {
116  this->GetOutput()->CopyInformation(this->GetInput());
117  this->GetOutput()->SetLargestPossibleRegion(this->GetInput()->GetLargestPossibleRegion());
118 
119  if (this->GetOutput()->GetRequestedRegion().GetNumberOfPixels() == 0)
120  {
121  this->GetOutput()->SetRequestedRegion(this->GetOutput()->GetLargestPossibleRegion());
122  }
123  }
124 }
125 
126 template <class TInputVectorImage, class TLabelImage>
128 {
129  // This is commented to prevent the streaming of the whole image for the first stream strip
130  // It shall not cause any problem because the output image of this filter is not intended to be used.
131  // InputImagePointer image = const_cast< TInputImage * >( this->GetInput() );
132  // this->GraftOutput( image );
133  // Nothing that needs to be allocated for the remaining outputs
134 }
135 
136 template <class TInputVectorImage, class TLabelImage>
138 {
139  // Update temporary accumulator
140  AccumulatorMapType outputAcc;
141  auto endAcc = outputAcc.end();
142 
143  for (auto const& threadAccMap : m_AccumulatorMaps)
144  {
145  for (auto const& it : threadAccMap)
146  {
147  auto label = it.first;
148  auto itAcc = outputAcc.find(label);
149  if (itAcc == endAcc)
150  {
151  outputAcc.emplace(label, it.second);
152  }
153  else
154  {
155  itAcc->second.Update(it.second);
156  }
157  }
158  }
159 
160  // Publish output maps
161  for (auto& it : outputAcc)
162  {
163  const LabelPixelType label = it.first;
164  const auto& bandCount = it.second.GetBandCount();
165  const auto& sum = it.second.GetSum();
166  const auto& sqSum = it.second.GetSqSum();
167 
168  // Count
169  m_LabelPopulation[label] = it.second.GetCount();
170 
171  // Mean & stdev
173  RealVectorPixelType std(sqSum);
174  RealVectorPixelType min(it.second.GetMin());
175  RealVectorPixelType max(it.second.GetMax());
176  for (unsigned int band = 0; band < mean.GetSize(); band++)
177  {
178  // Number of valid pixels in band
179  auto count = bandCount[band];
180  // Mean
181  mean[band] /= count;
182 
183  // Unbiased standard deviation (not sure unbiased is useful here)
184  const double variance = (sqSum[band] - (sum[band] * mean[band])) / (count - 1);
185  std[band] = std::sqrt(variance);
186 
187  // Use the no data value when no valid pixels were found
188  if (this->GetUseNoDataValue() && count == 0)
189  {
190  min[band] = this->GetNoDataValue();
191  max[band] = this->GetNoDataValue();
192  }
193  }
194  m_MeanRadiometricValue.emplace(label, std::move(mean));
195  m_StDevRadiometricValue.emplace(label, std::move(std));
196 
197  // Min & max
198  m_MinRadiometricValue.emplace(label, std::move(min));
199  m_MaxRadiometricValue.emplace(label, std::move(max));
200  }
201 }
202 
203 template <class TInputVectorImage, class TLabelImage>
205 {
206  m_AccumulatorMaps.clear();
207 
208  m_MeanRadiometricValue.clear();
209  m_StDevRadiometricValue.clear();
210  m_MinRadiometricValue.clear();
211  m_MaxRadiometricValue.clear();
212  m_LabelPopulation.clear();
213  m_AccumulatorMaps.resize(this->GetNumberOfWorkUnits());
214 }
215 
216 template <class TInputVectorImage, class TLabelImage>
218 {
219  // The Requested Regions of all the inputs are set to their Largest Possible Regions
220  this->itk::ProcessObject::GenerateInputRequestedRegion();
221 
222  // Iteration over all the inputs of the current filter (this)
223  for (itk::InputDataObjectIterator it(this); !it.IsAtEnd(); it++)
224  {
225  // Check whether the input is an image of the appropriate dimension
226  // dynamic_cast of all the input images as itk::ImageBase objects
227  // in order to pass the if ( input ) test whatever the inputImageType (vectorImage or labelImage)
228  ImageBaseType* input = dynamic_cast<ImageBaseType*>(it.GetInput());
229 
230  if (input)
231  {
232  // Use the function object RegionCopier to copy the output region
233  // to the input. The default region copier has default implementations
234  // to handle the cases where the input and output are the same
235  // dimension, the input a higher dimension than the output, and the
236  // input a lower dimension than the output.
237  InputImageRegionType inputRegion;
238  this->CallCopyOutputRegionToInputRegion(inputRegion, this->GetOutput()->GetRequestedRegion());
239  input->SetRequestedRegion(inputRegion);
240  }
241  }
242 }
243 
244 template <class TInputVectorImage, class TLabelImage>
246  itk::ThreadIdType threadId)
247 {
251  InputVectorImagePointer inputPtr = const_cast<TInputVectorImage*>(this->GetInput());
252  LabelImagePointer labelInputPtr = const_cast<TLabelImage*>(this->GetInputLabelImage());
254 
255  itk::ImageRegionConstIterator<TInputVectorImage> inIt(inputPtr, outputRegionForThread);
256  itk::ImageRegionConstIterator<TLabelImage> labelIt(labelInputPtr, outputRegionForThread);
257  itk::ProgressReporter progress(this, threadId, outputRegionForThread.GetNumberOfPixels());
258 
259  auto& acc = m_AccumulatorMaps[threadId];
260  auto endAcc = acc.end();
261 
262  // do the work
263  for (inIt.GoToBegin(), labelIt.GoToBegin(); !inIt.IsAtEnd() && !labelIt.IsAtEnd(); ++inIt, ++labelIt)
264  {
265  const auto& value = inIt.Get();
266  auto label = labelIt.Get();
267 
268  // Update the accumulator
269  auto itAcc = acc.find(label);
270  if (itAcc == endAcc)
271  {
272  acc.emplace(label, AccumulatorType(this->GetNoDataValue(), this->GetUseNoDataValue(), value));
273  }
274  else
275  {
276  itAcc->second.Update(value);
277  }
278 
279  progress.CompletedPixel();
280  }
281 }
282 
283 template <class TInputVectorImage, class TLabelImage>
285 {
286  Superclass::PrintSelf(os, indent);
287 }
288 
289 } // end namespace otb
290 #endif
std::unordered_map< LabelPixelType, RealVectorPixelType > PixelValueMapType
void ThreadedGenerateData(const RegionType &outputRegionForThread, itk::ThreadIdType threadId) override
DataObjectPointer MakeOutput(DataObjectPointerArraySizeType idx) override
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.