OTB  9.0.0
Orfeo Toolbox
otbScalarBufferToImageFileWriter.hxx
Go to the documentation of this file.
1 /*
2  * Copyright (C) 2005-2022 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 otbScalarBufferToImageFileWriter_hxx
22 #define otbScalarBufferToImageFileWriter_hxx
23 
24 #include "otbMacro.h"
25 
27 #include "itkImageRegionIterator.h"
28 #include "itkImageRegionIteratorWithIndex.h"
29 
30 namespace otb
31 {
32 template <class TBufferType, class TOutputPixelType>
34  : m_Buffer(nullptr), m_NumberOfChannels(0), m_InverseXSpacing(false)
35 {
37  m_ImageSize.Fill(0);
38 }
39 
40 
41 template <class TBufferType, class TOutputPixelType>
43 {
44  // Check image parameters
45  if ((m_ImageSize[0] == 0) || (m_ImageSize[1] == 0))
46  {
47  itkExceptionMacro("Invalid output image size, Size can't be null.");
48  }
49 
50  if (m_NumberOfChannels == 0)
51  {
52  itkExceptionMacro("Invalid output image number of channels.");
53  }
54 
55  RegionType lRegion;
56  IndexType lId;
57  lId.Fill(0);
58  lRegion.SetIndex(lId);
59  lRegion.SetSize(m_ImageSize);
60 
61  typename ImageType::Pointer lImage = ImageType::New();
62  lImage->SetRegions(lRegion);
63  lImage->SetNumberOfComponentsPerPixel(m_NumberOfChannels);
64  lImage->Allocate();
65  PixelType lPix;
66 
67  lPix.SetSize(m_NumberOfChannels);
68  lPix.Fill(itk::NumericTraits<OutputPixelType>::Zero);
69  lImage->FillBuffer(lPix);
70 
71  // 1 specific loop for each case to save time processing
72  if (m_InverseXSpacing == false)
73  {
74  itk::ImageRegionIterator<ImageType> it(lImage, lRegion);
75  it.GoToBegin();
76 
77  unsigned int cpt(0);
78  while (it.IsAtEnd() == false)
79  {
80  for (unsigned int i = 0; i < m_NumberOfChannels; ++i)
81  {
82  lPix[i] = static_cast<OutputPixelType>(m_Buffer[cpt]);
83  ++cpt;
84  }
85 
86  it.Set(lPix);
87  ++it;
88  }
89  }
90  else
91  {
92  itk::ImageRegionIteratorWithIndex<ImageType> it(lImage, lRegion);
93  it.GoToBegin();
94  // cpt is the first component of the last pixel
95  unsigned int cpt(0);
96  while (it.IsAtEnd() == false)
97  {
98  IndexType index = it.GetIndex();
99  cpt = (m_ImageSize[1] - 1 - index[1]) * m_NumberOfChannels * m_ImageSize[0] + m_NumberOfChannels * index[0];
100 
101  for (unsigned int i = 0; i < m_NumberOfChannels; ++i)
102  {
103  lPix[i] = static_cast<OutputPixelType>(m_Buffer[cpt + i]);
104  }
105 
106  it.Set(lPix);
107  ++it;
108  }
109  }
110 
111  m_Writer->SetInput(lImage);
112  m_Writer->Update();
113 }
114 
115 template <class TBufferType, class TOutputPixelType>
116 void ScalarBufferToImageFileWriter<TBufferType, TOutputPixelType>::PrintSelf(std::ostream& os, itk::Indent indent) const
117 {
118  Superclass::PrintSelf(os, indent);
119  os << indent << "FileName" << m_Writer->GetFileName() << std::endl;
120  os << indent << "Size" << m_ImageSize << std::endl;
121  os << indent << "NumberOfChannels" << m_NumberOfChannels << std::endl;
122 }
123 
124 } // end namespace otb
125 
126 #endif
otb::ScalarBufferToImageFileWriter::m_ImageSize
SizeType m_ImageSize
Definition: otbScalarBufferToImageFileWriter.h:127
otb::ScalarBufferToImageFileWriter::PixelType
ImageType::PixelType PixelType
Definition: otbScalarBufferToImageFileWriter.h:63
otbScalarBufferToImageFileWriter.h
otb::VectorImage::Pointer
itk::SmartPointer< Self > Pointer
Definition: otbVectorImage.h:53
otb::ScalarBufferToImageFileWriter::PrintSelf
void PrintSelf(std::ostream &os, itk::Indent indent) const override
Definition: otbScalarBufferToImageFileWriter.hxx:116
otb::ScalarBufferToImageFileWriter::OutputPixelType
TOutputPixelType OutputPixelType
Definition: otbScalarBufferToImageFileWriter.h:59
otb
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.
Definition: otbJoinContainer.h:32
otbMacro.h
otb::ScalarBufferToImageFileWriter::GenerateData
void GenerateData() override
Definition: otbScalarBufferToImageFileWriter.hxx:42
otb::ScalarBufferToImageFileWriter::ScalarBufferToImageFileWriter
ScalarBufferToImageFileWriter()
Definition: otbScalarBufferToImageFileWriter.hxx:33
otb::ScalarBufferToImageFileWriter::m_Writer
WriterPointer m_Writer
Definition: otbScalarBufferToImageFileWriter.h:118
otb::ImageFileWriter::New
static Pointer New()
otb::ScalarBufferToImageFileWriter::RegionType
ImageType::RegionType RegionType
Definition: otbScalarBufferToImageFileWriter.h:64
otb::ScalarBufferToImageFileWriter::IndexType
ImageType::IndexType IndexType
Definition: otbScalarBufferToImageFileWriter.h:66