OTB  9.0.0
Orfeo Toolbox
otbComplexMomentPathFunction.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 otbComplexMomentPathFunction_hxx
22 #define otbComplexMomentPathFunction_hxx
23 
25 #include "itkImageRegionIterator.h"
26 #include "itkConstNeighborhoodIterator.h"
27 #include "otbMacro.h"
28 
29 #include <complex>
30 namespace otb
31 {
32 
36 template <class TInputPath, class TOutput, class TPrecision>
38 {
39  m_P = 0;
40  m_Q = 0;
41 }
42 
46 template <class TInputPath, class TOutput, class TPrecision>
47 void ComplexMomentPathFunction<TInputPath, TOutput, TPrecision>::PrintSelf(std::ostream& os, itk::Indent indent) const
48 {
49  this->Superclass::PrintSelf(os, indent);
50  os << indent << " p indice value : " << m_P << std::endl;
51  os << indent << " q indice value : " << m_Q << std::endl;
52 }
54 
55 template <class TInputPath, class TOutput, class TPrecision>
58 {
61  ComplexPrecisionType Result;
62  PrecisionType PixelValue(1.0);
63 
64  ValP = ComplexPrecisionType(1.0, 0.0);
65  ValQ = ComplexPrecisionType(1.0, 0.0);
66  unsigned int p = m_P;
67  while (p > 0)
68  {
69  ValP *= ComplexPrecisionType(index[0], index[1]);
70  --p;
71  }
72  unsigned int q = m_Q;
73  while (q > 0)
74  {
75  ValQ *= ComplexPrecisionType(index[0], -index[1]);
76  --q;
77  }
78 
79  Result = ValP * ValQ * ComplexPrecisionType(static_cast<PrecisionType>(PixelValue), 0.0);
80  return Result;
81 }
82 
83 template <class TInputPath, class TOutput, class TPrecision>
86 {
87  // Retrieve the vertex list
88  VertexListPointer vertexList = path.GetVertexList();
89  // Get the number of vertices in the path
90  unsigned int pathSize = vertexList->Size();
91 
92  // value will store the result
93  ComplexPrecisionType value = static_cast<ComplexPrecisionType>(0.0);
94 
95  // Check if we there are enough vertices in the path to actually
96  // compute something
97  if (pathSize < 2)
98  {
99  return static_cast<OutputType>(value);
100  }
101 
102  // First, we compute the centroid of the path so as to center the moment
103  typename VertexListType::ConstIterator it = vertexList->Begin();
104  VertexType centroid = it.Value();
105  ++it;
106 
107  // Cumulate points
108  while (it != vertexList->End())
109  {
110  centroid[0] += it.Value()[0];
111  centroid[1] += it.Value()[1];
112  ++it;
113  }
114 
115  // Normalize
116  centroid[0] /= static_cast<PrecisionType>(pathSize);
117  centroid[1] /= static_cast<PrecisionType>(pathSize);
118 
119  // Second, we integrate along the edge
120  it = vertexList->Begin();
121 
122  VertexType source = it.Value();
123  source[0] -= centroid[0];
124  source[1] -= centroid[1];
125  ++it;
126 
127  PrecisionType ds;
128  VertexType dest;
129 
130  // This variable will be used to normalize the moment
131  PrecisionType norm = 0.;
132 
133  while (it != vertexList->End())
134  {
135  dest = it.Value();
136 
137  // Get source and destination coordinates
138  dest[0] -= centroid[0];
139  dest[1] -= centroid[1];
140 
141  // Don't forget the ds part of the integration process
142  ds = std::sqrt(std::pow(dest[0] - source[0], 2.) + std::pow(dest[1] - source[1], 2.));
143  norm += ds;
144  value += ds * EvaluateComplexMomentAtIndex(source);
145  source = dest;
146  ++it;
147  }
148  // Close the loop
149  dest = vertexList->Begin().Value();
150  dest[0] -= centroid[0];
151  dest[1] -= centroid[1];
152  ds = std::sqrt(std::pow(dest[0] - source[0], 2.) + std::pow(dest[1] - source[1], 2.));
153  norm += ds;
154  value += EvaluateComplexMomentAtIndex(source) * ds;
155  norm = std::pow(norm, ((PrecisionType)m_P + (PrecisionType)m_Q) / 2.);
156 
157  // Normalize with edge perimeter
158  value /= norm;
159 
160  return static_cast<OutputType>(value);
161 }
162 
163 template <class TInputPath, class TOutput, class TPrecision>
165 {
166  if (!this->GetInputPath())
167  {
168  otbMsgDevMacro(<< "Pb with GetInputPath");
169  return static_cast<OutputType>(ComplexPrecisionType(itk::NumericTraits<PrecisionType>::Zero, itk::NumericTraits<PrecisionType>::Zero));
170  }
171 
172  OutputType Result = Evaluate(*(this->GetInputPath()));
173 
174  return Result;
175 }
176 
177 } // namespace otb
178 
179 #endif
otb::ComplexMomentPathFunction::PrecisionType
Superclass::PrecisionType PrecisionType
Definition: otbComplexMomentPathFunction.h:85
otb::ComplexMomentPathFunction::EvaluateComplexMomentAtIndex
ComplexPrecisionType EvaluateComplexMomentAtIndex(VertexType index) const
Definition: otbComplexMomentPathFunction.hxx:57
otb::ComplexMomentPathFunction::VertexListPointer
VertexListType::ConstPointer VertexListPointer
Definition: otbComplexMomentPathFunction.h:75
otbComplexMomentPathFunction.h
otb::ComplexMomentPathFunction::PrintSelf
void PrintSelf(std::ostream &os, itk::Indent indent) const override
Definition: otbComplexMomentPathFunction.hxx:47
otb::ComplexMomentPathFunction::Evaluate
virtual OutputType Evaluate() const
Definition: otbComplexMomentPathFunction.hxx:164
otb::ComplexMomentPathFunction::OutputType
Superclass::OutputType OutputType
Definition: otbComplexMomentPathFunction.h:80
otb::ComplexMomentPathFunction::ComplexPrecisionType
std::complex< PrecisionType > ComplexPrecisionType
Definition: otbComplexMomentPathFunction.h:88
otb
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.
Definition: otbJoinContainer.h:32
otbMacro.h
otb::ComplexMomentPathFunction::ComplexMomentPathFunction
ComplexMomentPathFunction()
Definition: otbComplexMomentPathFunction.hxx:37
otbMsgDevMacro
#define otbMsgDevMacro(x)
Definition: otbMacro.h:64
otb::ComplexMomentPathFunction::VertexType
PathType::ContinuousIndexType VertexType
Definition: otbComplexMomentPathFunction.h:73
otb::ComplexMomentPathFunction::PathType
Superclass::PathType PathType
Definition: otbComplexMomentPathFunction.h:68