OTB  9.0.0
Orfeo Toolbox
otbRadiometricMomentsImageFilter.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 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 }
43 
44 template <class TInputImage, class TOutputImage>
46 {
47  // call the superclass' implementation of this method
48  Superclass::GenerateInputRequestedRegion();
49 
50  // get pointers to the input and output
51  typename Superclass::InputImagePointer inputPtr = const_cast<TInputImage*>(this->GetInput());
52  typename Superclass::OutputImagePointer outputPtr = this->GetOutput();
53 
54  if (!inputPtr || !outputPtr)
55  {
56  return;
57  }
58  // get a copy of the input requested region (should equal the output
59  // requested region)
60  typename TInputImage::RegionType inputRequestedRegion;
61  inputRequestedRegion = inputPtr->GetRequestedRegion();
62 
63  // pad the input requested region by the operator radius
64  inputRequestedRegion.PadByRadius(m_Radius);
65 
66  // crop the input requested region at the input's largest possible region
67  if (inputRequestedRegion.Crop(inputPtr->GetLargestPossibleRegion()))
68  {
69  inputPtr->SetRequestedRegion(inputRequestedRegion);
70  return;
71  }
72  else
73  {
74  // Couldn't crop the region (requested region is outside the largest
75  // possible region). Throw an exception.
76 
77  // store what we tried to request (prior to trying to crop)
78  inputPtr->SetRequestedRegion(inputRequestedRegion);
79 
80  // build an exception
81  itk::InvalidRequestedRegionError e(__FILE__, __LINE__);
82  std::ostringstream msg;
83  msg << this->GetNameOfClass() << "::GenerateInputRequestedRegion()";
84  e.SetLocation(msg.str());
85  e.SetDescription("Requested region is (at least partially) outside the largest possible region.");
86  e.SetDataObject(inputPtr);
87  throw e;
88  }
89 }
90 
94 template <class TInputImage, class TOutputImage>
96 {
97  Superclass::GenerateOutputInformation();
98  this->GetOutput()->SetNumberOfComponentsPerPixel(4);
99 }
100 
104 template <class TInputImage, class TOutputImage>
106  itk::ThreadIdType threadId)
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  // support progress methods/callbacks
133  itk::ProgressReporter progress(this, threadId, outputRegionForThread.GetNumberOfPixels());
134 
135  // Process each of the boundary faces. These are N-d regions which border
136  // the edge of the buffer.
137  for (fit = faceList.begin(); fit != faceList.end(); ++fit)
138  {
139  neighInputIt = itk::ConstNeighborhoodIterator<TInputImage>(r, inputPtr, *fit);
140 
141  outputIt = itk::ImageRegionIterator<TOutputImage>(outputPtr, *fit);
142  neighInputIt.OverrideBoundaryCondition(&nbc);
143  neighInputIt.GoToBegin();
144 
145  while (!outputIt.IsAtEnd())
146  {
147  outputIt.Set(m_Functor(neighInputIt));
148 
149  ++neighInputIt;
150  ++outputIt;
151  progress.CompletedPixel();
152  }
153  }
154 }
155 
156 } // end namespace otb
157 
158 #endif
otb::RadiometricMomentsImageFilter::InputImagePointer
InputImageType::ConstPointer InputImagePointer
Definition: otbRadiometricMomentsImageFilter.h:63
otb::RadiometricMomentsImageFilter::GenerateOutputInformation
void GenerateOutputInformation(void) override
Definition: otbRadiometricMomentsImageFilter.hxx:95
otb::RadiometricMomentsImageFilter::OutputImagePointer
OutputImageType::Pointer OutputImagePointer
Definition: otbRadiometricMomentsImageFilter.h:68
otb::RadiometricMomentsImageFilter::ThreadedGenerateData
void ThreadedGenerateData(const OutputImageRegionType &outputRegionForThread, itk::ThreadIdType threadId) override
Definition: otbRadiometricMomentsImageFilter.hxx:105
otbRadiometricMomentsImageFilter.h
otb::RadiometricMomentsImageFilter::NeighborhoodIteratorType
itk::ConstNeighborhoodIterator< TInputImage > NeighborhoodIteratorType
Definition: otbRadiometricMomentsImageFilter.h:88
otb
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.
Definition: otbJoinContainer.h:32
otb::RadiometricMomentsImageFilter::RadiometricMomentsImageFilter
RadiometricMomentsImageFilter()
Definition: otbRadiometricMomentsImageFilter.hxx:37
otb::RadiometricMomentsImageFilter::GenerateInputRequestedRegion
void GenerateInputRequestedRegion(void) override
Definition: otbRadiometricMomentsImageFilter.hxx:45
otb::RadiometricMomentsImageFilter::RadiusType
NeighborhoodIteratorType::RadiusType RadiusType
Definition: otbRadiometricMomentsImageFilter.h:89
otb::RadiometricMomentsImageFilter::OutputImageRegionType
OutputImageType::RegionType OutputImageRegionType
Definition: otbRadiometricMomentsImageFilter.h:69