OTB  9.0.0
Orfeo Toolbox
otbFourierMellinDescriptorsImageFunction.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 otbFourierMellinDescriptorsImageFunction_hxx
22 #define otbFourierMellinDescriptorsImageFunction_hxx
23 
25 #include "itkConstNeighborhoodIterator.h"
26 #include "itkNumericTraits.h"
27 
28 namespace otb
29 {
30 
34 template <class TInputImage, class TCoordRep>
36 {
37  m_NeighborhoodRadius = 1;
38  m_Pmax = 3;
39  m_Qmax = 3;
40  m_Sigma = 0.5;
41 }
42 
43 template <class TInputImage, class TCoordRep>
44 void FourierMellinDescriptorsImageFunction<TInputImage, TCoordRep>::PrintSelf(std::ostream& os, itk::Indent indent) const
45 {
46  this->Superclass::PrintSelf(os, indent);
47  os << indent << " p indice maximum value : " << m_Pmax << std::endl;
48  os << indent << " q indice maximum value : " << m_Qmax << std::endl;
49  os << indent << " Neighborhood radius value : " << m_NeighborhoodRadius << std::endl;
50 }
51 
52 template <class TInputImage, class TCoordRep>
55 {
56  // Build Fourier-Mellin Harmonics Matrix
57  ComplexType coefs;
58  coefs.resize(m_Pmax + 1);
59  OutputType descriptors;
60  descriptors.resize(m_Pmax + 1);
61 
62  // Initialize moments
63  for (unsigned int p = 0; p <= m_Pmax; p++)
64  {
65  coefs.at(p).resize(m_Qmax + 1);
66  descriptors.at(p).resize(m_Qmax + 1);
67  for (unsigned int q = 0; q <= m_Qmax; q++)
68  {
69  coefs.at(p).at(q) = itk::NumericTraits<ScalarComplexType>::Zero;
70  descriptors.at(p).at(q) = itk::NumericTraits<ScalarRealType>::Zero;
71  }
72  }
73 
74  // Check for input image
75  if (!this->GetInputImage())
76  {
77  return descriptors;
78  }
79 
80  // Check for out of buffer
81  if (!this->IsInsideBuffer(index))
82  {
83  return descriptors;
84  }
85 
86  // Create an N-d neighborhood kernel, using a zeroflux boundary condition
87  typename InputImageType::SizeType kernelSize;
88  kernelSize.Fill(m_NeighborhoodRadius);
89 
90  itk::ConstNeighborhoodIterator<InputImageType> it(kernelSize, this->GetInputImage(), this->GetInputImage()->GetBufferedRegion());
91 
92  // Set the iterator at the desired location
93  it.SetLocation(index);
94 
95  // Walk the neighborhood
96  const unsigned int size = it.Size();
97  for (unsigned int i = 0; i < size; ++i)
98  {
99  // Retrieve value, and centered-reduced position
100  ScalarRealType value = static_cast<ScalarRealType>(it.GetPixel(i));
101  ScalarRealType x = static_cast<ScalarRealType>(it.GetOffset(i)[0]) / (2 * m_NeighborhoodRadius + 1);
102  ScalarRealType y = static_cast<ScalarRealType>(it.GetOffset(i)[1]) / (2 * m_NeighborhoodRadius + 1);
103 
104  // Build complex value
105  ScalarComplexType xplusiy(x, y), x2plusy2(x * x + y * y, 0.0);
106 
107  // Update cumulants
108  for (unsigned int p = 0; p <= m_Pmax; p++)
109  {
110  for (unsigned int q = 0; q <= m_Qmax; q++)
111  {
112  ScalarComplexType power(double(p - 2.0 + m_Sigma) / 2.0, -double(q) / 2.0);
113 
114  if (x != 0 || y != 0) // std::pow limitation
115  {
116  coefs.at(p).at(q) += std::pow(xplusiy, -static_cast<double>(p)) * std::pow(x2plusy2, power) * value;
117  }
118  }
119  }
120  }
121 
122  // Normalisation
123 
124  for (int p = m_Pmax; p >= 0; p--)
125  {
126  for (int q = m_Qmax; q >= 0; q--)
127  {
128  coefs.at(p).at(q) /= 2 * CONST_PI * coefs.at(0).at(0);
129 
130  descriptors.at(p).at(q) = std::sqrt((coefs.at(p).at(q).real() * coefs.at(p).at(q).real() + coefs.at(p).at(q).imag() * coefs.at(p).at(q).imag()));
131  }
132  }
133 
134  // Return result
135  return descriptors;
136 }
137 
138 } // namespace otb
139 
140 #endif
otb::CONST_PI
constexpr double CONST_PI
Definition: otbMath.h:49
otb::FourierMellinDescriptorsImageFunction::ScalarComplexType
std::complex< ScalarRealType > ScalarComplexType
Definition: otbFourierMellinDescriptorsImageFunction.h:88
otb
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.
Definition: otbJoinContainer.h:32
otbFourierMellinDescriptorsImageFunction.h
otb::FourierMellinDescriptorsImageFunction::ComplexType
std::vector< std::vector< ScalarComplexType > > ComplexType
Definition: otbFourierMellinDescriptorsImageFunction.h:89
otb::FourierMellinDescriptorsImageFunction::IndexType
Superclass::IndexType IndexType
Definition: otbFourierMellinDescriptorsImageFunction.h:82
otb::FourierMellinDescriptorsImageFunction::PrintSelf
void PrintSelf(std::ostream &os, itk::Indent indent) const override
Definition: otbFourierMellinDescriptorsImageFunction.hxx:44
otb::FourierMellinDescriptorsImageFunction::ScalarRealType
double ScalarRealType
Definition: otbFourierMellinDescriptorsImageFunction.h:87
otb::FourierMellinDescriptorsImageFunction::EvaluateAtIndex
OutputType EvaluateAtIndex(const IndexType &index) const override
Definition: otbFourierMellinDescriptorsImageFunction.hxx:54
otb::FourierMellinDescriptorsImageFunction::OutputType
Superclass::OutputType OutputType
Definition: otbFourierMellinDescriptorsImageFunction.h:85
otb::FourierMellinDescriptorsImageFunction::FourierMellinDescriptorsImageFunction
FourierMellinDescriptorsImageFunction()
Definition: otbFourierMellinDescriptorsImageFunction.hxx:35