OTB  10.0.0
Orfeo Toolbox
otbMultiToMonoChannelExtractROI.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 otbMultiToMonoChannelExtractROI_hxx
22 #define otbMultiToMonoChannelExtractROI_hxx
23 
25 
26 namespace otb
27 {
28 
32 template <class TInputPixelType, class TOutputPixelType>
34  : ExtractROIBase<VectorImage<TInputPixelType, 2>, Image<TOutputPixelType, 2>>(), m_Channel(1)
35 {
36  this->DynamicMultiThreadingOn();
37 }
38 
42 template <class TInputPixelType, class TOutputPixelType>
44 {
45  Superclass::PrintSelf(os, indent);
46 }
47 
57 template <class TInputPixelType, class TOutputPixelType>
59 {
60  typename Superclass::InputImageConstPointer inputPtr = this->GetInput();
61  // Bounds checking for the channel to process
62  if ((m_Channel <= 0) || (m_Channel > inputPtr->GetVectorLength()))
63  {
64  itkExceptionMacro(<< "otb::MultiToMonoChannelExtractROI::GenerateOutputInformation "
65  << "The selected channel must in the range [1;" << inputPtr->GetVectorLength() << "] "
66  << typeid(itk::ImageBase<InputImageDimension>*).name());
67  }
69 
70  // Calling the superclass method
71  Superclass::GenerateOutputInformation();
72 }
73 
74 template <class TInputPixelType, class TOutputPixelType>
76 {
77  itkDebugMacro(<< "Actually executing");
78 
79  // Get the input and output pointers
80  typename Superclass::InputImageConstPointer inputPtr = this->GetInput();
81  typename Superclass::OutputImagePointer outputPtr = this->GetOutput();
82 
83  // Define the portion of the input to walk for this thread
84  InputImageRegionType inputRegionForThread;
85  this->CallCopyOutputRegionToInputRegion(inputRegionForThread, outputRegionForThread);
86 
87  // Define the iterators.
88  typedef itk::ImageRegionIterator<OutputImageType> OutputIterator;
89  typedef itk::ImageRegionConstIterator<InputImageType> InputIterator;
90 
91  OutputIterator outIt(outputPtr, outputRegionForThread);
92  InputIterator inIt(inputPtr, inputRegionForThread);
93 
94  // Go through channels to process
95  const unsigned int channelIn(m_Channel - 1);
96 
97  while (!outIt.IsAtEnd())
98  {
99 
100  InputImagePixelType const& pixelInput = inIt.Get();
101  OutputImagePixelType const pixelOutput = static_cast<OutputValueType>(pixelInput[channelIn]);
102  outIt.Set(pixelOutput);
103  ++outIt;
104  ++inIt;
105  }
106 }
107 
108 } // end namespace otb
109 
110 #endif
Base class to extract area of images.
Creation of an "otb" image which contains metadata.
Definition: otbImage.h:92
void DynamicThreadedGenerateData(const OutputImageRegionType &outputRegionForThread) override
void PrintSelf(std::ostream &os, itk::Indent indent) const override
Creation of an "otb" vector image which contains metadata.
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.