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