21 #ifndef otbLeeImageFilter_hxx
22 #define otbLeeImageFilter_hxx
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"
40 template <
class TInputImage,
class TOutputImage>
48 template <
class TInputImage,
class TOutputImage>
52 Superclass::GenerateInputRequestedRegion();
55 typename Superclass::InputImagePointer inputPtr =
const_cast<TInputImage*
>(this->GetInput());
56 typename Superclass::OutputImagePointer outputPtr = this->GetOutput();
58 if (!inputPtr || !outputPtr)
65 typename TInputImage::RegionType inputRequestedRegion;
66 inputRequestedRegion = inputPtr->GetRequestedRegion();
69 inputRequestedRegion.PadByRadius(m_Radius);
72 if (inputRequestedRegion.Crop(inputPtr->GetLargestPossibleRegion()))
74 inputPtr->SetRequestedRegion(inputRequestedRegion);
83 inputPtr->SetRequestedRegion(inputRequestedRegion);
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);
96 template <
class TInputImage,
class TOutputImage>
100 itk::ZeroFluxNeumannBoundaryCondition<InputImageType> nbc;
102 itk::ConstNeighborhoodIterator<InputImageType> bit;
103 itk::ImageRegionIterator<OutputImageType> it;
105 typename OutputImageType::Pointer output = this->GetOutput();
106 typename InputImageType::ConstPointer input = this->GetInput();
109 typename itk::NeighborhoodAlgorithm::ImageBoundaryFacesCalculator<InputImageType>::FaceListType faceList;
110 typename itk::NeighborhoodAlgorithm::ImageBoundaryFacesCalculator<InputImageType>::FaceListType::iterator fit;
112 itk::NeighborhoodAlgorithm::ImageBoundaryFacesCalculator<InputImageType> bC;
113 faceList = bC(input, outputRegionForThread, m_Radius);
116 itk::ProgressReporter progress(
this, threadId, outputRegionForThread.GetNumberOfPixels());
122 double Ci2, Cu2, w, E_I, I, Var_I, dPixel;
125 Cu2 = 1.0 / m_NbLooks;
129 for (fit = faceList.begin(); fit != faceList.end(); ++fit)
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);
140 while (!bit.IsAtEnd())
142 sum = itk::NumericTraits<InputRealType>::Zero;
143 sum2 = itk::NumericTraits<InputRealType>::Zero;
146 for (i = 0; i < neighborhoodSize; ++i)
148 dPixel =
static_cast<double>(bit.GetPixel(i));
151 E_I = sum /
static_cast<double>(neighborhoodSize);
153 for (i = 0; i < neighborhoodSize; ++i)
155 dPixel =
static_cast<double>(bit.GetPixel(i));
156 sum2 += (dPixel - E_I) * (dPixel - E_I);
158 Var_I = sum2 /
static_cast<double>(neighborhoodSize - 1);
160 I =
static_cast<double>(bit.GetCenterPixel());
162 Ci2 = Var_I / (E_I * E_I);
164 const double epsilon = 0.0000000001;
165 if (std::abs(E_I) < epsilon)
167 dPixel = itk::NumericTraits<OutputPixelType>::Zero;
169 else if (std::abs(Var_I) < epsilon)
180 dPixel = I * w + E_I * (1 - w);
189 progress.CompletedPixel();
197 template <
class TInputImage,
class TOutput>
200 Superclass::PrintSelf(os, indent);
201 os << indent <<
"Radius: " << m_Radius << std::endl;