OTB  9.0.0
Orfeo Toolbox
otbKuanImageFilter.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 otbKuanImageFilter_hxx
22 #define otbKuanImageFilter_hxx
23 
24 #include "otbKuanImageFilter.h"
25 
26 #include "itkDataObject.h"
27 #include "itkConstNeighborhoodIterator.h"
28 #include "itkNeighborhoodInnerProduct.h"
29 #include "itkImageRegionIterator.h"
30 #include "itkNeighborhoodAlgorithm.h"
31 #include "itkOffset.h"
32 #include "itkProgressReporter.h"
33 
34 namespace otb
35 {
36 
40 template <class TInputImage, class TOutputImage>
42 {
43  m_Radius.Fill(1);
44  SetNbLooks(1.0);
45 }
47 
48 template <class TInputImage, class TOutputImage>
50 {
51  // call the superclass' implementation of this method
52  Superclass::GenerateInputRequestedRegion();
53 
54  // get pointers to the input and output
55  typename Superclass::InputImagePointer inputPtr = const_cast<TInputImage*>(this->GetInput());
56  typename Superclass::OutputImagePointer outputPtr = this->GetOutput();
57 
58  if (!inputPtr || !outputPtr)
59  {
60  return;
61  }
62 
63  // get a copy of the input requested region (should equal the output
64  // requested region)
65  typename TInputImage::RegionType inputRequestedRegion;
66  inputRequestedRegion = inputPtr->GetRequestedRegion();
67 
68  // pad the input requested region by the operator radius
69  inputRequestedRegion.PadByRadius(m_Radius);
70 
71  // crop the input requested region at the input's largest possible region
72  if (inputRequestedRegion.Crop(inputPtr->GetLargestPossibleRegion()))
73  {
74  inputPtr->SetRequestedRegion(inputRequestedRegion);
75  return;
76  }
77  else
78  {
79  // Couldn't crop the region (requested region is outside the largest
80  // possible region). Throw an exception.
81 
82  // store what we tried to request (prior to trying to crop)
83  inputPtr->SetRequestedRegion(inputRequestedRegion);
84 
85  // build an exception
86  itk::InvalidRequestedRegionError e(__FILE__, __LINE__);
87  std::ostringstream msg;
88  msg << static_cast<const char*>(this->GetNameOfClass()) << "::GenerateInputRequestedRegion()";
89  e.SetLocation(msg.str());
90  e.SetDescription("Requested region is (at least partially) outside the largest possible region.");
91  e.SetDataObject(inputPtr);
92  throw e;
93  }
94 }
95 
96 template <class TInputImage, class TOutputImage>
97 void KuanImageFilter<TInputImage, TOutputImage>::ThreadedGenerateData(const OutputImageRegionType& outputRegionForThread, itk::ThreadIdType threadId)
98 {
99  unsigned int i;
100  itk::ZeroFluxNeumannBoundaryCondition<InputImageType> nbc;
101 
102  itk::ConstNeighborhoodIterator<InputImageType> bit;
103  itk::ImageRegionIterator<OutputImageType> it;
104 
105  typename OutputImageType::Pointer output = this->GetOutput();
106  typename InputImageType::ConstPointer input = this->GetInput();
107 
108  // Find the data-set boundary "faces"
109  typename itk::NeighborhoodAlgorithm::ImageBoundaryFacesCalculator<InputImageType>::FaceListType faceList;
110  typename itk::NeighborhoodAlgorithm::ImageBoundaryFacesCalculator<InputImageType>::FaceListType::iterator fit;
111 
112  itk::NeighborhoodAlgorithm::ImageBoundaryFacesCalculator<InputImageType> bC;
113  faceList = bC(input, outputRegionForThread, m_Radius);
114 
115  // support progress methods/callbacks
116  itk::ProgressReporter progress(this, threadId, outputRegionForThread.GetNumberOfPixels());
117 
118  // InputRealType pixel;
119  InputRealType sum;
120  InputRealType sum2;
121 
122  double Ci2, Cu2, w, E_I, I, Var_I, dPixel;
123 
124  // Compute the ratio using the number of looks
125  Cu2 = 1.0 / m_NbLooks;
126 
127  // Process each of the boundary faces. These are N-d regions which border
128  // the edge of the buffer.
129  for (fit = faceList.begin(); fit != faceList.end(); ++fit)
130  {
131  bit = itk::ConstNeighborhoodIterator<InputImageType>(m_Radius, input, *fit);
132  const unsigned int neighborhoodSize = bit.Size();
133  it = itk::ImageRegionIterator<OutputImageType>(output, *fit);
134  bit.OverrideBoundaryCondition(&nbc);
135 
136  bit.GoToBegin();
137  it.GoToBegin();
138 
139 
140  while (!bit.IsAtEnd())
141  {
142  sum = itk::NumericTraits<InputRealType>::Zero;
143  sum2 = itk::NumericTraits<InputRealType>::Zero;
144 
145  // Parcours du voisinage
146  for (i = 0; i < neighborhoodSize; ++i)
147  {
148  dPixel = static_cast<double>(bit.GetPixel(i));
149  sum += dPixel;
150  }
151  E_I = sum / static_cast<double>(neighborhoodSize);
152 
153  for (i = 0; i < neighborhoodSize; ++i)
154  {
155  dPixel = static_cast<double>(bit.GetPixel(i));
156  sum2 += (dPixel - E_I) * (dPixel - E_I);
157  }
158  Var_I = sum2 / static_cast<double>(neighborhoodSize - 1);
159 
160  I = static_cast<double>(bit.GetCenterPixel());
161 
162  Ci2 = Var_I / (E_I * E_I);
163 
164  const double epsilon = 0.0000000001;
165  if (std::abs(E_I) < epsilon)
166  {
167  dPixel = itk::NumericTraits<OutputPixelType>::Zero;
168  }
169  else if (std::abs(Var_I) < epsilon)
170  {
171  dPixel = E_I;
172  }
173  else if (Ci2 < Cu2)
174  {
175  dPixel = E_I;
176  }
177  else
178  {
179  w = (1 - Cu2 / Ci2) / (1 + Cu2);
180  dPixel = I * w + E_I * (1 - w);
181  }
182 
183  // set the weighted value
184  it.Set(static_cast<OutputPixelType>(dPixel));
185 
186  ++bit;
187  ++it;
188 
189  progress.CompletedPixel();
190  }
191  }
192 }
193 
197 template <class TInputImage, class TOutput>
198 void KuanImageFilter<TInputImage, TOutput>::PrintSelf(std::ostream& os, itk::Indent indent) const
199 {
200  Superclass::PrintSelf(os, indent);
201  os << indent << "Radius: " << m_Radius << std::endl;
202 }
204 
205 } // end namespace otb
206 
207 #endif
otb::KuanImageFilter::ThreadedGenerateData
void ThreadedGenerateData(const OutputImageRegionType &outputRegionForThread, itk::ThreadIdType threadId) override
Definition: otbKuanImageFilter.hxx:97
otb
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.
Definition: otbJoinContainer.h:32
otbKuanImageFilter.h
otb::KuanImageFilter::PrintSelf
void PrintSelf(std::ostream &os, itk::Indent indent) const override
Definition: otbKuanImageFilter.hxx:198
otb::KuanImageFilter::OutputImageRegionType
OutputImageType::RegionType OutputImageRegionType
Definition: otbKuanImageFilter.h:70
otb::KuanImageFilter::InputRealType
itk::NumericTraits< InputPixelType >::RealType InputRealType
Definition: otbKuanImageFilter.h:68
otb::KuanImageFilter::OutputPixelType
OutputImageType::PixelType OutputPixelType
Definition: otbKuanImageFilter.h:67
otb::KuanImageFilter::KuanImageFilter
KuanImageFilter()
Definition: otbKuanImageFilter.hxx:41
otb::KuanImageFilter::GenerateInputRequestedRegion
void GenerateInputRequestedRegion() override
Definition: otbKuanImageFilter.hxx:49