OTB  10.0.0
Orfeo Toolbox
otbTileImageFilter.hxx
Go to the documentation of this file.
1 /*
2  * Copyright (C) 2005-2024 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 otbTileImageFilter_hxx
22 #define otbTileImageFilter_hxx
23 
24 #include "otbTileImageFilter.h"
25 #include "itkImageRegionIterator.h"
26 
27 namespace otb
28 {
29 template <class TImage>
31 {
32  this->DynamicMultiThreadingOn();
33 }
34 
35 template <class TImage>
37 {
38 }
39 
40 template <class TImage>
41 void TileImageFilter<TImage>::PrintSelf(std::ostream& os, itk::Indent indent) const
42 {
43  Superclass::PrintSelf(os, indent);
44  os << indent << "Layout: " << m_Layout << std::endl;
45 }
46 
47 template <class TImage>
49 {
50  // First, call superclass implementation
51  Superclass::GenerateOutputInformation();
52 
53  // Get the number of images
54  unsigned int numberOfImages = m_Layout[0] * m_Layout[1];
55 
56  // Check we have enough inputs
57  if (numberOfImages != this->GetNumberOfInputs())
58  {
59  itkExceptionMacro(<< "Layout has " << numberOfImages << " tiles, but only " << this->GetNumberOfInputs() << " inputs are found.");
60  }
61 
62  typename ImageType::SpacingType spacing = this->GetInput()->GetSignedSpacing();
63  unsigned int nbComp = this->GetInput()->GetNumberOfComponentsPerPixel();
64 
65  m_ColumnsSizes.clear();
66  m_RowsSizes.clear();
67 
68  // Loop on the layout
69  for (unsigned int col = 0; col < m_Layout[0]; ++col)
70  {
71  for (unsigned int row = 0; row < m_Layout[1]; ++row)
72  {
73  // First, get current tile
74  const ImageType* currentTile = this->GetInput(col + row * m_Layout[0]);
75  typename ImageType::SizeType currentSize = currentTile->GetLargestPossibleRegion().GetSize();
76 
77  // Retrieve row and column sizes
78  if (col == 0)
79  {
80  m_RowsSizes.push_back(currentSize[1]);
81  }
82  if (row == 0)
83  {
84  m_ColumnsSizes.push_back(currentSize[0]);
85  }
86 
87  // Check for consistent layout
88  if (currentSize[1] != m_RowsSizes[row] || currentSize[0] != m_ColumnsSizes[col])
89  {
90  itkExceptionMacro(<< "Inconsistent sizes in layout detected!");
91  }
92 
93  if (spacing != currentTile->GetSignedSpacing())
94  {
95  itkExceptionMacro(<< "Inconsistent spacings in layout detected!");
96  }
97 
98  if (nbComp != currentTile->GetNumberOfComponentsPerPixel())
99  {
100  itkExceptionMacro(<< "Inconsistent number of components in layout detected!");
101  }
102  }
103  }
104 
105  // Compute total size
106  typename ImageType::SizeType totalSize;
107  totalSize.Fill(0);
108 
109  for (unsigned int i = 0; i < m_Layout[0]; ++i)
110  {
111  totalSize[0] += m_ColumnsSizes[i];
112  }
113 
114  for (unsigned int i = 0; i < m_Layout[1]; ++i)
115  {
116  totalSize[1] += m_RowsSizes[i];
117  }
118 
119  // Fill out output image settings
120  typename ImageType::RegionType outRegion;
121  outRegion.SetSize(totalSize);
122 
123  // std::cout<<"Output largest region: "<<outRegion<<std::endl;
124 
125  ImageType* outputPtr = this->GetOutput();
126 
127  // Copy information from first tile
128  outputPtr->CopyInformation(this->GetInput());
129 
130  // Set region
131  outputPtr->SetLargestPossibleRegion(outRegion);
132 }
133 
134 template <class TImage>
136 {
137  // Get output requested region
138  RegionType outRegion = this->GetOutput()->GetRequestedRegion();
139 
140  // Loop on al input images
141  unsigned int numberOfImages = m_Layout[0] * m_Layout[1];
142 
143  // std::cout<<"Number of input images: "<<this->GetNumberOfInputs()<<std::endl;
144 
145  // std::cout<<"Out region: "<<outRegion<<std::endl;
146 
147  for (unsigned int i = 0; i < numberOfImages; ++i)
148  {
149  ImageType* inputTile = const_cast<ImageType*>(this->GetInput(i));
150 
151  RegionType inRegion = OutputRegionToInputRegion(i, outRegion);
152 
153  // std::cout<<"In region: "<<i<<", "<<inRegion<<std::endl;
154  // std::cout<<"tile pointer: "<<inputTile<<std::endl;
155 
156  inputTile->SetRequestedRegion(inRegion);
157  }
158 }
159 
160 template <class TImage>
162 {
163  // Retrieve output image pointer
164  ImageType* outputPtr = this->GetOutput();
165 
166  // Loop on all input tiles
167  unsigned int numberOfImages = m_Layout[0] * m_Layout[1];
168 
169  for (unsigned int i = 0; i < numberOfImages; ++i)
170  {
171  const ImageType* inputTile = this->GetInput(i);
172 
173  RegionType inRegion = OutputRegionToInputRegion(i, outputRegionForThread);
174  RegionType outRegion = InputRegionToOutputRegion(i, inRegion);
175 
176  // std::cout<<"[thread] outRegion: "<<outRegion<<std::endl;
177 
178  if (inRegion.GetNumberOfPixels() > 0)
179  {
180  // std::cout<<"[thread] "<<i<<" inRegion: "<<inRegion<<std::endl;
181 
182  // std::cout<<"[thread] "<<i<< inputTile->GetBufferedRegion()<<std::endl;
183 
184  itk::ImageRegionConstIterator<ImageType> inIt(inputTile, inRegion);
185  itk::ImageRegionIterator<ImageType> outIt(outputPtr, outRegion);
186 
187  inIt.GoToBegin();
188  outIt.GoToBegin();
189 
190  while (!inIt.IsAtEnd() && !outIt.IsAtEnd())
191  {
192  outIt.Set(inIt.Get());
193  ++inIt;
194  ++outIt;
195  }
196  }
197  }
198 }
199 
200 template <class TImage>
202 {
203  // Retrieve the nth tile pointer
204  const ImageType* tile = this->GetInput(tileIndex);
205 
206  unsigned int tileYIndex = tileIndex / m_Layout[0];
207  unsigned int tileXIndex = tileIndex % m_Layout[0];
208 
209  RegionType out2inRegion = requestedRegion;
210  typename RegionType::IndexType regionIndex = out2inRegion.GetIndex();
211 
212  // Compute tile offsets
213  for (unsigned int i = 0; i < tileXIndex; ++i)
214  {
215  regionIndex[0] -= m_ColumnsSizes.at(i);
216  }
217 
218  for (unsigned int i = 0; i < tileYIndex; ++i)
219  {
220  regionIndex[1] -= m_RowsSizes.at(i);
221  }
222 
223  out2inRegion.SetIndex(regionIndex);
224 
225 
226  // Crop input region
227  if (!out2inRegion.Crop(tile->GetLargestPossibleRegion()))
228  {
229  typename RegionType::IndexType nullIndex;
230  nullIndex.Fill(0);
231  out2inRegion.SetIndex(nullIndex);
232 
233  SizeType nullSize;
234  nullSize.Fill(0);
235  out2inRegion.SetSize(nullSize);
236  }
237 
238  return out2inRegion;
239 }
240 
241 template <class TImage>
243 {
244  // Retrieve the nth tile pointer
245  // const ImageType * tile = this->GetInput(tileIndex);
246 
247  unsigned int tileYIndex = tileIndex / m_Layout[0];
248  unsigned int tileXIndex = tileIndex % m_Layout[0];
249 
250  RegionType out2inRegion = requestedRegion;
251  typename RegionType::IndexType regionIndex = out2inRegion.GetIndex();
252 
253 
254  // Compute tile offsets
255  for (unsigned int i = 0; i < tileXIndex; ++i)
256  {
257  regionIndex[0] += m_ColumnsSizes.at(i);
258  }
259 
260  for (unsigned int i = 0; i < tileYIndex; ++i)
261  {
262  regionIndex[1] += m_RowsSizes.at(i);
263  }
264 
265  out2inRegion.SetIndex(regionIndex);
266 
267  // Crop input region
268  if (!out2inRegion.Crop(this->GetOutput()->GetLargestPossibleRegion()))
269  {
270  SizeType nullSize;
271  nullSize.Fill(0);
272  out2inRegion.SetSize(nullSize);
273  }
274 
275  return out2inRegion;
276 }
277 
278 
279 } // end namespace otb
280 
281 #endif
void GenerateInputRequestedRegion() override
ImageType::SizeType SizeType
void GenerateOutputInformation() override
void DynamicThreadedGenerateData(const RegionType &outputRegionForThread) override
RegionType InputRegionToOutputRegion(unsigned int tileIndex, const RegionType &requestedRegion)
ImageType::RegionType RegionType
void PrintSelf(std::ostream &os, itk::Indent indent) const override
RegionType OutputRegionToInputRegion(unsigned int tileIndex, const RegionType &requestedRegion)
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.