OTB  10.0.0
Orfeo Toolbox
otbImageFileWriter.hxx
Go to the documentation of this file.
1 /*
2  * Copyright (C) 2005-2024 Centre National d'Etudes Spatiales (CNES)
3  * Copyright (C) 2018-2020 CS Systemes d'Information (CS SI)
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 otbImageFileWriter_hxx
23 #define otbImageFileWriter_hxx
24 
25 #include "otbImageFileWriter.h"
26 #include "itkImageFileWriter.h"
27 
28 #include "itkObjectFactoryBase.h"
29 
30 #include "itkImageRegionMultidimensionalSplitter.h"
31 #include "otbImageIOFactory.h"
32 
33 #include "itkImageRegionIterator.h"
34 
35 #include "itkMetaDataObject.h"
36 #include "otbMetaDataKey.h"
37 #include "otbImageCommons.h"
38 
39 #include "otbConfigure.h"
40 
48 
50 
51 #include "otbStringUtils.h"
52 #include "otbUtils.h"
53 #include "itkProgressTransformer.h"
54 
55 namespace otb
56 {
57 
61 template <class TInputImage>
63  : m_NumberOfDivisions(0),
64  m_CurrentDivision(0),
65  m_DivisionProgress(0.0),
66  m_UserSpecifiedImageIO(true),
67  m_UserSpecifiedIORegion(false),
68  m_FactorySpecifiedImageIO(false),
69  m_UseCompression(false),
70  m_UseInputMetaDataDictionary(false),
71  m_WriteGeomFile(false),
72  m_FilenameHelper(),
73  m_IsObserving(true),
74  m_ObserverID(0),
75  m_IOComponents(0)
76 {
77  // Init output index shift
78  m_ShiftOutputIndex.Fill(0);
79 
80  // By default, we use tiled streaming, with automatic tile size
81  // We don't set any parameter, so the memory size is retrieved from the OTB configuration options
83 
85 }
86 
90 template <class TInputImage>
92 {
93 }
94 
95 template <class TInputImage>
97 {
98  typedef NumberOfDivisionsStrippedStreamingManager<TInputImage> NumberOfDivisionsStrippedStreamingManagerType;
99  typename NumberOfDivisionsStrippedStreamingManagerType::Pointer streamingManager = NumberOfDivisionsStrippedStreamingManagerType::New();
100  streamingManager->SetNumberOfDivisions(nbDivisions);
101  m_StreamingManager = streamingManager;
102 }
103 
104 template <class TInputImage>
106 {
107  typedef NumberOfDivisionsTiledStreamingManager<TInputImage> NumberOfDivisionsTiledStreamingManagerType;
108  typename NumberOfDivisionsTiledStreamingManagerType::Pointer streamingManager = NumberOfDivisionsTiledStreamingManagerType::New();
109  streamingManager->SetNumberOfDivisions(nbDivisions);
110  m_StreamingManager = streamingManager;
111 }
112 
113 template <class TInputImage>
115 {
116  typedef NumberOfLinesStrippedStreamingManager<TInputImage> NumberOfLinesStrippedStreamingManagerType;
117  typename NumberOfLinesStrippedStreamingManagerType::Pointer streamingManager = NumberOfLinesStrippedStreamingManagerType::New();
118  streamingManager->SetNumberOfLinesPerStrip(nbLinesPerStrip);
119  m_StreamingManager = streamingManager;
120 }
121 
122 template <class TInputImage>
123 void ImageFileWriter<TInputImage>::SetAutomaticStrippedStreaming(unsigned int availableRAM, double bias)
124 {
125  typedef RAMDrivenStrippedStreamingManager<TInputImage> RAMDrivenStrippedStreamingManagerType;
126  typename RAMDrivenStrippedStreamingManagerType::Pointer streamingManager = RAMDrivenStrippedStreamingManagerType::New();
127  streamingManager->SetAvailableRAMInMB(availableRAM);
128  streamingManager->SetBias(bias);
129  m_StreamingManager = streamingManager;
130 }
131 
132 template <class TInputImage>
134 {
135  typedef TileDimensionTiledStreamingManager<TInputImage> TileDimensionTiledStreamingManagerType;
136  typename TileDimensionTiledStreamingManagerType::Pointer streamingManager = TileDimensionTiledStreamingManagerType::New();
137  streamingManager->SetTileDimension(tileDimension);
138  m_StreamingManager = streamingManager;
139 }
140 
141 template <class TInputImage>
142 void ImageFileWriter<TInputImage>::SetAutomaticTiledStreaming(unsigned int availableRAM, double bias)
143 {
144  typedef RAMDrivenTiledStreamingManager<TInputImage> RAMDrivenTiledStreamingManagerType;
145  typename RAMDrivenTiledStreamingManagerType::Pointer streamingManager = RAMDrivenTiledStreamingManagerType::New();
146  streamingManager->SetAvailableRAMInMB(availableRAM);
147  streamingManager->SetBias(bias);
148  m_StreamingManager = streamingManager;
149 }
150 
151 template <class TInputImage>
152 void ImageFileWriter<TInputImage>::SetAutomaticAdaptativeStreaming(unsigned int availableRAM, double bias)
153 {
154  typedef RAMDrivenAdaptativeStreamingManager<TInputImage> RAMDrivenAdaptativeStreamingManagerType;
155  typename RAMDrivenAdaptativeStreamingManagerType::Pointer streamingManager = RAMDrivenAdaptativeStreamingManagerType::New();
156  streamingManager->SetAvailableRAMInMB(availableRAM);
157  streamingManager->SetBias(bias);
158  m_StreamingManager = streamingManager;
159 }
160 
164 template <class TInputImage>
165 void ImageFileWriter<TInputImage>::PrintSelf(std::ostream& os, itk::Indent indent) const
166 {
167  Superclass::PrintSelf(os, indent);
168 
169  os << indent << "File Name: " << (m_FileName.data() ? m_FileName.data() : "(none)") << std::endl;
170 
171  os << indent << "Image IO: ";
172  if (m_ImageIO.IsNull())
173  {
174  os << "(none)\n";
175  }
176  else
177  {
178  os << m_ImageIO << "\n";
179  }
180 
181  os << indent << "IO Region: " << m_IORegion << "\n";
182 
183  if (m_UseCompression)
184  {
185  os << indent << "Compression: On\n";
186  }
187  else
188  {
189  os << indent << "Compression: Off\n";
190  }
191 
192  if (m_UseInputMetaDataDictionary)
193  {
194  os << indent << "UseInputMetaDataDictionary: On\n";
195  }
196  else
197  {
198  os << indent << "UseInputMetaDataDictionary: Off\n";
199  }
200 
201  if (m_FactorySpecifiedImageIO)
202  {
203  os << indent << "FactorySpecifiedmageIO: On\n";
204  }
205  else
206  {
207  os << indent << "FactorySpecifiedmageIO: Off\n";
208  }
209 }
210 
211 //---------------------------------------------------------
212 template <class TInputImage>
213 void ImageFileWriter<TInputImage>::SetIORegion(const itk::ImageIORegion& region)
214 {
215  if (m_IORegion != region)
216  {
217  m_IORegion = region;
218  this->Modified();
219  m_UserSpecifiedIORegion = true;
220  }
221 }
222 
223 template <class TInputImage>
225 {
226  this->ProcessObject::SetNthInput(0, const_cast<InputImageType*>(input));
227 }
228 
229 template <class TInputImage>
231 {
232  if (this->GetNumberOfInputs() < 1)
233  {
234  return nullptr;
235  }
236 
237  return static_cast<const InputImageType*>(this->ProcessObject::GetInput(0));
238 }
239 
241 template <class TInputImage>
243 {
244  // Update output information on input image
245  InputImagePointer inputPtr = const_cast<InputImageType*>(this->GetInput());
246 
247  // Make sure input is available
248  if (inputPtr.IsNull())
249  {
250  itkExceptionMacro(<< "No input to writer");
251  }
252 
253  otb::Logger::Instance()->LogSetupInformation();
254 
256  if (m_FilenameHelper->StreamingTypeIsSet())
257  {
258  otbLogMacro(
259  Warning,
260  << "Streaming configuration through extended filename is used. Any previous streaming configuration (ram value, streaming mode ...) will be ignored.");
261 
262  std::string type = m_FilenameHelper->GetStreamingType();
263 
264  std::string sizemode = "auto";
265 
266  if (m_FilenameHelper->StreamingSizeModeIsSet())
267  {
268  sizemode = m_FilenameHelper->GetStreamingSizeMode();
269  }
270 
271  unsigned int sizevalue = 0;
272  // Save the DefaultRAM value for later
273  unsigned int oldDefaultRAM = m_StreamingManager->GetDefaultRAM();
274  if (sizemode == "auto")
275  {
276  sizevalue = oldDefaultRAM;
277  }
278 
279  if (m_FilenameHelper->StreamingSizeValueIsSet())
280  {
281  sizevalue = static_cast<unsigned int>(m_FilenameHelper->GetStreamingSizeValue());
282  }
283 
284  if (type == "auto")
285  {
286  if (sizemode != "auto")
287  {
288  otbLogMacro(Warning, << "In auto streaming type, the sizemode option will be ignored.");
289  }
290  if (sizevalue == 0)
291  {
292  otbLogMacro(Warning, << "sizemode is auto but sizevalue is 0. Value will be fetched from the OTB_MAX_RAM_HINT environment variable if set, or else use "
293  "the default value");
294  }
295  this->SetAutomaticAdaptativeStreaming(sizevalue);
296  }
297  else if (type == "tiled")
298  {
299  if (sizemode == "auto")
300  {
301  if (sizevalue == 0)
302  {
303  otbLogMacro(Warning, << "sizemode is auto but sizevalue is 0. Value will be fetched from the OTB_MAX_RAM_HINT environment variable if set, or else "
304  "use the default value");
305  }
306  this->SetAutomaticTiledStreaming(sizevalue);
307  }
308  else if (sizemode == "nbsplits")
309  {
310  if (sizevalue == 0)
311  {
312  otbLogMacro(Warning, << "Streaming sizemode is set to nbsplits but sizevalue is 0. This will result in undefined behaviour. Please consider setting "
313  "the sizevalue by using &streaming:sizevalue=x.");
314  }
315  this->SetNumberOfDivisionsTiledStreaming(sizevalue);
316  }
317  else if (sizemode == "height")
318  {
319  if (sizevalue == 0)
320  {
321  otbLogMacro(Warning, << "Streaming sizemode is set to height but sizevalue is 0. This will result in undefined behaviour. Please consider setting "
322  "the sizevalue by using &streaming:sizevalue=x.");
323  }
324 
325  this->SetTileDimensionTiledStreaming(sizevalue);
326  }
327  }
328  else if (type == "stripped")
329  {
330  if (sizemode == "auto")
331  {
332  if (sizevalue == 0)
333  {
334  otbLogMacro(
335  Warning, << "sizemode is auto but sizevalue is 0. Value will be fetched from configuration file if any, or from cmake configuration otherwise.");
336  }
337 
338  this->SetAutomaticStrippedStreaming(sizevalue);
339  }
340  else if (sizemode == "nbsplits")
341  {
342  if (sizevalue == 0)
343  {
344  otbLogMacro(Warning, << "Streaming sizemode is set to nbsplits but sizevalue is 0. This will result in undefined behaviour. Please consider setting "
345  "the sizevalue by using &streaming:sizevalue=x.");
346  }
347  this->SetNumberOfDivisionsStrippedStreaming(sizevalue);
348  }
349  else if (sizemode == "height")
350  {
351  if (sizevalue == 0)
352  {
353  otbLogMacro(Warning, << "Streaming sizemode is set to height but sizevalue is 0. This will result in undefined behaviour. Please consider setting "
354  "the sizevalue by using &streaming:sizevalue=x.");
355  }
356  this->SetNumberOfLinesStrippedStreaming(sizevalue);
357  }
358  }
359  else if (type == "none")
360  {
361  if (sizemode != "" || sizevalue != 0)
362  {
363  otbLogMacro(Warning, << "Streaming is explicitly disabled, sizemode and sizevalue will be ignored.");
364  }
365  this->SetNumberOfDivisionsTiledStreaming(0);
366  }
367 
368  // since we change the m_StreamingManager under the hood, we copy the DefaultRAM
369  // value to the new streamingManager.
370  m_StreamingManager->SetDefaultRAM(oldDefaultRAM);
371  }
372  else
373  {
374  if (m_FilenameHelper->StreamingSizeValueIsSet() || m_FilenameHelper->StreamingSizeModeIsSet())
375  {
376  otbLogMacro(Warning, << "No streaming type is set, streaming sizemode and sizevalue will be ignored.");
377  }
378  }
379 
382  if (m_FileName == "")
383  {
384  // Make sure that we can write the file given the name
385  itkExceptionMacro(<< "No filename was specified");
386  }
387 
388 
389  // Make sure CanWriteFile is called
390  // either from ImageIOFactory::CreateImageIO or directly in this file
391  // GDALImageIO uses it to store the filename
392  // and later answer to CanStreamWrite()
393  // This is a needed workaround to a defect in the itk::ImageIO interface
394 
395  if (m_ImageIO.IsNull()) // try creating via factory
396  {
397  this->SetImageIO(ImageIOFactory::CreateImageIO(m_FileName.c_str(), otb::ImageIOFactory::WriteMode));
398 
399  m_FactorySpecifiedImageIO = true;
400  }
401  else
402  {
403  if (!m_ImageIO->CanWriteFile(m_FileName.c_str()))
404  {
405  if (m_FactorySpecifiedImageIO)
406  {
407  m_ImageIO = ImageIOFactory::CreateImageIO(m_FileName.c_str(), otb::ImageIOFactory::WriteMode);
408  m_FactorySpecifiedImageIO = true;
409  }
410  }
411  }
412 
413  if (m_ImageIO.IsNull())
414  {
415  itk::ImageFileWriterException e(__FILE__, __LINE__);
416  std::ostringstream msg;
417  msg << "Cannot write image " << m_FileName << ". Probably unsupported format or incorrect filename extension.";
418  e.SetDescription(msg.str());
419  e.SetLocation(ITK_LOCATION);
420  throw e;
421  }
422 
423  // Manage extended filename
424  if ((strcmp(m_ImageIO->GetNameOfClass(), "GDALImageIO") == 0) &&
425  (m_FilenameHelper->gdalCreationOptionsIsSet() || m_FilenameHelper->WriteRPCTagsIsSet() || m_FilenameHelper->NoDataValueIsSet() || m_FilenameHelper->SrsValueIsSet()))
426  {
427  typename GDALImageIO::Pointer imageIO = dynamic_cast<GDALImageIO*>(m_ImageIO.GetPointer());
428 
429  if (imageIO.IsNull())
430  {
431  itk::ImageFileWriterException e(__FILE__, __LINE__);
432  std::ostringstream msg;
433  msg << " ImageIO is of kind GDALImageIO, but fails to dynamic_cast (this should never happen)." << std::endl;
434  e.SetDescription(msg.str());
435  throw e;
436  }
437 
438  imageIO->SetOptions(m_FilenameHelper->GetgdalCreationOptions());
439  imageIO->SetWriteRPCTags(m_FilenameHelper->GetWriteRPCTags());
440  if (m_FilenameHelper->NoDataValueIsSet())
441  imageIO->SetNoDataList(m_FilenameHelper->GetNoDataList());
442  if (m_FilenameHelper->SrsValueIsSet())
443  imageIO->SetEpsgCode(m_FilenameHelper->GetSrsValue());
444  }
445 
446 
452  InputImageRegionType inputRegion = inputPtr->GetLargestPossibleRegion();
453 
455  if (m_FilenameHelper->BoxIsSet())
456  {
457  std::vector<unsigned int> boxVector;
458  Utils::ConvertStringToVector(m_FilenameHelper->GetBox(), boxVector, "ExtendedFileName:box", ":");
460 
461  if (boxVector.size() != 4)
462  {
463  itk::ImageFileWriterException e(__FILE__, __LINE__);
464  std::ostringstream msg;
465  msg << "Invalid box option " << m_FilenameHelper->GetBox() << ". The box should contains four elements: startx:starty:sizex:sizey";
466  e.SetDescription(msg.str());
467  e.SetLocation(ITK_LOCATION);
468  throw e;
469  }
470 
471  typename InputImageRegionType::IndexType start;
472  typename InputImageRegionType::SizeType size;
473  start[0] = boxVector[0]; // first index on X
474  start[1] = boxVector[1]; // first index on Y
475  size[0] = boxVector[2]; // size along X
476  size[1] = boxVector[3]; // size along Y
477 
478  inputRegion.SetSize(size);
479  inputRegion.SetIndex(start);
480 
481  if (!inputRegion.Crop(inputPtr->GetLargestPossibleRegion()))
482  {
483  // Couldn't crop the region (requested region is outside the largest
484  // possible region). Throw an exception.
485 
486  // build an exception
487  itk::InvalidRequestedRegionError e(__FILE__, __LINE__);
488  e.SetLocation(ITK_LOCATION);
489  e.SetDescription("Requested box region is (at least partially) outside the largest possible region.");
490  e.SetDataObject(inputPtr);
491  throw e;
492  }
493  otbLogMacro(Info, << "Writing user defined region [" << start[0] << ", " << start[0] + size[0] - 1 << "]x[" << start[1] << ", " << start[1] + size[1]
494  << "]");
495  }
496  m_ShiftOutputIndex = inputRegion.GetIndex();
497 
505  if (m_ImageIO->CanStreamWrite() == false)
506  {
507  otbLogMacro(Warning, << "The file format of " << m_FileName << " does not support streaming. All data will be loaded to memory");
508  this->SetNumberOfDivisionsStrippedStreaming(1);
509  }
511 
515  else if (inputPtr->GetBufferedRegion() == inputRegion)
516  {
517  otbLogMacro(Debug, << "Buffered region is the largest possible region, there is no need for streaming.");
518  this->SetNumberOfDivisionsStrippedStreaming(1);
519  }
520  m_StreamingManager->PrepareStreaming(inputPtr, inputRegion);
521  m_NumberOfDivisions = m_StreamingManager->GetNumberOfSplits();
523 
524  const auto firstSplitSize = m_StreamingManager->GetSplit(0).GetSize();
525  otbLogMacro(Info, << "File " << m_FileName << " will be written in " << m_NumberOfDivisions << " blocks of " << firstSplitSize[0] << "x" << firstSplitSize[1]
526  << " pixels");
527 
528  //
529  // Setup the ImageIO with information from inputPtr
530  //
531  typename TInputImage::PointType origin;
532  inputPtr->TransformIndexToPhysicalPoint(inputRegion.GetIndex(), origin);
533  const typename TInputImage::SpacingType& spacing = inputPtr->GetSpacing();
534  const typename TInputImage::DirectionType& direction = inputPtr->GetDirection();
535  m_ImageIO->SetNumberOfDimensions(TInputImage::ImageDimension);
536  int direction_sign(0);
537  for (unsigned int i = 0; i < TInputImage::ImageDimension; ++i)
538  {
539  if (direction[i][i] < 0)
540  direction_sign = -1;
541  else
542  direction_sign = 1;
543  // Final image size
544  m_ImageIO->SetDimensions(i, inputRegion.GetSize(i));
545  m_ImageIO->SetSpacing(i, direction_sign * spacing[i]);
546  m_ImageIO->SetOrigin(i, origin[i]);
547 
548  vnl_vector<double> axisDirection(TInputImage::ImageDimension);
549  // Please note: direction cosines are stored as columns of the
550  // direction matrix
551  for (unsigned int j = 0; j < TInputImage::ImageDimension; ++j)
552  {
553  axisDirection[j] = direction_sign * direction[j][i];
554  }
555  m_ImageIO->SetDirection(i, axisDirection);
556  }
557 
558  m_ImageIO->SetUseCompression(m_UseCompression);
559  m_ImageIO->SetMetaDataDictionary(inputPtr->GetMetaDataDictionary());
560 
561  const ImageCommons* img_common = dynamic_cast<const ImageCommons*>(inputPtr.GetPointer());
562  if (img_common != nullptr)
563  {
564  m_ImageIO->SetImageMetadata(img_common->GetImageMetadata());
565  }
566 
568  // Setup the image IO for writing.
569  //
570  m_ImageIO->SetFileName(m_FileName);
571 
572  m_ImageIO->WriteImageInformation();
573 }
574 
578 template <class TInputImage>
580 {
581  this->UpdateOutputInformation();
582 
583  this->SetAbortGenerateData(0);
584  this->SetProgress(0.0);
585 
586  itk::ProgressTransformer pt( 0.0f, 1.0f, this );
590  this->InvokeEvent(itk::StartEvent());
591 
592  pt.GetProcessObject();
593  m_CurrentDivision = 0;
594  m_DivisionProgress = 0;
595 
596  // Get the source process object
597  InputImagePointer inputPtr = const_cast<InputImageType*>(this->GetInput());
598  itk::ProcessObject* source = inputPtr->GetSource();
599  m_IsObserving = false;
600  m_ObserverID = 0;
601 
602  // Check if source exists
603  if (source)
604  {
605  typedef itk::MemberCommand<Self> CommandType;
606  typedef typename CommandType::Pointer CommandPointerType;
607 
608  CommandPointerType command = CommandType::New();
609  command->SetCallbackFunction(this, &Self::ObserveSourceFilterProgress);
610 
611  m_ObserverID = source->AddObserver(itk::ProgressEvent(), command);
612  m_IsObserving = true;
613  }
614  else
615  {
616  otbLogMacro(Warning, << "Could not get the source process object. Progress report might be buggy");
617  }
618 
623  InputImageRegionType streamRegion;
624 
625  for (m_CurrentDivision = 0; m_CurrentDivision < m_NumberOfDivisions && !this->GetAbortGenerateData();
626  m_CurrentDivision++, m_DivisionProgress = 0, pt.GetProcessObject())
627  {
628  streamRegion = m_StreamingManager->GetSplit(m_CurrentDivision);
629 
630  inputPtr->SetRequestedRegion(streamRegion);
631  inputPtr->PropagateRequestedRegion();
632  inputPtr->UpdateOutputData();
633 
634  // Write the whole image
635  itk::ImageIORegion ioRegion(TInputImage::ImageDimension);
636  for (unsigned int i = 0; i < TInputImage::ImageDimension; ++i)
637  {
638  ioRegion.SetSize(i, streamRegion.GetSize(i));
639  // Set the ioRegion index using the shifted index ( (0,0 without box parameter))
640  ioRegion.SetIndex(i, streamRegion.GetIndex(i) - m_ShiftOutputIndex[i]);
641  }
642  this->SetIORegion(ioRegion);
643  m_ImageIO->SetIORegion(m_IORegion);
644 
645  // Start writing stream region in the image file
646  this->GenerateData();
647  }
648 
653  if (!this->GetAbortGenerateData())
654  {
655  this->UpdateProgress(1.0);
656  }
657  else
658  {
659  itk::ProcessAborted e(__FILE__, __LINE__);
660  e.SetLocation(ITK_LOCATION);
661  e.SetDescription("Image writing has been aborted");
662  throw e;
663  }
665 
666  // Notify end event observers
667  this->InvokeEvent(itk::EndEvent());
668 
669  if (m_IsObserving)
670  {
671  m_IsObserving = false;
672  source->RemoveObserver(m_ObserverID);
673  }
674 
678  this->ReleaseInputs();
679 
680  // Reset global shift on input region (box parameter)
681  // It allows calling multiple update over the writer
682  m_ShiftOutputIndex = inputPtr->GetLargestPossibleRegion().GetIndex();
683 }
684 
685 
689 template <class TInputImage>
691 {
692  const InputImageType* input = this->GetInput();
693  InputImagePointer cacheImage;
695 
696  // Make sure that the image is the right type and no more than
697  // four components.
698  typedef typename InputImageType::PixelType ImagePixelType;
699 
700  if (strcmp(input->GetNameOfClass(), "VectorImage") == 0)
701  {
702  typedef typename InputImageType::InternalPixelType VectorImagePixelType;
703  m_ImageIO->SetPixelTypeInfo(typeid(VectorImagePixelType));
704 
705  typedef typename InputImageType::AccessorFunctorType AccessorFunctorType;
706  m_ImageIO->SetNumberOfComponents(AccessorFunctorType::GetVectorLength(input));
707 
708  m_IOComponents = m_ImageIO->GetNumberOfComponents();
709  m_BandList.clear();
710  if (m_FilenameHelper->BandRangeIsSet())
711  {
712  // get band range
713  bool retBandRange = m_FilenameHelper->ResolveBandRange(m_FilenameHelper->GetBandRange(), m_IOComponents, m_BandList);
714  if (retBandRange == false || m_BandList.empty())
715  {
716  // invalid range
717  itkGenericExceptionMacro("The given band range is either empty or invalid for a " << m_IOComponents << " bands input image!");
718  }
719  }
720  }
721  else
722  {
723  // Set the pixel and component type; the number of components.
724  m_ImageIO->SetPixelTypeInfo(typeid(ImagePixelType));
725  }
726 
727  // Setup the image IO for writing.
728  //
729  // okay, now extract the data as a raw buffer pointer
730  const void* dataPtr = (const void*)input->GetBufferPointer();
731 
732  // check that the image's buffered region is the same as
733  // ImageIO is expecting and we requested
734  InputImageRegionType ioRegion;
735 
736  // No shift of the ioRegion from the buffered region is expected
737  itk::ImageIORegionAdaptor<TInputImage::ImageDimension>::Convert(m_ImageIO->GetIORegion(), ioRegion, m_ShiftOutputIndex);
738  InputImageRegionType bufferedRegion = input->GetBufferedRegion();
739 
740  // before this test, bad stuff would happened when they don't match.
741  // In case of the buffer has not enough components, adapt the region.
742  if ((bufferedRegion != ioRegion) || (m_FilenameHelper->BandRangeIsSet() && (m_IOComponents < m_BandList.size())))
743  {
744  if (m_NumberOfDivisions > 1 || m_UserSpecifiedIORegion)
745  {
746  cacheImage = InputImageType::New();
747  cacheImage->CopyInformation(input);
748 
749  // set number of components at the band range size
750  if (m_FilenameHelper->BandRangeIsSet() && (m_IOComponents < m_BandList.size()))
751  {
752  cacheImage->SetNumberOfComponentsPerPixel(m_BandList.size());
753  }
754 
755  cacheImage->SetBufferedRegion(ioRegion);
756  cacheImage->Allocate();
757 
758  // set number of components at the initial size
759  if (m_FilenameHelper->BandRangeIsSet() && (m_IOComponents < m_BandList.size()))
760  {
761  cacheImage->SetNumberOfComponentsPerPixel(m_IOComponents);
762  }
763 
764  typedef itk::ImageRegionConstIterator<TInputImage> ConstIteratorType;
765  typedef itk::ImageRegionIterator<TInputImage> IteratorType;
766 
767  ConstIteratorType in(input, ioRegion);
768  IteratorType out(cacheImage, ioRegion);
769 
770  // copy the data into a buffer to match the ioregion
771  for (in.GoToBegin(), out.GoToBegin(); !in.IsAtEnd(); ++in, ++out)
772  {
773  out.Set(in.Get());
774  }
775 
776  dataPtr = (const void*)cacheImage->GetBufferPointer();
777  }
778  else
779  {
780  itk::ImageFileWriterException e(__FILE__, __LINE__);
781  std::ostringstream msg;
782  msg << "Did not get requested region!" << std::endl;
783  msg << "Requested:" << std::endl;
784  msg << ioRegion;
785  msg << "Actual:" << std::endl;
786  msg << bufferedRegion;
787  e.SetDescription(msg.str());
788  e.SetLocation(ITK_LOCATION);
789  throw e;
790  }
791  }
792 
793  if (m_FilenameHelper->BandRangeIsSet() && (!m_BandList.empty()))
794  {
795  // Adapt the image size with the region and take into account a potential
796  // remapping of the components. m_BandList is empty if no band range is set
797  m_ImageIO->DoMapBuffer(const_cast<void*>(dataPtr), bufferedRegion.GetNumberOfPixels(), this->m_BandList);
798  m_ImageIO->SetNumberOfComponents(m_BandList.size());
799  }
800 
801  m_ImageIO->Write(dataPtr);
802 }
803 
804 template <class TInputImage>
805 void ImageFileWriter<TInputImage>::SetFileName(const std::string& extendedFileName)
806 {
807  this->m_FilenameHelper->SetExtendedFileName(extendedFileName);
808  m_FileName = this->m_FilenameHelper->GetSimpleFileName();
809  m_ImageIO = nullptr;
810  this->Modified();
811 }
812 
813 template <class TInputImage>
815 {
816  return this->m_FilenameHelper->GetSimpleFileName();
817 }
818 
819 template <class TInputImage>
821 {
822  std::lock_guard<std::mutex> mutexHolder(m_Lock);
823  // protected read here
824  bool ret = Superclass::GetAbortGenerateData();
825  if (ret)
828 }
829 
830 template <class TInputImage>
832 {
833  std::lock_guard<std::mutex> mutexHolder(m_Lock);
834  Superclass::SetAbortGenerateData(val);
835 }
836 
837 } // end namespace otb
838 
839 #endif
ImageIO object for reading and writing images with GDAL.
itk::SmartPointer< Self > Pointer
const ImageMetadata & GetImageMetadata() const
void SetImageMetadata(ImageMetadata imd)
void SetIORegion(const itk::ImageIORegion &region)
void SetAbortGenerateData(const bool val) override
void PrintSelf(std::ostream &os, itk::Indent indent) const override
InputImageType::RegionType InputImageRegionType
void SetNumberOfLinesStrippedStreaming(unsigned int nbLinesPerStrip)
virtual void SetInput(const InputImageType *input)
InputImageType::Pointer InputImagePointer
void SetAutomaticTiledStreaming(unsigned int availableRAM=0, double bias=1.0)
InputIndexType m_ShiftOutputIndex
const bool & GetAbortGenerateData() const override
void SetAutomaticStrippedStreaming(unsigned int availableRAM=0, double bias=1.0)
void SetNumberOfDivisionsStrippedStreaming(unsigned int nbDivisions)
virtual void SetFileName(const std::string &extendedFileName)
FNameHelperType::Pointer m_FilenameHelper
void SetTileDimensionTiledStreaming(unsigned int tileDimension)
void SetNumberOfDivisionsTiledStreaming(unsigned int nbDivisions)
void GenerateOutputInformation(void) override
virtual const char * GetFileName() const
void SetAutomaticAdaptativeStreaming(unsigned int availableRAM=0, double bias=1.0)
const InputImageType * GetInput()
void GenerateData(void) override
static ImageIOBasePointer CreateImageIO(const char *path, FileModeType mode)
static Logger * Instance()
This class computes the divisions needed to stream an image by strips, driven by a user-defined numbe...
This class computes the divisions needed to stream an image by strips, driven by a user-defined numbe...
This class computes the divisions needed to stream an image by strips, driven by a user-defined desir...
This class computes the divisions needed to stream an image according to the input image tiling schem...
This class computes the divisions needed to stream an image by strips, according to a user-defined av...
This class computes the divisions needed to stream an image in square tiles, according to a user-defi...
This class computes the divisions needed to stream an image in square tiles, driven by a user-defined...
OTBCommon_EXPORT bool const FalseConstant
void ConvertStringToVector(std::string const &str, T &ret, std::string const &errmsg, const char *delims=" ")
OTBCommon_EXPORT bool const TrueConstant
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.
#define otbLogMacro(level, msg)
Definition: otbMacro.h:104