OTB  9.0.0
Orfeo Toolbox
otbFrostImageFilter.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 otbFrostImageFilter_hxx
22 #define otbFrostImageFilter_hxx
23 
24 #include "otbFrostImageFilter.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  m_Deramp = 2;
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 FrostImageFilter<TInputImage, TOutputImage>::ThreadedGenerateData(const OutputImageRegionType& outputRegionForThread, itk::ThreadIdType threadId)
98 {
99  unsigned int i;
100  itk::ZeroFluxNeumannBoundaryCondition<InputImageType> nbc;
101  itk::ConstNeighborhoodIterator<InputImageType> bit;
102  typename itk::ConstNeighborhoodIterator<InputImageType>::OffsetType off;
103  itk::ImageRegionIterator<OutputImageType> it;
104 
105  // Allocate output
106  typename OutputImageType::Pointer output = this->GetOutput();
107  typename InputImageType::ConstPointer input = this->GetInput();
108 
109  // Find the data-set boundary "faces"
110  typename itk::NeighborhoodAlgorithm::ImageBoundaryFacesCalculator<InputImageType>::FaceListType faceList;
111  typename itk::NeighborhoodAlgorithm::ImageBoundaryFacesCalculator<InputImageType>::FaceListType::iterator fit;
112 
113  itk::NeighborhoodAlgorithm::ImageBoundaryFacesCalculator<InputImageType> bC;
114  faceList = bC(input, outputRegionForThread, m_Radius);
115 
116  // support progress methods/callbacks
117  itk::ProgressReporter progress(this, threadId, outputRegionForThread.GetNumberOfPixels());
118 
119  InputRealType sum;
120  InputRealType sum2;
121 
122  double Mean, Variance;
123  double Alpha;
124  double NormFilter;
125  double FrostFilter;
126  double CoefFilter;
127  double dPixel;
128 
129  // Process each of the boundary faces. These are N-d regions which border
130  // the edge of the buffer.
131  for (fit = faceList.begin(); fit != faceList.end(); ++fit)
132  {
133  bit = itk::ConstNeighborhoodIterator<InputImageType>(m_Radius, input, *fit);
134  unsigned int neighborhoodSize = bit.Size();
135  it = itk::ImageRegionIterator<OutputImageType>(output, *fit);
136  bit.OverrideBoundaryCondition(&nbc);
137 
138  bit.GoToBegin();
139  it.GoToBegin();
140 
141  while (!bit.IsAtEnd())
142  {
143  sum = itk::NumericTraits<InputRealType>::Zero;
144  sum2 = itk::NumericTraits<InputRealType>::Zero;
145 
146  for (i = 0; i < neighborhoodSize; ++i)
147  {
148  dPixel = static_cast<double>(bit.GetPixel(i));
149  sum += dPixel;
150  }
151  Mean = 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 - Mean) * (dPixel - Mean);
157  }
158  Variance = sum2 / double(neighborhoodSize - 1);
159 
160  const double epsilon = 0.0000000001;
161  if (std::abs(Mean) < epsilon)
162  {
163  dPixel = itk::NumericTraits<OutputPixelType>::Zero;
164  }
165  else if (std::abs(Variance) < epsilon)
166  {
167  dPixel = Mean;
168  }
169  else
170  {
171  Alpha = m_Deramp * Variance / (Mean * Mean);
172 
173  NormFilter = 0.0;
174  FrostFilter = 0.0;
175 
176  const int rad_x = m_Radius[0];
177  const int rad_y = m_Radius[1];
178 
179  for (int x = -rad_x; x <= rad_x; ++x)
180  {
181  for (int y = -rad_y; y <= rad_y; ++y)
182  {
183  double Dist = std::sqrt(static_cast<double>(x * x + y * y));
184  off[0] = x;
185  off[1] = y;
186 
187  dPixel = static_cast<double>(bit.GetPixel(off));
188 
189  CoefFilter = std::exp(-Alpha * Dist);
190  NormFilter += CoefFilter;
191  FrostFilter += (CoefFilter * dPixel);
192  }
193  }
194 
195  dPixel = FrostFilter / NormFilter;
196  }
197 
198  it.Set(static_cast<OutputPixelType>(dPixel));
199 
200  ++bit;
201  ++it;
202  progress.CompletedPixel();
203  }
204  }
205 }
206 
210 template <class TInputImage, class TOutput>
211 void FrostImageFilter<TInputImage, TOutput>::PrintSelf(std::ostream& os, itk::Indent indent) const
212 {
213  Superclass::PrintSelf(os, indent);
214  os << indent << "Radius: " << m_Radius << std::endl;
215 }
217 
218 } // end namespace otb
219 
220 #endif
otb::FrostImageFilter::OutputImageRegionType
OutputImageType::RegionType OutputImageRegionType
Definition: otbFrostImageFilter.h:81
otb::FrostImageFilter::OutputPixelType
OutputImageType::PixelType OutputPixelType
Definition: otbFrostImageFilter.h:75
otb
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.
Definition: otbJoinContainer.h:32
otb::FrostImageFilter::ThreadedGenerateData
void ThreadedGenerateData(const OutputImageRegionType &outputRegionForThread, itk::ThreadIdType threadId) override
Definition: otbFrostImageFilter.hxx:97
otbFrostImageFilter.h
otb::FrostImageFilter::PrintSelf
void PrintSelf(std::ostream &os, itk::Indent indent) const override
Definition: otbFrostImageFilter.hxx:211
otb::FrostImageFilter::GenerateInputRequestedRegion
void GenerateInputRequestedRegion() override
Definition: otbFrostImageFilter.hxx:49
otb::FrostImageFilter::FrostImageFilter
FrostImageFilter()
Definition: otbFrostImageFilter.hxx:41
otb::FrostImageFilter::InputRealType
itk::NumericTraits< InputPixelType >::RealType InputRealType
Definition: otbFrostImageFilter.h:78