OTB  10.0.0
Orfeo Toolbox
otbUnaryFunctorNeighborhoodWithOffsetImageFilter.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 otbUnaryFunctorNeighborhoodWithOffsetImageFilter_hxx
22 #define otbUnaryFunctorNeighborhoodWithOffsetImageFilter_hxx
23 
24 #include "otbMacro.h" //for
26 #include "itkImageRegionIterator.h"
27 #include "itkNeighborhoodAlgorithm.h"
28 #include "itkProgressReporter.h"
29 #include "itkZeroFluxNeumannBoundaryCondition.h"
30 #include "itkNeighborhoodAlgorithm.h"
31 
32 namespace otb
33 {
37 template <class TInputImage, class TOutputImage, class TFunction>
39 {
40  this->DynamicMultiThreadingOff();
41  this->SetNumberOfRequiredInputs(1);
42  m_Radius.Fill(1);
43  m_Offset.Fill(1);
44  m_FunctorList.clear();
45 }
47 
48 template <class TInputImage, class TOutputImage, class TFunction>
50 {
51  Superclass::BeforeThreadedGenerateData();
52 
53  for (itk::ThreadIdType i = 0; i < this->GetNumberOfWorkUnits(); ++i)
54  {
55  m_FunctorList.push_back(m_Functor);
56  }
57 }
58 
59 template <class TInputImage, class TOutputImage, class TFunction>
61 {
62  // call the superclass' implementation of this method
63  Superclass::GenerateInputRequestedRegion();
64 
65  // get pointers to the input and output
66  typename Superclass::InputImagePointer inputPtr = const_cast<TInputImage*>(this->GetInput());
67  typename Superclass::OutputImagePointer outputPtr = this->GetOutput();
68 
69  if (!inputPtr || !outputPtr)
70  {
71  return;
72  }
73  // get a copy of the input requested region (should equal the output
74  // requested region)
75  typename TInputImage::RegionType inputRequestedRegion;
76  inputRequestedRegion = inputPtr->GetRequestedRegion();
77 
78  // pad the input requested region by the operator radius
79  InputImageSizeType maxRad;
80  maxRad[0] = m_Radius[0] + std::abs(m_Offset[0]);
81  maxRad[1] = m_Radius[1] + std::abs(m_Offset[1]);
82  inputRequestedRegion.PadByRadius(maxRad);
83 
84  // crop the input requested region at the input's largest possible region
85  if (inputRequestedRegion.Crop(inputPtr->GetLargestPossibleRegion()))
86  {
87  inputPtr->SetRequestedRegion(inputRequestedRegion);
88  return;
89  }
90  else
91  {
92  // Couldn't crop the region (requested region is outside the largest
93  // possible region). Throw an exception.
94 
95  // store what we tried to request (prior to trying to crop)
96  inputPtr->SetRequestedRegion(inputRequestedRegion);
97 
98  // build an exception
99  itk::InvalidRequestedRegionError e(__FILE__, __LINE__);
100  std::ostringstream msg;
101  msg << this->GetNameOfClass() << "::GenerateInputRequestedRegion()";
102  e.SetLocation(msg.str());
103  e.SetDescription("Requested region is (at least partially) outside the largest possible region.");
104  e.SetDataObject(inputPtr);
105  throw e;
106  }
107 }
108 
112 template <class TInputImage, class TOutputImage, class TFunction>
114  const OutputImageRegionType& outputRegionForThread, itk::ThreadIdType threadId)
115 {
116  itk::ZeroFluxNeumannBoundaryCondition<TInputImage> nbc;
117  itk::ZeroFluxNeumannBoundaryCondition<TInputImage> nbcOff;
118  // We use dynamic_cast since inputs are stored as DataObjects. The
119  // ImageToImageFilter::GetInput(int) always returns a pointer to a
120  // TInputImage so it cannot be used for the second input.
121  InputImagePointer inputPtr = dynamic_cast<const TInputImage*>(ProcessObjectType::GetInput(0));
122  OutputImagePointer outputPtr = this->GetOutput(0);
124 
125  itk::ImageRegionIterator<TOutputImage> outputIt;
126 
127  // Neighborhood+offset iterator
128  RadiusType rOff;
129  rOff[0] = m_Radius[0] + std::abs(m_Offset[0]);
130  rOff[1] = m_Radius[1] + std::abs(m_Offset[1]);
131  NeighborhoodIteratorType neighInputOffIt;
132  // Find the data-set boundary "faces"
133  typename itk::NeighborhoodAlgorithm::ImageBoundaryFacesCalculator<TInputImage>::FaceListType faceListOff;
134  typename itk::NeighborhoodAlgorithm::ImageBoundaryFacesCalculator<TInputImage> bCOff;
135  faceListOff = bCOff(inputPtr, outputRegionForThread, rOff);
136  typename itk::NeighborhoodAlgorithm::ImageBoundaryFacesCalculator<TInputImage>::FaceListType::iterator fitOff;
137 
138  // support progress methods/callbacks
139  itk::ProgressReporter progress(this, threadId, outputRegionForThread.GetNumberOfPixels());
140 
141  // Process each of the boundary faces. These are N-d regions which border
142  // the edge of the buffer.
143  fitOff = faceListOff.begin();
144 
145  while (fitOff != faceListOff.end())
146  {
147  // Neighborhood+offset iterator
148  neighInputOffIt = itk::ConstNeighborhoodIterator<TInputImage>(rOff, inputPtr, *fitOff);
149  neighInputOffIt.OverrideBoundaryCondition(&nbcOff);
150  neighInputOffIt.GoToBegin();
151 
152  outputIt = itk::ImageRegionIterator<TOutputImage>(outputPtr, *fitOff);
153 
154  while (!outputIt.IsAtEnd())
155  {
156 
157  outputIt.Set(m_FunctorList[threadId](neighInputOffIt.GetNeighborhood()));
158 
159  ++neighInputOffIt;
160  ++outputIt;
161  progress.CompletedPixel();
162  }
163  ++fitOff;
164  }
165 }
166 
167 } // end namespace otb
168 
169 #endif
void ThreadedGenerateData(const OutputImageRegionType &outputRegionForThread, itk::ThreadIdType threadId) override
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.