OTB  7.2.0
Orfeo Toolbox
otbStreamingMinMaxVectorImageFilter.h
Go to the documentation of this file.
1 /*
2  * Copyright (C) 1999-2011 Insight Software Consortium
3  * Copyright (C) 2005-2020 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 otbStreamingMinMaxVectorImageFilter_h
23 #define otbStreamingMinMaxVectorImageFilter_h
24 
27 #include "itkSimpleDataObjectDecorator.h"
28 #include "itkImageRegionSplitter.h"
29 #include "itkVariableSizeMatrix.h"
30 #include "itkVariableLengthVector.h"
31 
32 namespace otb
33 {
34 
53 template <class TInputImage>
54 class ITK_EXPORT PersistentMinMaxVectorImageFilter : public PersistentImageFilter<TInputImage, TInputImage>
55 {
56 public:
60  typedef itk::SmartPointer<Self> Pointer;
61  typedef itk::SmartPointer<const Self> ConstPointer;
62 
64  itkNewMacro(Self);
65 
68 
70  typedef TInputImage ImageType;
71  typedef typename TInputImage::Pointer InputImagePointer;
72  typedef typename TInputImage::RegionType RegionType;
73  typedef typename TInputImage::SizeType SizeType;
75  typedef typename TInputImage::PixelType PixelType;
76  typedef typename TInputImage::InternalPixelType InternalPixelType;
77 
78  itkStaticConstMacro(InputImageDimension, unsigned int, TInputImage::ImageDimension);
79 
81  itkStaticConstMacro(ImageDimension, unsigned int, TInputImage::ImageDimension);
82 
85  typedef itk::VariableLengthVector<RealType> RealPixelType;
86 
88  typedef typename itk::DataObject::Pointer DataObjectPointer;
89  typedef itk::ProcessObject::DataObjectPointerArraySizeType DataObjectPointerArraySizeType;
90 
92  typedef typename itk::VariableSizeMatrix<RealType> MatrixType;
93  typedef typename std::vector<MatrixType> ArrayMatrixType;
94  typedef typename itk::Array<long> ArrayLongPixelType;
95  typedef typename std::vector<RealPixelType> ArrayRealPixelType;
96  typedef typename std::vector<PixelType> ArrayPixelType;
97  typedef itk::SimpleDataObjectDecorator<RealPixelType> RealPixelObjectType;
98  typedef itk::SimpleDataObjectDecorator<PixelType> PixelObjectType;
99  typedef itk::SimpleDataObjectDecorator<MatrixType> MatrixObjectType;
100 
104  itkSetMacro(NoDataValue, InternalPixelType);
105 
109  itkGetConstReferenceMacro(NoDataValue, InternalPixelType);
110 
114  itkSetMacro(NoDataFlag, bool);
115 
119  itkGetMacro(NoDataFlag, bool);
120 
124  itkBooleanMacro(NoDataFlag);
125 
126 
128  PixelType GetMinimum() const
129  {
130  return this->GetMinimumOutput()->Get();
131  }
132  PixelObjectType* GetMinimumOutput();
133  const PixelObjectType* GetMinimumOutput() const;
135 
137  PixelType GetMaximum() const
138  {
139  return this->GetMaximumOutput()->Get();
140  }
141  PixelObjectType* GetMaximumOutput();
142  const PixelObjectType* GetMaximumOutput() const;
144 
148  DataObjectPointer MakeOutput(DataObjectPointerArraySizeType idx) override;
149  using Superclass::MakeOutput;
150 
154  void AllocateOutputs() override;
155  void GenerateOutputInformation() override;
156  void Synthetize(void) override;
157  void Reset(void) override;
159 
160 protected:
163  {
164  }
165  void PrintSelf(std::ostream& os, itk::Indent indent) const override;
167  void ThreadedGenerateData(const RegionType& outputRegionForThread, itk::ThreadIdType threadId) override;
168 
169 private:
170  PersistentMinMaxVectorImageFilter(const Self&) = delete;
171  void operator=(const Self&) = delete;
172 
173  ArrayPixelType m_ThreadMin;
174  ArrayPixelType m_ThreadMax;
176  InternalPixelType m_NoDataValue;
177 
178 }; // end of class PersistentStatisticsVectorImageFilter
179 
202 template <class TInputImage>
203 class ITK_EXPORT StreamingMinMaxVectorImageFilter : public PersistentFilterStreamingDecorator<PersistentMinMaxVectorImageFilter<TInputImage>>
204 {
205 public:
209  typedef itk::SmartPointer<Self> Pointer;
210  typedef itk::SmartPointer<const Self> ConstPointer;
211 
213  itkNewMacro(Self);
214 
217 
218  typedef TInputImage InputImageType;
228 
233 
234  using Superclass::SetInput;
235  void SetInput(InputImageType* input)
236  {
237  this->GetFilter()->SetInput(input);
238  }
239  const InputImageType* GetInput()
240  {
241  return this->GetFilter()->GetInput();
242  }
243 
245  PixelType GetMinimum() const
246  {
247  return this->GetFilter()->GetMinimumOutput()->Get();
248  }
249  PixelObjectType* GetMinimumOutput()
250  {
251  return this->GetFilter()->GetMinimumOutput();
252  }
253  const PixelObjectType* GetMinimumOutput() const
254  {
255  return this->GetFilter()->GetMinimumOutput();
256  }
258 
260  PixelType GetMaximum() const
261  {
262  return this->GetFilter()->GetMaximumOutput()->Get();
263  }
264  PixelObjectType* GetMaximumOutput()
265  {
266  return this->GetFilter()->GetMaximumOutput();
267  }
268  const PixelObjectType* GetMaximumOutput() const
269  {
270  return this->GetFilter()->GetMaximumOutput();
271  }
273 
274 protected:
277 
280  {
281  }
282 
283 private:
284  StreamingMinMaxVectorImageFilter(const Self&) = delete;
285  void operator=(const Self&) = delete;
286 };
287 
288 } // end namespace otb
289 
290 #ifndef OTB_MANUAL_INSTANTIATION
292 #endif
293 
294 #endif
itk::Size< Monteverdi_DIMENSION > SizeType
Definition: mvdTypes.h:137
itk::Index< Monteverdi_DIMENSION > IndexType
Definition: mvdTypes.h:133
StatFilterType::RealPixelObjectType RealPixelObjectType
PersistentFilterStreamingDecorator< PersistentMinMaxVectorImageFilter< TInputImage > > Superclass
PersistentImageFilter< TInputImage, TInputImage > Superclass
itk::SimpleDataObjectDecorator< MatrixType > MatrixObjectType
itk::NumericTraits< InternalPixelType >::RealType RealType
This class streams the whole input image through the PersistentMinMaxVectorImageFilter.
itk::SimpleDataObjectDecorator< PixelType > PixelObjectType
itk::ProcessObject::DataObjectPointerArraySizeType DataObjectPointerArraySizeType
Monteverdi_FLOATING_TYPE RealType
Definition: mvdTypes.h:84
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.
Compute min. max of a large image using streaming.
OTBMetadata_EXPORT char const * NoDataValue
itk::VariableLengthVector< RealType > RealPixelType
itk::SimpleDataObjectDecorator< RealPixelType > RealPixelObjectType
This filter is the base class for all filter persisting data through multiple update. For instance, a filter computing global statistics on an image with streaming capabilities will have to keep the temporary results for each streamed piece of the image in order to synthesize the global statistics at the end. This filter is an itk::ImageToImageFilter, providing two additional methods. The first one, Synthetize(), allows the user to synthesize temporary data produced by the multiple updates on different pieces of the image to the global result. The second one, Reset(), allows the user to reset the temporary data for a new input image for instance.
This filter link a persistent filter with a StreamingImageVirtualWriter.