OTB  10.0.0
Orfeo Toolbox
otbStreamingImageVirtualWriter.h
Go to the documentation of this file.
1 /*
2  * Copyright (C) 2005-2024 Centre National d'Etudes Spatiales (CNES)
3  *
4  * This file is part of Orfeo Toolbox
5  *
6  * https://www.orfeo-toolbox.org/
7  *
8  * Licensed under the Apache License, Version 2.0 (the "License");
9  * you may not use this file except in compliance with the License.
10  * You may obtain a copy of the License at
11  *
12  * http://www.apache.org/licenses/LICENSE-2.0
13  *
14  * Unless required by applicable law or agreed to in writing, software
15  * distributed under the License is distributed on an "AS IS" BASIS,
16  * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
17  * See the License for the specific language governing permissions and
18  * limitations under the License.
19  */
20 
21 #ifndef otbStreamingImageVirtualWriter_h
22 #define otbStreamingImageVirtualWriter_h
23 
24 #include "itkMacro.h"
25 #include "itkImageToImageFilter.h"
26 #include "otbStreamingManager.h"
27 #include <mutex>
28 
29 namespace otb
30 {
31 
59 template <class TInputImage>
60 class ITK_EXPORT StreamingImageVirtualWriter : public itk::ImageToImageFilter<TInputImage, TInputImage>
61 {
62 public:
65  typedef itk::ImageToImageFilter<TInputImage, TInputImage> Superclass;
66  typedef itk::SmartPointer<Self> Pointer;
67  typedef itk::SmartPointer<const Self> ConstPointer;
68 
70  itkNewMacro(Self);
71 
73  itkTypeMacro(StreamingImageVirtualWriter, itk::ImageToImageFilter);
74 
76  typedef TInputImage InputImageType;
77  typedef typename InputImageType::Pointer InputImagePointer;
78  typedef typename InputImageType::RegionType InputImageRegionType;
79  typedef typename InputImageType::PixelType InputImagePixelType;
80 
84 
86  itkStaticConstMacro(InputImageDimension, unsigned int, InputImageType::ImageDimension);
87 
91  {
92  return m_StreamingManager;
93  }
94 
97  void SetStreamingManager(StreamingManagerType* streamingManager)
98  {
99  m_StreamingManager = streamingManager;
100  }
101 
104  void SetNumberOfDivisionsStrippedStreaming(unsigned int nbDivisions);
105 
108  void SetNumberOfDivisionsTiledStreaming(unsigned int nbDivisions);
109 
112  void SetNumberOfLinesStrippedStreaming(unsigned int nbLinesPerStrip);
113 
119  void SetAutomaticStrippedStreaming(unsigned int availableRAM = 0, double bias = 1.0);
120 
123  void SetTileDimensionTiledStreaming(unsigned int tileDimension);
124 
131  void SetAutomaticTiledStreaming(unsigned int availableRAM = 0, double bias = 1.0);
132 
139  void SetAutomaticAdaptativeStreaming(unsigned int availableRAM = 0, double bias = 1.0);
140 
143  void Update() override;
144 
146  const bool& GetAbortGenerateData() const override;
147 
148  void SetAbortGenerateData(const bool val) override;
149 
150 protected:
152 
153  ~StreamingImageVirtualWriter() override;
154 
155  void PrintSelf(std::ostream& os, itk::Indent indent) const override;
156 
157  void GenerateData(void) override;
158 
159  void GenerateInputRequestedRegion(void) override;
160 
161 private:
163  void operator=(const StreamingImageVirtualWriter&) = delete;
164 
165  void ObserveSourceFilterProgress(itk::Object* object, const itk::EventObject& event)
166  {
167  if (typeid(event) != typeid(itk::ProgressEvent))
168  {
169  return;
170  }
171 
172  itk::ProcessObject* processObject = dynamic_cast<itk::ProcessObject*>(object);
173  if (processObject)
174  {
175  m_DivisionProgress = processObject->GetProgress();
176  }
177 
178  this->UpdateFilterProgress();
179  }
180 
182  {
183  this->UpdateProgress((m_DivisionProgress + m_CurrentDivision) / m_NumberOfDivisions);
184  }
185 
186  unsigned int m_NumberOfDivisions;
187  unsigned int m_CurrentDivision;
189 
191 
193  unsigned long m_ObserverID;
194 
196  mutable std::mutex m_Lock;
197 
198 };
199 
200 } // end namespace otb
201 
202 #ifndef OTB_MANUAL_INSTANTIATION
204 #endif
205 
206 #endif
This class acts like a StreamingImageFileWriter, but without actually writing data to the disk.
void SetStreamingManager(StreamingManagerType *streamingManager)
StreamingManager< InputImageType > StreamingManagerType
itk::ImageToImageFilter< TInputImage, TInputImage > Superclass
itk::SmartPointer< const Self > ConstPointer
StreamingManagerType * GetStreamingManager(void)
StreamingImageVirtualWriter(const StreamingImageVirtualWriter &)=delete
StreamingManagerType::Pointer StreamingManagerPointerType
void operator=(const StreamingImageVirtualWriter &)=delete
void ObserveSourceFilterProgress(itk::Object *object, const itk::EventObject &event)
This class handles the streaming process used in the writers implementation.
itk::SmartPointer< Self > Pointer
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.