OTB  6.7.0
Orfeo Toolbox
otbStreamingImageVirtualWriter.hxx
Go to the documentation of this file.
1 /*
2  * Copyright (C) 2005-2019 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_hxx
22 #define otbStreamingImageVirtualWriter_hxx
24 
25 #include "otbMacro.h"
26 #include "itkCommand.h"
27 
35 #include "otbUtils.h"
36 
37 namespace otb
38 {
39 
40 template <class TInputImage>
43  : m_NumberOfDivisions(0),
44  m_CurrentDivision(0),
45  m_DivisionProgress(0.0),
46  m_IsObserving(true),
47  m_ObserverID(0)
48 {
49  // By default, we use tiled streaming, with automatic tile size
50  // We don't set any parameter, so the memory size is retrieved from the OTB configuration options
51  this->SetAutomaticAdaptativeStreaming();
52 }
53 
54 template <class TInputImage>
57 {
58 }
59 
60 template <class TInputImage>
61 void
63 ::PrintSelf(std::ostream& os, itk::Indent indent) const
64 {
65  Superclass::PrintSelf(os, indent);
66 }
67 
68 template <class TInputImage>
69 void
71 ::SetNumberOfDivisionsStrippedStreaming(unsigned int nbDivisions)
72 {
73  typedef NumberOfDivisionsStrippedStreamingManager<TInputImage> NumberOfDivisionsStrippedStreamingManagerType;
74  typename NumberOfDivisionsStrippedStreamingManagerType::Pointer streamingManager = NumberOfDivisionsStrippedStreamingManagerType::New();
75  streamingManager->SetNumberOfDivisions(nbDivisions);
76 
77  m_StreamingManager = streamingManager;
78 }
79 
80 template <class TInputImage>
81 void
83 ::SetNumberOfDivisionsTiledStreaming(unsigned int nbDivisions)
84 {
85  typedef NumberOfDivisionsTiledStreamingManager<TInputImage> NumberOfDivisionsTiledStreamingManagerType;
86  typename NumberOfDivisionsTiledStreamingManagerType::Pointer streamingManager = NumberOfDivisionsTiledStreamingManagerType::New();
87  streamingManager->SetNumberOfDivisions(nbDivisions);
88 
89  m_StreamingManager = streamingManager;
90 }
91 
92 template <class TInputImage>
93 void
95 ::SetNumberOfLinesStrippedStreaming(unsigned int nbLinesPerStrip)
96 {
97  typedef NumberOfLinesStrippedStreamingManager<TInputImage> NumberOfLinesStrippedStreamingManagerType;
98  typename NumberOfLinesStrippedStreamingManagerType::Pointer streamingManager = NumberOfLinesStrippedStreamingManagerType::New();
99  streamingManager->SetNumberOfLinesPerStrip(nbLinesPerStrip);
100 
101  m_StreamingManager = streamingManager;
102 }
103 
104 template <class TInputImage>
105 void
107 ::SetAutomaticStrippedStreaming(unsigned int availableRAM, double bias)
108 {
109  typedef RAMDrivenStrippedStreamingManager<TInputImage> RAMDrivenStrippedStreamingManagerType;
110  typename RAMDrivenStrippedStreamingManagerType::Pointer streamingManager = RAMDrivenStrippedStreamingManagerType::New();
111  streamingManager->SetAvailableRAMInMB(availableRAM);
112  streamingManager->SetBias(bias);
113  m_StreamingManager = streamingManager;
114 }
115 
116 template <class TInputImage>
117 void
119 ::SetTileDimensionTiledStreaming(unsigned int tileDimension)
120 {
121  typedef TileDimensionTiledStreamingManager<TInputImage> TileDimensionTiledStreamingManagerType;
122  typename TileDimensionTiledStreamingManagerType::Pointer streamingManager = TileDimensionTiledStreamingManagerType::New();
123  streamingManager->SetTileDimension(tileDimension);
124 
125  m_StreamingManager = streamingManager;
126 }
127 
128 template <class TInputImage>
129 void
131 ::SetAutomaticTiledStreaming(unsigned int availableRAM, double bias)
132 {
133  typedef RAMDrivenTiledStreamingManager<TInputImage> RAMDrivenTiledStreamingManagerType;
134  typename RAMDrivenTiledStreamingManagerType::Pointer streamingManager = RAMDrivenTiledStreamingManagerType::New();
135  streamingManager->SetAvailableRAMInMB(availableRAM);
136  streamingManager->SetBias(bias);
137  m_StreamingManager = streamingManager;
138 }
139 
140 template <class TInputImage>
141 void
143 ::SetAutomaticAdaptativeStreaming(unsigned int availableRAM, double bias)
144 {
145  typedef RAMDrivenAdaptativeStreamingManager<TInputImage> RAMDrivenAdaptativeStreamingManagerType;
146  typename RAMDrivenAdaptativeStreamingManagerType::Pointer streamingManager = RAMDrivenAdaptativeStreamingManagerType::New();
147  streamingManager->SetAvailableRAMInMB(availableRAM);
148  streamingManager->SetBias(bias);
149  m_StreamingManager = streamingManager;
150 }
151 
152 template <class TInputImage>
153 void
156 {
157  InputImagePointer inputPtr = const_cast<InputImageType *>(this->GetInput(0));
158  inputPtr->UpdateOutputInformation();
159 
160  this->GenerateData();
161 }
162 
163 template <class TInputImage>
164 void
167 {
168  InputImagePointer inputPtr = const_cast<InputImageType *>(this->GetInput(0));
169 
170  InputImageRegionType region;
171  typename InputImageRegionType::SizeType size;
172  typename InputImageRegionType::IndexType index;
173 
174  index.Fill(0);
175  size.Fill(0);
176  region.SetSize(size);
177  region.SetIndex(index);
178  inputPtr->SetRequestedRegion(region);
179 }
180 
181 template<class TInputImage>
182 void
185 {
186  otb::Logger::Instance()->LogSetupInformation();
187 
191  this->PrepareOutputs();
192  this->SetAbortGenerateData(0);
193  this->SetProgress(0.0);
194  this->m_Updating = true;
195 
199  this->InvokeEvent(itk::StartEvent());
200 
204  InputImagePointer inputPtr = const_cast<InputImageType *>(this->GetInput(0));
205  InputImageRegionType outputRegion = inputPtr->GetLargestPossibleRegion();
206 
212  m_StreamingManager->PrepareStreaming(inputPtr, outputRegion);
213  m_NumberOfDivisions = m_StreamingManager->GetNumberOfSplits();
215 
219  // Get the source process object
220  itk::ProcessObject* source = inputPtr->GetSource();
221  m_IsObserving = false;
222  m_ObserverID = 0;
223 
224  // Check if source exists
225  if(source)
226  {
227  typedef itk::MemberCommand<Self> CommandType;
228  typedef typename CommandType::Pointer CommandPointerType;
229 
230  CommandPointerType command = CommandType::New();
231  command->SetCallbackFunction(this, &Self::ObserveSourceFilterProgress);
232 
233  m_ObserverID = source->AddObserver(itk::ProgressEvent(), command);
234  m_IsObserving = true;
235  }
236 
237  const auto firstSplitSize = m_StreamingManager->GetSplit(0).GetSize();
238  otbLogMacro(Info,<<"Estimation will be performed in "<<m_NumberOfDivisions<<" blocks of "<<firstSplitSize[0]<<"x"<<firstSplitSize[1]<<" pixels");
239 
240 
245  InputImageRegionType streamRegion;
246  for (m_CurrentDivision = 0;
247  m_CurrentDivision < m_NumberOfDivisions && !this->GetAbortGenerateData();
248  m_CurrentDivision++, m_DivisionProgress = 0, this->UpdateFilterProgress())
249  {
250  streamRegion = m_StreamingManager->GetSplit(m_CurrentDivision);
251  //inputPtr->ReleaseData();
252  //inputPtr->SetRequestedRegion(streamRegion);
253  //inputPtr->Update();
254  inputPtr->SetRequestedRegion(streamRegion);
255  inputPtr->PropagateRequestedRegion();
256  inputPtr->UpdateOutputData();
257  }
259 
264  if (!this->GetAbortGenerateData())
265  {
266  this->UpdateProgress(1.0);
267  }
268  else
269  {
270  itk::ProcessAborted e(__FILE__, __LINE__);
271  e.SetLocation(ITK_LOCATION);
272  e.SetDescription("Image streaming has been aborted");
273  throw e;
274  }
276 
277  // Notify end event observers
278  this->InvokeEvent(itk::EndEvent());
279 
280  if (m_IsObserving)
281  {
282  m_IsObserving = false;
283  source->RemoveObserver(m_ObserverID);
284  }
285 
289  for (unsigned int idx = 0; idx < this->GetNumberOfOutputs(); ++idx)
290  {
291  if (this->GetOutput(idx))
292  {
293  this->GetOutput(idx)->DataHasBeenGenerated();
294  }
295  }
297 
301  this->ReleaseInputs();
302 }
303 
304 template <class TInputImage>
305 const bool &
308 {
309  m_Lock.Lock();
310  bool ret = Superclass::GetAbortGenerateData();
311  m_Lock.Unlock();
312  if (ret) return otb::Utils::TrueConstant;
314 }
315 
316 template <class TInputImage>
317 void
320 {
321  m_Lock.Lock();
322  Superclass::SetAbortGenerateData(val);
323  m_Lock.Unlock();
324 }
325 
326 } // end namespace otb
327 
328 #endif
itk::Size< Monteverdi_DIMENSION > SizeType
Definition: mvdTypes.h:137
This class computes the divisions needed to stream an image according to the input image tiling schem...
void SetNumberOfDivisionsStrippedStreaming(unsigned int nbDivisions)
This class computes the divisions needed to stream an image by strips, driven by a user-defined numbe...
void PrintSelf(std::ostream &os, itk::Indent indent) const override
virtual void SetDescription(const std::string &s)
virtual void SetLocation(const std::string &s)
void SetNumberOfDivisionsTiledStreaming(unsigned int nbDivisions)
void SetNumberOfLinesStrippedStreaming(unsigned int nbLinesPerStrip)
This class computes the divisions needed to stream an image by strips, driven by a user-defined desir...
void RemoveObserver(unsigned long tag)
OTBCommon_EXPORT bool const TrueConstant
This class computes the divisions needed to stream an image by strips, according to a user-defined av...
void SetAutomaticTiledStreaming(unsigned int availableRAM=0, double bias=1.0)
OTBCommon_EXPORT bool const FalseConstant
unsigned long AddObserver(const EventObject &event, Command *)
This class computes the divisions needed to stream an image in square tiles, according to a user-defi...
itk::Index< Monteverdi_DIMENSION > IndexType
Definition: mvdTypes.h:133
void SetAutomaticStrippedStreaming(unsigned int availableRAM=0, double bias=1.0)
#define otbLogMacro(level, msg)
Definition: otbMacro.h:54
void SetTileDimensionTiledStreaming(unsigned int tileDimension)
void SetAutomaticAdaptativeStreaming(unsigned int availableRAM=0, double bias=1.0)
This class computes the divisions needed to stream an image by strips, driven by a user-defined numbe...
void SetAbortGenerateData(const bool val) override
static Logger * Instance()
This class computes the divisions needed to stream an image in square tiles, driven by a user-defined...
const bool & GetAbortGenerateData() const override