OTB  9.0.0
Orfeo Toolbox
otbEigenvalueLikelihoodMaximisation.hxx
Go to the documentation of this file.
1 /*
2  * Copyright (C) 1999-2011 Insight Software Consortium
3  * Copyright (C) 2005-2022 Centre National d'Etudes Spatiales (CNES)
4  *
5  * This file is part of Orfeo Toolbox
6  *
7  * https://www.orfeo-toolbox.org/
8  *
9  * Licensed under the Apache License, Version 2.0 (the "License");
10  * you may not use this file except in compliance with the License.
11  * You may obtain a copy of the License at
12  *
13  * http://www.apache.org/licenses/LICENSE-2.0
14  *
15  * Unless required by applicable law or agreed to in writing, software
16  * distributed under the License is distributed on an "AS IS" BASIS,
17  * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
18  * See the License for the specific language governing permissions and
19  * limitations under the License.
20  */
21 
22 #ifndef otbEigenvalueLikelihoodMaximisation_hxx
23 #define otbEigenvalueLikelihoodMaximisation_hxx
24 
26 
27 #include <algorithm>
28 
29 namespace otb
30 {
31 
32 template <class TPrecision>
34 {
35 }
36 
37 template <class TInputImage>
39 {
40  // TODO check size
41  const unsigned int nbBands = m_Covariance.rows();
42 
43  // Compute diagonalisation of covariance and correlation
44  vnl_symmetric_eigensystem<PrecisionType> eigenK(m_Covariance);
45  VectorType eigenCovariance = eigenK.D.diagonal();
46  std::sort(eigenCovariance.begin(), eigenCovariance.end());
47  eigenCovariance.flip();
48 
49  vnl_symmetric_eigensystem<PrecisionType> eigenR(m_Correlation);
50  VectorType eigenCorrelation = eigenR.D.diagonal();
51  std::sort(eigenCorrelation.begin(), eigenCorrelation.end());
52  eigenCorrelation.flip();
53 
54  // Compute likelihood log
55  m_Likelihood.set_size(nbBands);
56  const double coef = 2.0 / m_NumberOfPixels;
57 
58  for (unsigned int i = 0; i < nbBands; ++i)
59  {
60  const unsigned int nl = nbBands - i;
61  VectorType sigma(nl), t(nl);
62 
63  for (unsigned int j = 0; j < nl; ++j)
64  {
65  PrecisionType r = eigenCorrelation[j + i];
66  PrecisionType k = eigenCovariance[j + i];
67 
68  sigma[j] = coef * (r * r + k * k);
69  // std::cout << "sigma[" << j << "]=" << sigma[j] << std::endl;
70  t[j] = (r - k) * (r - k) / sigma[j];
71  // std::cout << "t[" << j <<"]=" << t[j] << std::endl;
72  sigma[j] = std::log(sigma[j]);
73  }
74  m_Likelihood(i) = -0.5 * t.sum() - 0.5 * sigma.sum();
75  }
76 
77  // Extract first local maximum
78  // double max = m_Likelihood[0];
79  unsigned int iMax = 0;
80  for (unsigned int i = 1; i < m_Likelihood.size() - 1; ++i)
81  {
82  if ((m_Likelihood[i] > m_Likelihood[i - 1]) && (m_Likelihood[i] > m_Likelihood[i + 1]))
83  {
84  // max = m_Likelihood[i];
85  iMax = i;
86  break;
87  }
88  }
89  m_NumberOfEndmembers = iMax;
90 }
91 
92 template <class TImage>
93 void EigenvalueLikelihoodMaximisation<TImage>::PrintSelf(std::ostream& os, itk::Indent indent) const
94 {
95  Superclass::PrintSelf(os, indent);
96 
97  os << indent << "Covariance: " << m_Covariance << std::endl;
98  os << indent << "Correlation: " << m_Correlation << std::endl;
99  os << indent << "NumberOfEndmembers: " << m_NumberOfEndmembers << std::endl;
100  os << indent << "Likelihood: " << m_Likelihood << std::endl;
101 }
102 
103 } // end namespace otb
104 #endif
otb::EigenvalueLikelihoodMaximisation::PrintSelf
void PrintSelf(std::ostream &os, itk::Indent indent) const override
Definition: otbEigenvalueLikelihoodMaximisation.hxx:93
otb
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.
Definition: otbJoinContainer.h:32
otbEigenvalueLikelihoodMaximisation.h
otb::EigenvalueLikelihoodMaximisation::PrecisionType
TPrecision PrecisionType
Definition: otbEigenvalueLikelihoodMaximisation.h:70
otb::EigenvalueLikelihoodMaximisation::VectorType
vnl_vector< PrecisionType > VectorType
Definition: otbEigenvalueLikelihoodMaximisation.h:75
otb::EigenvalueLikelihoodMaximisation::Compute
void Compute()
Definition: otbEigenvalueLikelihoodMaximisation.hxx:38
otb::EigenvalueLikelihoodMaximisation::EigenvalueLikelihoodMaximisation
EigenvalueLikelihoodMaximisation()
Definition: otbEigenvalueLikelihoodMaximisation.hxx:33