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