OTB  10.0.0
Orfeo Toolbox
otbImageOfVectorsToMonoChannelExtractROI.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 otbImageOfVectorsToMonoChannelExtractROI_hxx
22 #define otbImageOfVectorsToMonoChannelExtractROI_hxx
23 
25 
26 namespace otb
27 {
28 
32 template <class TInputImage, class TOutputImage>
34 {
35  this->DynamicMultiThreadingOn();
36 }
37 
41 template <class TInputImage, class TOutputImage>
43 {
44  Superclass::PrintSelf(os, indent);
45 }
46 
56 template <class TInputImage, class TOutputImage>
58 {
59  typename Superclass::InputImageConstPointer inputPtr = this->GetInput();
60  // Analysis of processed channel
61  if ((m_Channel <= 0) || (m_Channel > InputImagePixelType::Dimension))
62  {
63  itkExceptionMacro(<< "otb::ExtractImImageOfVectorsToMonoChannelExtractROIageFilter::GenerateOutputInformation "
64  << " Channel must be in the following range: [1;" << InputImagePixelType::Dimension << "] "
65  << typeid(itk::ImageBase<InputImageDimension>*).name());
66  }
68 
69  // Call base class implementation
70  Superclass::GenerateOutputInformation();
71 }
72 
73 template <class TInputImage, class TOutputImage>
75 {
76  itkDebugMacro(<< "Actually executing");
77 
78  // Get the input and output pointers
79  typename Superclass::InputImageConstPointer inputPtr = this->GetInput();
80  typename Superclass::OutputImagePointer outputPtr = this->GetOutput();
81 
82  // Define the portion of the input to walk for this thread
83  InputImageRegionType inputRegionForThread;
84  this->CallCopyOutputRegionToInputRegion(inputRegionForThread, outputRegionForThread);
85 
86  // Define the iterators.
87  typedef itk::ImageRegionIterator<OutputImageType> OutputIterator;
88  typedef itk::ImageRegionConstIterator<InputImageType> InputIterator;
89 
90  OutputIterator outIt(outputPtr, outputRegionForThread);
91  InputIterator inIt(inputPtr, inputRegionForThread);
92 
93  // Loop through the processed channels
94  unsigned int channelIn(m_Channel - 1);
95 
96  InputImagePixelType pixelInput;
97  while (!outIt.IsAtEnd())
98  {
99  OutputImagePixelType pixelOutput;
100  pixelInput = inIt.Get();
101  pixelOutput = static_cast<OutputValueType>(pixelInput[channelIn]);
102  outIt.Set(pixelOutput);
103  ++outIt;
104  ++inIt;
105  }
106 }
107 
108 } // end namespace otb
109 
110 #endif
TOutputImage::RegionType OutputImageRegionType
TOutputImage::PixelType OutputImagePixelType
TInputImage::PixelType InputImagePixelType
TInputImage::RegionType InputImageRegionType
void PrintSelf(std::ostream &os, itk::Indent indent) const override
void DynamicThreadedGenerateData(const OutputImageRegionType &outputRegionForThread) override
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.