OTB  9.0.0
Orfeo Toolbox
otbHistogramOfOrientedGradientCovariantImageFunction.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 otbHistogramOfOrientedGradientCovariantImageFunction_hxx
22 #define otbHistogramOfOrientedGradientCovariantImageFunction_hxx
23 
25 #include "itkConstNeighborhoodIterator.h"
26 #include "itkNumericTraits.h"
27 #include "otbMath.h"
28 
29 namespace otb
30 {
31 
35 template <class TInputImage, class TOutputPrecision, class TCoordRep>
37  : m_NeighborhoodRadius(8), m_NumberOfOrientationBins(18)
38 {
39 }
40 
41 template <class TInputImage, class TOutputPrecision, class TCoordRep>
43 {
44  this->Superclass::PrintSelf(os, indent);
45  os << indent << " Neighborhood radius value : " << m_NeighborhoodRadius << std::endl;
46 }
47 
48 template <class TInputImage, class TOutputPrecision, class TCoordRep>
51 {
52  // Build output vector
53  OutputType hog;
54 
55  // Check for input image
56  if (!this->GetInputImage())
57  {
58  return hog;
59  }
60 
61  // Check for out of buffer
62  if (!this->IsInsideBuffer(index))
63  {
64  return hog;
65  }
66 
67  // Create an N-d neighborhood kernel, using a zeroflux boundary condition
68  typename InputImageType::SizeType kernelSize;
69  kernelSize.Fill(m_NeighborhoodRadius);
70 
71  itk::ConstNeighborhoodIterator<InputImageType> it(kernelSize, this->GetInputImage(), this->GetInputImage()->GetBufferedRegion());
72 
73  // Set the iterator at the desired location
74  it.SetLocation(index);
75 
76  // Offset to be used in the loops
77  typename InputImageType::OffsetType offset;
78 
79  // Compute the center bin radius
80  double centerBinRadius = static_cast<double>(m_NeighborhoodRadius) / 2;
81 
82  // Define a gaussian kernel around the center location
83  double squaredRadius = m_NeighborhoodRadius * m_NeighborhoodRadius;
84  double squaredCenterBinRadius = centerBinRadius * centerBinRadius;
85  double squaredSigma = 0.25 * squaredRadius;
86 
87  // Build a global orientation histogram
88  std::vector<TOutputPrecision> globalOrientationHistogram(m_NumberOfOrientationBins, 0.);
89 
90  // Compute the orientation bin width
91  double orientationBinWidth = otb::CONST_2PI / m_NumberOfOrientationBins;
92 
93  // Build the global orientation histogram
94  for (int i = -(int)m_NeighborhoodRadius; i < (int)m_NeighborhoodRadius; ++i)
95  {
96  for (int j = -(int)m_NeighborhoodRadius; j < (int)m_NeighborhoodRadius; ++j)
97  {
98  // Check if the current pixel lies within a disc of radius m_NeighborhoodRadius
99  double currentSquaredRadius = i * i + j * j;
100  if (currentSquaredRadius < squaredRadius)
101  {
102  // If so, compute the gaussian weighting (this could be
103  // computed once for all for the sake of optimisation)
104  double gWeight = (1 / std::sqrt(otb::CONST_2PI * squaredSigma)) * std::exp(-currentSquaredRadius / (2 * squaredSigma));
105 
106  // Compute pixel location
107  offset[0] = i;
108  offset[1] = j;
109 
110  // Get the current gradient covariant value
111  InputPixelType gradient = it.GetPixel(offset);
112 
113  // Then, compute the gradient orientation
114  double angle = std::atan2(gradient[1], gradient[0]);
115 
116  // Also compute its magnitude
117  TOutputPrecision magnitude = std::sqrt(gradient[0] * gradient[0] + gradient[1] * gradient[1]);
118 
119  // Determine the bin index (shift of otb::CONST_PI since atan2 values
120  // lies in [-pi, pi]
121  unsigned int binIndex = std::floor((otb::CONST_PI + angle) / orientationBinWidth);
122 
123  // Handle special case where angle = pi, and binIndex is out-of-bound
124  if (binIndex == m_NumberOfOrientationBins)
125  binIndex = m_NumberOfOrientationBins - 1;
126 
127  // Cumulate values
128  globalOrientationHistogram[binIndex] += magnitude * gWeight;
129  }
130  }
131  }
132 
133  // Compute principal orientation
134  // TODO: Replace this with a stl algorithm
135  double maxOrientationHistogramValue = globalOrientationHistogram[0];
136  unsigned int maxOrientationHistogramBin = 0;
137 
138  for (unsigned int i = 1; i < m_NumberOfOrientationBins; ++i)
139  {
140  // Retrieve current value
141  double currentValue = globalOrientationHistogram[i];
142 
143  // Look for new maximum
144  if (maxOrientationHistogramValue < currentValue)
145  {
146  maxOrientationHistogramValue = currentValue;
147  maxOrientationHistogramBin = i;
148  }
149  }
150 
151  // Derive principal orientation
152  double principalOrientation = maxOrientationHistogramBin * orientationBinWidth - otb::CONST_PI;
153 
154  // Initialize the five spatial bins
155  std::vector<TOutputPrecision> centerHistogram(m_NumberOfOrientationBins, 0.);
156  std::vector<TOutputPrecision> upperLeftHistogram(m_NumberOfOrientationBins, 0.);
157  std::vector<TOutputPrecision> upperRightHistogram(m_NumberOfOrientationBins, 0.);
158  std::vector<TOutputPrecision> lowerLeftHistogram(m_NumberOfOrientationBins, 0.);
159  std::vector<TOutputPrecision> lowerRightHistogram(m_NumberOfOrientationBins, 0.);
160 
161  // Walk the image once more to fill the spatial bins
162  for (int i = -(int)m_NeighborhoodRadius; i < (int)m_NeighborhoodRadius; ++i)
163  {
164  for (int j = -(int)m_NeighborhoodRadius; j < (int)m_NeighborhoodRadius; ++j)
165  {
166  // Check if the current pixel lies within a disc of radius m_NeighborhoodRadius
167  double currentSquaredRadius = i * i + j * j;
168  if (currentSquaredRadius < squaredRadius)
169  {
170  // If so, compute the gaussian weighting (this could be
171  // computed once for all for the sake of optimisation)
172  double gWeight = (1 / std::sqrt(otb::CONST_2PI * squaredSigma)) * std::exp(-currentSquaredRadius / (2 * squaredSigma));
173 
174  // Compute pixel location
175  offset[0] = i;
176  offset[1] = j;
177 
178  // Get the current gradient covariant value
179  InputPixelType gradient = it.GetPixel(offset);
180 
181  // Then, compute the compensated gradient orientation
182  double angle = std::atan2(gradient[1], gradient[0]) - principalOrientation;
183 
184  // Angle is supposed to lie with [-pi, pi], so we ensure that
185  // compenstation did not introduce out-of-range values
186  if (angle > otb::CONST_PI)
187  {
188  angle -= otb::CONST_2PI;
189  }
190  else if (angle < -otb::CONST_PI)
191  {
192  angle += otb::CONST_2PI;
193  }
194 
195  // Also compute its magnitude
196  TOutputPrecision magnitude = std::sqrt(gradient[0] * gradient[0] + gradient[1] * gradient[1]);
197 
198  // Determine the bin index (shift of otb::CONST_PI since atan2 values
199  // lies in [-pi, pi]
200  unsigned int binIndex = std::floor((otb::CONST_PI + angle) / orientationBinWidth);
201 
202  if (binIndex == m_NumberOfOrientationBins)
203  binIndex = m_NumberOfOrientationBins - 1;
204 
205  // Compute the angular position
206  double angularPosition = std::atan2((double)j, (double)i) - principalOrientation;
207 
208  // Angle is supposed to lie within [-pi, pi], so we ensure that
209  // the compensation did not introduce out-of-range values
210  if (angularPosition > otb::CONST_PI)
211  {
212  angularPosition -= otb::CONST_2PI;
213  }
214  else if (angularPosition < -otb::CONST_PI)
215  {
216  angularPosition += otb::CONST_2PI;
217  }
218 
219  // Check if we lie in center bin
220  if (currentSquaredRadius < squaredCenterBinRadius)
221  {
222  centerHistogram[binIndex] += magnitude * gWeight;
223  }
224  else if (angularPosition > 0)
225  {
226  if (angularPosition < otb::CONST_PI_2)
227  {
228  upperRightHistogram[binIndex] += magnitude * gWeight;
229  }
230  else
231  {
232  upperLeftHistogram[binIndex] += magnitude * gWeight;
233  }
234  }
235  else
236  {
237  if (angularPosition > -otb::CONST_PI_2)
238  {
239  lowerRightHistogram[binIndex] += magnitude * gWeight;
240  }
241  else
242  {
243  lowerLeftHistogram[binIndex] += magnitude * gWeight;
244  }
245  }
246 
247  // Cumulate values
248  globalOrientationHistogram[binIndex] += magnitude * gWeight;
249  }
250  }
251  }
252 
253  // Build the final output
254  hog.push_back(centerHistogram);
255  hog.push_back(upperLeftHistogram);
256  hog.push_back(upperRightHistogram);
257  hog.push_back(lowerRightHistogram);
258  hog.push_back(lowerLeftHistogram);
259 
260  // Normalize each histogram
261  for (typename OutputType::iterator oIt = hog.begin(); oIt != hog.end(); ++oIt)
262  {
263  // Compute L2 norm
264  double squaredCumul = 1e-10;
265  for (typename std::vector<TOutputPrecision>::const_iterator vIt = oIt->begin(); vIt != oIt->end(); ++vIt)
266  {
267  squaredCumul += (*vIt) * (*vIt);
268  }
269  double scale = 1 / std::sqrt(squaredCumul);
270  // Apply normalisation factor
271  for (typename std::vector<TOutputPrecision>::iterator vIt = oIt->begin(); vIt != oIt->end(); ++vIt)
272  {
273  (*vIt) *= scale;
274  }
275  }
276 
277 
278  // Return result
279  return hog;
280 }
281 
282 } // namespace otb
283 
284 #endif
otb::CONST_PI
constexpr double CONST_PI
Definition: otbMath.h:49
otbMath.h
otb
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.
Definition: otbJoinContainer.h:32
otb::HistogramOfOrientedGradientCovariantImageFunction::OutputType
Superclass::OutputType OutputType
Definition: otbHistogramOfOrientedGradientCovariantImageFunction.h:95
otbHistogramOfOrientedGradientCovariantImageFunction.h
otb::CONST_2PI
constexpr double CONST_2PI
Definition: otbMath.h:55
otb::HistogramOfOrientedGradientCovariantImageFunction::InputPixelType
InputImageType::PixelType InputPixelType
Definition: otbHistogramOfOrientedGradientCovariantImageFunction.h:89
otb::HistogramOfOrientedGradientCovariantImageFunction::HistogramOfOrientedGradientCovariantImageFunction
HistogramOfOrientedGradientCovariantImageFunction()
Definition: otbHistogramOfOrientedGradientCovariantImageFunction.hxx:36
otb::HistogramOfOrientedGradientCovariantImageFunction::IndexType
Superclass::IndexType IndexType
Definition: otbHistogramOfOrientedGradientCovariantImageFunction.h:90
otb::CONST_PI_2
constexpr double CONST_PI_2
Definition: otbMath.h:50
otb::HistogramOfOrientedGradientCovariantImageFunction::PrintSelf
void PrintSelf(std::ostream &os, itk::Indent indent) const override
Definition: otbHistogramOfOrientedGradientCovariantImageFunction.hxx:42
otb::HistogramOfOrientedGradientCovariantImageFunction::EvaluateAtIndex
OutputType EvaluateAtIndex(const IndexType &index) const override
Definition: otbHistogramOfOrientedGradientCovariantImageFunction.hxx:50