OTB  10.0.0
Orfeo Toolbox
otbRadiometricMomentsImageFilter.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 otbRadiometricMomentsImageFilter_hxx
22 #define otbRadiometricMomentsImageFilter_hxx
23 
25 #include "itkImageRegionIterator.h"
26 #include "itkNeighborhoodAlgorithm.h"
27 #include "itkProgressReporter.h"
28 #include "itkZeroFluxNeumannBoundaryCondition.h"
29 #include "itkNeighborhoodAlgorithm.h"
30 
31 namespace otb
32 {
36 template <class TInputImage, class TOutputImage>
38 {
39  this->SetNumberOfRequiredInputs(1);
40  m_Radius.Fill(1);
41  this->DynamicMultiThreadingOn();
42 }
44 
45 template <class TInputImage, class TOutputImage>
47 {
48  // call the superclass' implementation of this method
49  Superclass::GenerateInputRequestedRegion();
50 
51  // get pointers to the input and output
52  typename Superclass::InputImagePointer inputPtr = const_cast<TInputImage*>(this->GetInput());
53  typename Superclass::OutputImagePointer outputPtr = this->GetOutput();
54 
55  if (!inputPtr || !outputPtr)
56  {
57  return;
58  }
59  // get a copy of the input requested region (should equal the output
60  // requested region)
61  typename TInputImage::RegionType inputRequestedRegion;
62  inputRequestedRegion = inputPtr->GetRequestedRegion();
63 
64  // pad the input requested region by the operator radius
65  inputRequestedRegion.PadByRadius(m_Radius);
66 
67  // crop the input requested region at the input's largest possible region
68  if (inputRequestedRegion.Crop(inputPtr->GetLargestPossibleRegion()))
69  {
70  inputPtr->SetRequestedRegion(inputRequestedRegion);
71  return;
72  }
73  else
74  {
75  // Couldn't crop the region (requested region is outside the largest
76  // possible region). Throw an exception.
77 
78  // store what we tried to request (prior to trying to crop)
79  inputPtr->SetRequestedRegion(inputRequestedRegion);
80 
81  // build an exception
82  itk::InvalidRequestedRegionError e(__FILE__, __LINE__);
83  std::ostringstream msg;
84  msg << this->GetNameOfClass() << "::GenerateInputRequestedRegion()";
85  e.SetLocation(msg.str());
86  e.SetDescription("Requested region is (at least partially) outside the largest possible region.");
87  e.SetDataObject(inputPtr);
88  throw e;
89  }
90 }
91 
95 template <class TInputImage, class TOutputImage>
97 {
98  Superclass::GenerateOutputInformation();
99  this->GetOutput()->SetNumberOfComponentsPerPixel(4);
100 }
101 
105 template <class TInputImage, class TOutputImage>
107 {
108  itk::ZeroFluxNeumannBoundaryCondition<TInputImage> nbc;
109 
110  // We use dynamic_cast since inputs are stored as DataObjects. The
111  // ImageToImageFilter::GetInput(int) always returns a pointer to a
112  // TInputImage so it cannot be used for the second input.
113  InputImagePointer inputPtr = dynamic_cast<const TInputImage*>(ProcessObjectType::GetInput(0));
114  OutputImagePointer outputPtr = this->GetOutput(0);
115 
116  RadiusType r;
117 
118  r[0] = m_Radius[0];
119  r[1] = m_Radius[1];
120 
121  NeighborhoodIteratorType neighInputIt;
122 
123  itk::ImageRegionIterator<TOutputImage> outputIt;
124 
125  // Find the data-set boundary "faces"
126  typename itk::NeighborhoodAlgorithm::ImageBoundaryFacesCalculator<TInputImage>::FaceListType faceList;
127  typename itk::NeighborhoodAlgorithm::ImageBoundaryFacesCalculator<TInputImage> bC;
128  faceList = bC(inputPtr, outputRegionForThread, r);
129 
130  typename itk::NeighborhoodAlgorithm::ImageBoundaryFacesCalculator<TInputImage>::FaceListType::iterator fit;
131 
132  // Process each of the boundary faces. These are N-d regions which border
133  // the edge of the buffer.
134  for (fit = faceList.begin(); fit != faceList.end(); ++fit)
135  {
136  neighInputIt = itk::ConstNeighborhoodIterator<TInputImage>(r, inputPtr, *fit);
137 
138  outputIt = itk::ImageRegionIterator<TOutputImage>(outputPtr, *fit);
139  neighInputIt.OverrideBoundaryCondition(&nbc);
140  neighInputIt.GoToBegin();
141 
142  while (!outputIt.IsAtEnd())
143  {
144  outputIt.Set(m_Functor(neighInputIt));
145 
146  ++neighInputIt;
147  ++outputIt;
148  }
149  }
150 }
151 
152 } // end namespace otb
153 
154 #endif
NeighborhoodIteratorType::RadiusType RadiusType
itk::ConstNeighborhoodIterator< TInputImage > NeighborhoodIteratorType
void DynamicThreadedGenerateData(const OutputImageRegionType &outputRegionForThread) override
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.