OTB  10.0.0
Orfeo Toolbox
otbStreamingHistogramVectorImageFilter.h
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_h
23 #define otbStreamingHistogramVectorImageFilter_h
24 
27 
28 #include "otbObjectList.h"
29 #include "itkStatisticsAlgorithm.h"
30 #include "itkNumericTraits.h"
31 #include "itkHistogram.h"
32 
33 namespace otb
34 {
35 
54 template <class TInputImage>
55 class ITK_EXPORT PersistentHistogramVectorImageFilter : public PersistentImageFilter<TInputImage, TInputImage>
56 {
57 public:
61  typedef itk::SmartPointer<Self> Pointer;
62  typedef itk::SmartPointer<const Self> ConstPointer;
63 
65  itkNewMacro(Self);
66 
69 
71  typedef TInputImage ImageType;
72  typedef typename TInputImage::Pointer InputImagePointer;
73  typedef typename TInputImage::RegionType RegionType;
74  typedef typename TInputImage::SizeType SizeType;
75  typedef typename TInputImage::IndexType IndexType;
76  typedef typename TInputImage::PixelType PixelType;
77  typedef typename TInputImage::InternalPixelType InternalPixelType;
78 
79  itkStaticConstMacro(InputImageDimension, unsigned int, TInputImage::ImageDimension);
80 
82  itkStaticConstMacro(ImageDimension, unsigned int, TInputImage::ImageDimension);
83 
85  typedef typename itk::NumericTraits<InternalPixelType>::RealType RealType;
86  typedef itk::VariableLengthVector<RealType> RealPixelType;
87 
89  typedef typename itk::DataObject::Pointer DataObjectPointer;
90  typedef itk::ProcessObject::DataObjectPointerArraySizeType DataObjectPointerArraySizeType;
91 
93  typedef itk::Statistics::DenseFrequencyContainer2 DFContainerType;
94 
95  typedef typename itk::NumericTraits<InternalPixelType>::RealType HistogramMeasurementRealType;
96 
97  typedef itk::Statistics::Histogram<HistogramMeasurementRealType, DFContainerType> HistogramType;
98 
99  typedef itk::VariableLengthVector<unsigned int> CountVectorType;
100 
104  typedef typename std::vector<HistogramListPointerType> ArrayHistogramListType;
105 
106 
111 
115  itkGetConstReferenceMacro(NoDataValue, InternalPixelType);
116 
120  itkSetMacro(NoDataFlag, bool);
121 
125  itkGetMacro(NoDataFlag, bool);
126 
130  itkBooleanMacro(NoDataFlag);
131 
132  inline void SetNumberOfBins(unsigned int i, CountVectorType::ValueType size)
133  {
134  m_Size[i] = size;
135  }
136 
137  inline void SetNumberOfBins(const CountVectorType& size)
138  {
139  m_Size = size;
140  }
141 
143  HistogramListType* GetHistogramListOutput();
144  const HistogramListType* GetHistogramListOutput() const;
146 
148  itkGetConstReferenceMacro(HistogramMin, MeasurementVectorType);
149 
151  itkSetMacro(HistogramMin, MeasurementVectorType);
152 
154  itkGetConstReferenceMacro(HistogramMax, MeasurementVectorType);
155 
157  itkSetMacro(HistogramMax, MeasurementVectorType);
158 
160  itkSetMacro(SubSamplingRate, unsigned int);
161 
163  itkGetMacro(SubSamplingRate, unsigned int);
164 
168  DataObjectPointer MakeOutput(DataObjectPointerArraySizeType idx) override;
169  using Superclass::MakeOutput;
170 
174  void AllocateOutputs() override;
175  void GenerateOutputInformation() override;
176  void Synthetize(void) override;
177  void Reset(void) override;
179 
180 protected:
183  {
184  }
185  void PrintSelf(std::ostream& os, itk::Indent indent) const override;
187  void ThreadedGenerateData(const RegionType& outputRegionForThread, itk::ThreadIdType threadId) override;
188 
189 private:
191  void operator=(const Self&) = delete;
192 
199 
201  unsigned int m_SubSamplingRate;
202 
203 }; // end of class PersistentStatisticsVectorImageFilter
204 
227 template <class TInputImage>
228 class ITK_EXPORT StreamingHistogramVectorImageFilter : public PersistentFilterStreamingDecorator<PersistentHistogramVectorImageFilter<TInputImage>>
229 {
230 public:
234  typedef itk::SmartPointer<Self> Pointer;
235  typedef itk::SmartPointer<const Self> ConstPointer;
236 
238  itkNewMacro(Self);
239 
242 
243  typedef TInputImage InputImageType;
245 
249 
250  using Superclass::SetInput;
252  {
253  this->GetFilter()->SetInput(input);
254  }
256  {
257  return this->GetFilter()->GetInput();
258  }
259 
262  {
263  return this->GetFilter()->GetHistogramListOutput();
264  }
265 
266 
267 protected:
270 
273  {
274  }
275 
276 private:
278  void operator=(const Self&) = delete;
279 };
280 
281 } // end namespace otb
282 
283 #ifndef OTB_MANUAL_INSTANTIATION
285 #endif
286 
287 #endif
This class is a generic all-purpose wrapping around an std::vector<itk::SmartPointer<ObjectType> >.
Definition: otbObjectList.h:41
itk::SmartPointer< Self > Pointer
Definition: otbObjectList.h:46
This filter link a persistent filter with a StreamingImageVirtualWriter.
Compute the histogram of a large image using streaming.
itk::NumericTraits< InternalPixelType >::RealType RealType
itk::Statistics::DenseFrequencyContainer2 DFContainerType
itk::ProcessObject::DataObjectPointerArraySizeType DataObjectPointerArraySizeType
PersistentImageFilter< TInputImage, TInputImage > Superclass
itk::VariableLengthVector< unsigned int > CountVectorType
void operator=(const Self &)=delete
itk::Statistics::Histogram< HistogramMeasurementRealType, DFContainerType > HistogramType
std::vector< HistogramListPointerType > ArrayHistogramListType
itk::NumericTraits< InternalPixelType >::RealType HistogramMeasurementRealType
void SetNumberOfBins(unsigned int i, CountVectorType::ValueType size)
PersistentHistogramVectorImageFilter(const Self &)=delete
This filter is the base class for all filter persisting data through multiple update....
This class streams the whole input image through the PersistentHistogramVectorImageFilter.
PersistentFilterStreamingDecorator< PersistentHistogramVectorImageFilter< TInputImage > > Superclass
StreamingHistogramVectorImageFilter(const Self &)=delete
void operator=(const Self &)=delete
OTBMetadata_EXPORT char const * NoDataValue
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.