OTB  10.0.0
Orfeo Toolbox
otbMultiImageFileWriter.h
Go to the documentation of this file.
1 /*
2  * Copyright (C) 2017-2019 CS Systemes d'Information (CS SI)
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 otbMultiImageFileWriter_h
22 #define otbMultiImageFileWriter_h
23 
24 #include "otbImageFileWriter.h"
25 #include "otbImage.h"
26 #include "itkImageBase.h"
27 #include "itkProcessObject.h"
28 #include "itkImageIOBase.h"
29 #include "OTBImageIOExport.h"
30 
31 #include <memory>
32 
33 namespace otb
34 {
35 
47 class OTBImageIO_EXPORT MultiImageFileWriter : public itk::ProcessObject
48 {
49 public:
52  typedef itk::ProcessObject Superclass;
53  typedef itk::SmartPointer<Self> Pointer;
54  typedef itk::SmartPointer<const Self> ConstPointer;
55 
56  itkNewMacro(Self);
57 
58  itkTypeMacro(MultiImageFileWriter, itk::ProcessObject);
59 
61  typedef itk::ImageBase<2> ImageBaseType;
62  typedef ImageBaseType::RegionType RegionType;
63  typedef ImageBaseType::IndexType IndexType;
64  typedef ImageBaseType::IndexValueType IndexValueType;
65  typedef ImageBaseType::SizeType SizeType;
66  typedef ImageBaseType::SizeValueType SizeValueType;
67 
70  void SetNumberOfDivisionsStrippedStreaming(unsigned int nbDivisions);
71 
74  void SetNumberOfDivisionsTiledStreaming(unsigned int nbDivisions);
75 
79  void SetNumberOfLinesStrippedStreaming(unsigned int nbLinesPerStrip);
80 
90  void SetAutomaticStrippedStreaming(unsigned int availableRAM = 0, double bias = 1.0);
91 
94  void SetTileDimensionTiledStreaming(unsigned int tileDimension);
95 
106  void SetAutomaticTiledStreaming(unsigned int availableRAM = 0, double bias = 1.0);
107 
114  void SetAutomaticAdaptativeStreaming(unsigned int availableRAM = 0, double bias = 1.0);
115 
116  virtual void UpdateOutputData(itk::DataObject* itkNotUsed(output)) override;
117 
128  template <class TImage>
129  void AddInputImage(const TImage* inputPtr, const std::string& fileName);
130 
131 
135  template <class TWriter>
136  void AddInputWriter(typename TWriter::Pointer writer);
137 
138  virtual void UpdateOutputInformation() override;
139 
140  virtual void Update() override
141  {
142  this->UpdateOutputInformation();
143  this->UpdateOutputData(NULL);
144  }
145 
146 protected:
149  {
150  }
151 
153  virtual void GenerateData(void) override;
154 
159  virtual void GenerateInputRequestedRegion() override;
160 
162  virtual void InitializeStreaming();
163 
168 
170  virtual RegionType GetStreamRegion(int inputIndex);
171 
172  void operator=(const MultiImageFileWriter&) = delete;
173 
174  void ObserveSourceFilterProgress(itk::Object* object, const itk::EventObject& event)
175  {
176  if (typeid(event) != typeid(itk::ProgressEvent))
177  {
178  return;
179  }
180 
181  itk::ProcessObject* processObject = dynamic_cast<itk::ProcessObject*>(object);
182  if (processObject)
183  {
184  m_DivisionProgress = processObject->GetProgress();
185  }
186 
187  this->UpdateFilterProgress();
188  }
189 
191  {
192  this->UpdateProgress((m_DivisionProgress + m_CurrentDivision) / m_NumberOfDivisions);
193  }
194 
195 private:
200 
201  // Division parameters
202  unsigned int m_NumberOfDivisions;
203  unsigned int m_CurrentDivision;
205 
207  unsigned long m_ObserverID;
208 
214  class SinkBase
215  {
216  public:
218  {
219  }
220  SinkBase(ImageBaseType::ConstPointer inputImage) : m_InputImage(inputImage)
221  {
222  }
223  virtual ~SinkBase()
224  {
225  }
226  virtual ImageBaseType::ConstPointer GetInput() const
227  {
228  return m_InputImage;
229  }
230  virtual ImageBaseType::Pointer GetInput()
231  {
232  return const_cast<ImageBaseType*>(m_InputImage.GetPointer());
233  }
234  virtual void WriteImageInformation() = 0;
235  virtual void Write(const RegionType& streamRegion) = 0;
236  virtual bool CanStreamWrite() const = 0;
237  typedef std::shared_ptr<SinkBase> Pointer;
238 
239 
240  virtual itk::ImageRegion<2> GetRegionToWrite() const = 0;
241 
242  protected:
244  ImageBaseType::ConstPointer m_InputImage;
245  };
246 
252  template <class TImage>
253  class Sink : public SinkBase
254  {
255  public:
257  {
258  }
259  Sink(typename TImage::ConstPointer inputImage, const std::string& filename);
261 
262  virtual ~Sink()
263  {
264  }
265 
266  void WriteImageInformation() override;
267  void Write(const RegionType& streamRegion) override;
268  bool CanStreamWrite() const override;
269  typedef std::shared_ptr<Sink> Pointer;
270 
274  itk::ImageRegion<2> GetRegionToWrite() const override;
275 
276  private:
279 
282  };
283 
285  typedef std::vector<std::shared_ptr<SinkBase>> SinkListType;
287 
288  std::vector<RegionType> m_StreamRegionList;
289 };
290 
291 } // end of namespace otb
292 
293 #ifndef OTB_MANUAL_INSTANTIATION
295 #endif
296 
297 #endif // otbMultiImageFileWriter_h
itk::SmartPointer< Self > Pointer
itk::SmartPointer< Self > Pointer
Creation of an "otb" image which contains metadata.
Definition: otbImage.h:92
virtual void Write(const RegionType &streamRegion)=0
ImageBaseType::ConstPointer m_InputImage
virtual itk::ImageRegion< 2 > GetRegionToWrite() const =0
virtual ImageBaseType::ConstPointer GetInput() const
virtual bool CanStreamWrite() const =0
virtual ImageBaseType::Pointer GetInput()
SinkBase(ImageBaseType::ConstPointer inputImage)
otb::ImageFileWriter< TImage >::Pointer m_Writer
Streams a pipeline with multiple outputs. It writes each output to a file. Inputs are connected to th...
virtual void GenerateData(void) override
virtual void GenerateInputRequestedRegion() override
void SetTileDimensionTiledStreaming(unsigned int tileDimension)
itk::SmartPointer< const Self > ConstPointer
virtual void UpdateOutputInformation() override
virtual void InitializeStreaming()
ImageBaseType::IndexValueType IndexValueType
void operator=(const MultiImageFileWriter &)=delete
ImageBaseType::RegionType RegionType
void ResetAllRequestedRegions(ImageBaseType *imagePtr)
void SetNumberOfDivisionsStrippedStreaming(unsigned int nbDivisions)
void SetAutomaticAdaptativeStreaming(unsigned int availableRAM=0, double bias=1.0)
virtual void Update() override
StreamingManager< FakeOutputType > StreamingManagerType
virtual RegionType GetStreamRegion(int inputIndex)
void SetAutomaticStrippedStreaming(unsigned int availableRAM=0, double bias=1.0)
void SetNumberOfLinesStrippedStreaming(unsigned int nbLinesPerStrip)
ImageBaseType::SizeType SizeType
std::vector< RegionType > m_StreamRegionList
ImageBaseType::IndexType IndexType
void ObserveSourceFilterProgress(itk::Object *object, const itk::EventObject &event)
void SetNumberOfDivisionsTiledStreaming(unsigned int nbDivisions)
otb::Image< unsigned char, 2 > FakeOutputType
void SetAutomaticTiledStreaming(unsigned int availableRAM=0, double bias=1.0)
StreamingManagerType::Pointer m_StreamingManager
virtual void UpdateOutputData(itk::DataObject *) override
itk::SmartPointer< Self > Pointer
std::vector< std::shared_ptr< SinkBase > > SinkListType
ImageBaseType::SizeValueType SizeValueType
This class handles the streaming process used in the writers implementation.
itk::SmartPointer< Self > Pointer
OTBApplicationEngine_EXPORT void Write(const std::string &filename, Application::Pointer application)
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.