OTB  10.0.0
Orfeo Toolbox
otbStreamingMinMaxImageFilter.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 otbStreamingMinMaxImageFilter_hxx
23 #define otbStreamingMinMaxImageFilter_hxx
25 
26 #include <algorithm>
27 #include "itkImageRegionIterator.h"
28 #include "itkProgressReporter.h"
29 #include "otbMacro.h"
30 
31 namespace otb
32 {
33 
34 template <class TInputImage>
36 {
37  this->DynamicMultiThreadingOff();
38  // TODO : SetNumberOfRequiredOutputs
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 & index types
45  for (int i = 1; i < 5; ++i)
46  {
47  this->itk::ProcessObject::SetNthOutput(i, this->MakeOutput(i));
48  }
49 
50  this->GetMinimumOutput()->Set(itk::NumericTraits<PixelType>::max());
51  this->GetMaximumOutput()->Set(itk::NumericTraits<PixelType>::NonpositiveMin());
52 
53  this->Reset();
54 }
55 
56 template <class TInputImage>
58 {
59  itk::DataObject::Pointer ret;
60  switch (output)
61  {
62  case 0:
63  ret = static_cast<itk::DataObject*>(TInputImage::New().GetPointer());
64  break;
65  case 1:
66  case 2:
67  ret = static_cast<itk::DataObject*>(PixelObjectType::New().GetPointer());
68  break;
69  case 3:
70  case 4:
71  ret = static_cast<itk::DataObject*>(IndexObjectType::New().GetPointer());
72  break;
73  }
74  return ret;
75 }
76 
77 template <class TInputImage>
79 {
80  return static_cast<PixelObjectType*>(this->itk::ProcessObject::GetOutput(1));
81 }
82 
83 template <class TInputImage>
85 {
86  return static_cast<const PixelObjectType*>(this->itk::ProcessObject::GetOutput(1));
87 }
88 
89 template <class TInputImage>
91 {
92  return static_cast<PixelObjectType*>(this->itk::ProcessObject::GetOutput(2));
93 }
94 
95 template <class TInputImage>
97 {
98  return static_cast<const PixelObjectType*>(this->itk::ProcessObject::GetOutput(2));
99 }
100 
101 
102 template <class TInputImage>
104 {
105  return static_cast<IndexObjectType*>(this->itk::ProcessObject::GetOutput(3));
106 }
107 
108 template <class TInputImage>
110 {
111  return static_cast<const IndexObjectType*>(this->itk::ProcessObject::GetOutput(3));
112 }
113 
114 template <class TInputImage>
116 {
117  return static_cast<IndexObjectType*>(this->itk::ProcessObject::GetOutput(4));
118 }
119 
120 template <class TInputImage>
122 {
123  return static_cast<const IndexObjectType*>(this->itk::ProcessObject::GetOutput(4));
124 }
125 
126 template <class TInputImage>
128 {
129  Superclass::GenerateOutputInformation();
130  if (this->GetInput())
131  {
132  this->GetOutput()->CopyInformation(this->GetInput());
133  this->GetOutput()->SetLargestPossibleRegion(this->GetInput()->GetLargestPossibleRegion());
134 
135  if (this->GetOutput()->GetRequestedRegion().GetNumberOfPixels() == 0)
136  {
137  this->GetOutput()->SetRequestedRegion(this->GetOutput()->GetLargestPossibleRegion());
138  }
139  }
140 }
141 template <class TInputImage>
143 {
144  // This is commented to prevent the streaming of the whole image for the first stream strip
145  // It shall not cause any problem because the output image of this filter is not intended to be used.
146  // InputImagePointer image = const_cast< TInputImage * >( this->GetInput() );
147  // this->GraftOutput( image );
148  // Nothing that needs to be allocated for the remaining outputs
149 }
150 
151 template <class TInputImage>
153 {
154  int i;
155  int numberOfThreads = this->GetNumberOfWorkUnits();
156 
157  PixelType minimum = itk::NumericTraits<PixelType>::max();
158  PixelType maximum = itk::NumericTraits<PixelType>::NonpositiveMin();
159  IndexType minimumIdx;
160  IndexType maximumIdx;
161 
162  for (i = 0; i < numberOfThreads; ++i)
163  {
164  if (m_ThreadMin[i] < minimum)
165  {
166  minimum = m_ThreadMin[i];
167  minimumIdx = m_ThreadMinIndex[i];
168  }
169  if (m_ThreadMax[i] > maximum)
170  {
171  maximum = m_ThreadMax[i];
172  maximumIdx = m_ThreadMaxIndex[i];
173  }
174  }
175 
176  // Set the outputs
177  this->GetMinimumOutput()->Set(minimum);
178  this->GetMaximumOutput()->Set(maximum);
179  this->GetMinimumIndexOutput()->Set(minimumIdx);
180  this->GetMaximumIndexOutput()->Set(maximumIdx);
181 }
182 
183 template <class TInputImage>
185 {
186  int numberOfThreads = this->GetNumberOfWorkUnits();
187 
188  m_ThreadMin.resize(numberOfThreads);
189  m_ThreadMax.resize(numberOfThreads);
190  std::fill(m_ThreadMin.begin(), m_ThreadMin.end(), itk::NumericTraits<PixelType>::max());
191  std::fill(m_ThreadMax.begin(), m_ThreadMax.end(), itk::NumericTraits<PixelType>::NonpositiveMin());
192 
193  IndexType zeroIdx;
194  zeroIdx.Fill(0);
195  m_ThreadMinIndex.resize(numberOfThreads);
196  m_ThreadMaxIndex.resize(numberOfThreads);
197  std::fill(m_ThreadMinIndex.begin(), m_ThreadMinIndex.end(), zeroIdx);
198  std::fill(m_ThreadMaxIndex.begin(), m_ThreadMaxIndex.end(), zeroIdx);
199 }
200 
201 template <class TInputImage>
202 void PersistentMinMaxImageFilter<TInputImage>::ThreadedGenerateData(const RegionType& outputRegionForThread, itk::ThreadIdType threadId)
203 {
204  // support progress methods/callbacks
205  itk::ProgressReporter progress(this, threadId, outputRegionForThread.GetNumberOfPixels());
206 
207  InputImagePointer inputPtr = const_cast<TInputImage*>(this->GetInput(0));
208  itk::ImageRegionConstIterator<TInputImage> it(inputPtr, outputRegionForThread);
209  it.GoToBegin();
210  // do the work
211  while (!it.IsAtEnd())
212  {
213  PixelType value = it.Get();
214  if (value < m_ThreadMin[threadId])
215  {
216  m_ThreadMin[threadId] = value;
217  m_ThreadMinIndex[threadId] = it.GetIndex();
218  }
219  if (value > m_ThreadMax[threadId])
220  {
221  m_ThreadMax[threadId] = value;
222  m_ThreadMaxIndex[threadId] = it.GetIndex();
223  }
224  ++it;
225  progress.CompletedPixel();
226  }
227 }
228 
229 template <class TImage>
230 void PersistentMinMaxImageFilter<TImage>::PrintSelf(std::ostream& os, itk::Indent indent) const
231 {
232  Superclass::PrintSelf(os, indent);
233 
234  os << indent << "Minimum: " << static_cast<typename itk::NumericTraits<PixelType>::PrintType>(this->GetMinimum()) << std::endl;
235  os << indent << "Maximum: " << static_cast<typename itk::NumericTraits<PixelType>::PrintType>(this->GetMaximum()) << std::endl;
236  os << indent << "Minimum Index: " << this->GetMinimumIndex() << std::endl;
237  os << indent << "Maximum Index: " << this->GetMaximumIndex() << std::endl;
238 }
239 
240 } // end namespace otb
241 #endif
void PrintSelf(std::ostream &os, itk::Indent indent) const override
itk::ProcessObject::DataObjectPointerArraySizeType DataObjectPointerArraySizeType
itk::SimpleDataObjectDecorator< PixelType > PixelObjectType
itk::SimpleDataObjectDecorator< IndexType > IndexObjectType
DataObjectPointer MakeOutput(DataObjectPointerArraySizeType idx) override
void ThreadedGenerateData(const RegionType &outputRegionForThread, itk::ThreadIdType threadId) override
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.