OTB  9.0.0
Orfeo Toolbox
otbExtractROI.hxx
Go to the documentation of this file.
1 /*
2  * Copyright (C) 2005-2022 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 otbExtractROI_hxx
22 #define otbExtractROI_hxx
23 
24 #include "otbExtractROI.h"
25 #include "itkImageRegionIterator.h"
26 #include "itkProgressReporter.h"
27 
28 namespace otb
29 {
30 
34 template <class TInputPixel, class TOutputPixel>
35 ExtractROI<TInputPixel, TOutputPixel>::ExtractROI() //: ExtractROIBase< itk::Image<TInputPixel, VImageDimension> , itk::Image<TOutputPixel, VImageDimension> >()
36 {
37 }
38 
42 template <class TInputPixel, class TOutputPixel>
43 void ExtractROI<TInputPixel, TOutputPixel>::PrintSelf(std::ostream& os, itk::Indent indent) const
44 {
45  Superclass::PrintSelf(os, indent);
46 }
47 
57 template <class TInputPixel, class TOutputPixel>
59 {
60  // Call to the base class method
61  Superclass::GenerateOutputInformation();
62 }
63 
64 template <class TInputPixel, class TOutputPixel>
65 void ExtractROI<TInputPixel, TOutputPixel>::ThreadedGenerateData(const OutputImageRegionType& outputRegionForThread, itk::ThreadIdType threadId)
66 {
67  itkDebugMacro(<< "Actually executing");
68 
69  // Get the input and output pointers
70  typename Superclass::InputImageConstPointer inputPtr = this->GetInput();
71  typename Superclass::OutputImagePointer outputPtr = this->GetOutput();
72 
73  // support progress methods/callbacks
74  itk::ProgressReporter progress(this, threadId, outputRegionForThread.GetNumberOfPixels());
75 
76  // Define the portion of the input to walk for this thread
77  InputImageRegionType inputRegionForThread;
78  this->CallCopyOutputRegionToInputRegion(inputRegionForThread, outputRegionForThread);
79 
80  // Define the iterators.
81  typedef itk::ImageRegionIterator<OutputImageType> OutputIterator;
82  typedef itk::ImageRegionConstIterator<InputImageType> InputIterator;
83 
84  OutputIterator outIt(outputPtr, outputRegionForThread);
85  InputIterator inIt(inputPtr, inputRegionForThread);
86 
87  // walk the output region, and sample the input image
88  while (!outIt.IsAtEnd())
89  {
90  // copy the input pixel to the output
91  outIt.Set(inIt.Get());
92  ++outIt;
93  ++inIt;
94  progress.CompletedPixel();
95  }
96 }
97 
98 } // end namespace otb
99 
100 #endif
otb::ExtractROI::ExtractROI
ExtractROI()
Definition: otbExtractROI.hxx:35
otb
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.
Definition: otbJoinContainer.h:32
otb::ExtractROI::PrintSelf
void PrintSelf(std::ostream &os, itk::Indent indent) const override
Definition: otbExtractROI.hxx:43
otb::ExtractROI::GenerateOutputInformation
void GenerateOutputInformation() override
Definition: otbExtractROI.hxx:58
otb::ExtractROI::ThreadedGenerateData
void ThreadedGenerateData(const OutputImageRegionType &outputRegionForThread, itk::ThreadIdType threadId) override
Definition: otbExtractROI.hxx:65
otb::ExtractROI::OutputImageRegionType
OutputImageType::RegionType OutputImageRegionType
Definition: otbExtractROI.h:63
otbExtractROI.h
otb::ExtractROI::InputImageRegionType
InputImageType::RegionType InputImageRegionType
Definition: otbExtractROI.h:64