OTB  10.0.0
Orfeo Toolbox
otbSubsampleImageFilter.hxx
Go to the documentation of this file.
1 /*
2  * Copyright (C) 2005-2024 Centre National d'Etudes Spatiales (CNES)
3  * Copyright (C) 2007-2012 Institut Mines Telecom / Telecom Bretagne
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 
23 #ifndef otbSubsampleImageFilter_hxx
24 #define otbSubsampleImageFilter_hxx
26 
27 #include "otbMacro.h"
29 #include "itkImageRegionIterator.h"
30 
31 namespace otb
32 {
33 
34 template <class TInputImage, class TOutputImage, Wavelet::WaveletDirection TDirectionOfTransformation>
36 {
37  Superclass::PrintSelf(os, indent);
38  os << indent << "SubsampleFactor = [" << m_SubsampleFactor[0];
39  for (unsigned int i = 1; i < InputImageDimension; ++i)
40  {
41  os << ", " << m_SubsampleFactor[i];
42  }
43  os << "]\n";
44 }
45 
46 template <class TInputImage, class TOutputImage, Wavelet::WaveletDirection TDirectionOfTransformation>
48 {
49  for (unsigned int i = 0; i < InputImageDimension; ++i)
50  {
51  if (m_SubsampleFactor[i] != 1)
52  return false;
53  }
54 
55  return true;
56 }
57 
58 template <class TInputImage, class TOutputImage, Wavelet::WaveletDirection TDirectionOfTransformation>
60 {
61  Superclass::GenerateOutputInformation();
62 
63  if (!IsSubsampleFactorOne())
64  {
65  OutputImageRegionType newRegion;
66  this->CallCopyInputRegionToOutputRegion(newRegion, this->GetInput()->GetLargestPossibleRegion());
67  this->GetOutput()->SetRegions(newRegion);
68  }
69 }
70 
71 template <class TInputImage, class TOutputImage, Wavelet::WaveletDirection TDirectionOfTransformation>
73  const OutputImageRegionType& srcRegion)
74 {
75  Superclass::CallCopyOutputRegionToInputRegion(destRegion, srcRegion);
76 
77  if (DirectionOfTransformation == Wavelet::INVERSE)
78  {
79  typename OutputImageRegionType::IndexType srcIndex = srcRegion.GetIndex();
80  typename OutputImageRegionType::SizeType srcSize = srcRegion.GetSize();
81 
82  typename InputImageRegionType::IndexType destIndex;
83  typename InputImageRegionType::SizeType destSize;
84 
85  for (unsigned int i = 0; i < InputImageDimension; ++i)
86  {
87  // TODO: This seems not right in odd index cases
88  destIndex[i] = srcIndex[i] / m_SubsampleFactor[i];
89  destSize[i] = srcSize[i] / m_SubsampleFactor[i];
90  }
91 
92  destRegion.SetIndex(destIndex);
93  destRegion.SetSize(destSize);
94  }
95 }
96 
97 template <class TInputImage, class TOutputImage, Wavelet::WaveletDirection TDirectionOfTransformation>
99  const InputImageRegionType& srcRegion)
100 {
101  Superclass::CallCopyInputRegionToOutputRegion(destRegion, srcRegion);
102 
103  if (DirectionOfTransformation == Wavelet::INVERSE)
104  {
105  typename InputImageRegionType::IndexType srcIndex = srcRegion.GetIndex();
106  typename InputImageRegionType::SizeType srcSize = srcRegion.GetSize();
107 
108  typename OutputImageRegionType::IndexType destIndex;
109  typename OutputImageRegionType::SizeType destSize;
110 
111  for (unsigned int i = 0; i < InputImageDimension; ++i)
112  {
113  destIndex[i] = srcIndex[i] * m_SubsampleFactor[i];
114  destSize[i] = srcSize[i] * m_SubsampleFactor[i];
115  }
116 
117  destRegion.SetIndex(destIndex);
118  destRegion.SetSize(destSize);
119  }
120 }
121 
122 template <class TInputImage, class TOutputImage, Wavelet::WaveletDirection TDirectionOfTransformation>
124 {
125  OutputImagePointerType output = this->GetOutput();
126  output->FillBuffer(0);
127 }
128 
129 template <class TInputImage, class TOutputImage, Wavelet::WaveletDirection TDirectionOfTransformation>
131 {
132  OutputImagePointerType output = this->GetOutput();
133 
134  SubsampledImageRegionIterator<OutputImageType> outputIter(this->GetOutput(), outputRegionForThread);
135  outputIter.SetSubsampleFactor(1);
136  outputIter.GoToBegin();
137 
138  InputImageRegionType inputRegionForThread;
139  this->CallCopyOutputRegionToInputRegion(inputRegionForThread, outputRegionForThread);
140 
141  SubsampledImageRegionConstIterator<InputImageType> inputIter(this->GetInput(), inputRegionForThread);
142 
143  if (DirectionOfTransformation == Wavelet::FORWARD)
144  {
145  inputIter.SetSubsampleFactor(GetSubsampleFactor());
146  inputIter.GoToBegin();
147 
148  while (!inputIter.IsAtEnd())
149  {
150  outputIter.SetOffset(static_cast<typename SubsampledImageRegionIterator<OutputImageType>::OffsetType>(inputIter.GetOffset()));
151  outputIter.Set(static_cast<OutputPixelType>(inputIter.Get()));
152  ++inputIter;
153  }
154  }
155  else
156  {
157  inputIter.SetSubsampleFactor(1);
158  inputIter.GoToBegin();
159 
160  while (!inputIter.IsAtEnd())
161  {
162  InputImageIndexType inputIndex = inputIter.GetIndex();
163  OutputImageIndexType outputIndex;
164  for (unsigned int i = 0; i < OutputImageDimension; ++i)
165  {
166  outputIndex[i] = inputIndex[i] * m_SubsampleFactor[i];
167  }
168  outputIter.SetIndex(outputIndex);
169  outputIter.Set(static_cast<OutputPixelType>(inputIter.Get()));
170  ++inputIter;
171  }
172  }
173 }
174 
175 } // end of namespace otb
176 
177 #endif
OutputImageType::IndexType OutputImageIndexType
OutputImageType::RegionType OutputImageRegionType
void DynamicThreadedGenerateData(const OutputImageRegionType &outputRegionForThread) override
InputImageType::RegionType InputImageRegionType
InputImageType::IndexType InputImageIndexType
void PrintSelf(std::ostream &os, itk::Indent indent) const override
void CallCopyOutputRegionToInputRegion(InputImageRegionType &destRegion, const OutputImageRegionType &srcRegion) override
void CallCopyInputRegionToOutputRegion(OutputImageRegionType &destRegion, const InputImageRegionType &srcRegion) override
OutputImageType::PixelType OutputPixelType
OutputImageType::Pointer OutputImagePointerType
Regular subsample iterator over an image.
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.