OTB  10.0.0
Orfeo Toolbox
otbStreamingSimpleMosaicFilter.hxx
Go to the documentation of this file.
1 /*
2  * Copyright (C) 1999-2011 Insight Software Consortium
3  * Copyright (C) 2005-2024 Centre National d'Etudes Spatiales (CNES)
4  * Copyright (C) 2016-2019 IRSTEA
5  *
6  * This file is part of Orfeo Toolbox
7  *
8  * https://www.orfeo-toolbox.org/
9  *
10  * Licensed under the Apache License, Version 2.0 (the "License");
11  * you may not use this file except in compliance with the License.
12  * You may obtain a copy of the License at
13  *
14  * http://www.apache.org/licenses/LICENSE-2.0
15  *
16  * Unless required by applicable law or agreed to in writing, software
17  * distributed under the License is distributed on an "AS IS" BASIS,
18  * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
19  * See the License for the specific language governing permissions and
20  * limitations under the License.
21  */
22 #ifndef __StreamingSimpleMosaicFilter_hxx
23 #define __StreamingSimpleMosaicFilter_hxx
24 
26 
27 namespace otb
28 {
29 
33 template <class TInputImage, class TOutputImage, class TInternalValueType>
35 {
36 
37  // Get output pointer
38  OutputImageType* mosaicImage = this->GetOutput();
39 
40  // Get number of used inputs
41  const unsigned int nbOfUsedInputImages = Superclass::GetNumberOfUsedInputImages();
42 
43  // Get number of bands
44  const unsigned int nBands = Superclass::GetNumberOfBands();
45 
46  // Iterate through the thread region
47  IteratorType outputIt(mosaicImage, outputRegionForThread);
48 
49  // Prepare interpolated pixel
50  InternalPixelType interpolatedMathPixel;
51  interpolatedMathPixel.SetSize(nBands);
52 
53  // Prepare input pointers, interpolators, and valid regions (input images)
54  typename std::vector<InputImageType*> currentImage;
55  typename std::vector<InterpolatorPointerType> interp;
56  Superclass::PrepareImageAccessors(currentImage, interp);
57 
58  // Container for geo coordinates
59  OutputImagePointType geoPoint;
60 
61  for (outputIt.GoToBegin(); !outputIt.IsAtEnd(); ++outputIt)
62  {
63  // Prepare output pixel
64  OutputImagePixelType outputPixel(Superclass::GetNoDataOutputPixel());
65 
66  // Current pixel --> Geographical point
67  mosaicImage->TransformIndexToPhysicalPoint(outputIt.GetIndex(), geoPoint);
68 
69  // Loop on used input images
70  for (unsigned int i = 0; i < nbOfUsedInputImages; i++)
71  {
72  // Get the input image pointer
73  unsigned int imgIndex = Superclass::GetUsedInputImageIndice(i);
74 
75  // Check if the point is inside the transformed thread region
76  // (i.e. the region in the current input image which match the thread
77  // region)
78  if (interp[i]->IsInsideBuffer(geoPoint))
79  {
80 
81  // Compute the interpolated pixel value
82  InputImagePixelType interpolatedPixel = interp[i]->Evaluate(geoPoint);
83 
84  // Check that interpolated pixel is not empty
85  if (Superclass::IsPixelNotEmpty(interpolatedPixel))
86  {
87  // Update the output pixel
88  for (unsigned int band = 0; band < nBands; band++)
89  {
90  if (this->GetShiftScaleInputImages())
91  {
92  InternalValueType pixelValue = static_cast<InternalValueType>(interpolatedPixel[band]);
93  this->ShiftScaleValue(pixelValue, imgIndex, band);
94  outputPixel[band] = static_cast<OutputImageInternalPixelType>(pixelValue);
95  }
96  else
97  {
98  outputPixel[band] = static_cast<OutputImageInternalPixelType>(interpolatedPixel[band]);
99  }
100  }
101  } // Interpolated pixel is not empty
102  } // point inside buffer
103  } // next image
104 
105  // Update output pixel value
106  outputIt.Set(outputPixel);
107 
108  } // next output pixel
109 }
110 
111 } // end namespace gtb
112 
113 #endif
Superclass::OutputImagePointType OutputImagePointType
Superclass::InputImagePixelType InputImagePixelType
Superclass::OutputImagePixelType OutputImagePixelType
Superclass::OutputImageRegionType OutputImageRegionType
Superclass::OutputImageInternalPixelType OutputImageInternalPixelType
void DynamicThreadedGenerateData(const OutputImageRegionType &outputRegionForThread) override
Superclass::InternalPixelType InternalPixelType
Superclass::InternalValueType InternalValueType
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.