OTB  10.0.0
Orfeo Toolbox
otbVirtualDimensionality.h
Go to the documentation of this file.
1 /*
2  * Copyright (C) 1999-2011 Insight Software Consortium
3  * Copyright (C) 2005-2024 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 otbVirtualDimensionality_h
23 #define otbVirtualDimensionality_h
24 
25 #include "itkObjectFactory.h"
26 #include "itkLightObject.h"
27 #include "vnl/vnl_vector.h"
28 #include "vnl/vnl_matrix.h"
29 #include "vnl/algo/vnl_symmetric_eigensystem.h"
30 
31 namespace otb
32 {
33 
58 template <class TPrecision>
59 class ITK_EXPORT VirtualDimensionality : public itk::LightObject
60 {
61 public:
64  typedef itk::LightObject Superclass;
65  typedef itk::SmartPointer<Self> Pointer;
66  typedef itk::SmartPointer<const Self> ConstPointer;
67 
69  itkNewMacro(Self);
70 
72  itkTypeMacro(VirtualDimensionality, itk::LightObject);
73 
75  typedef TPrecision PrecisionType;
76 
77  typedef vnl_vector<PrecisionType> VectorType;
78  typedef vnl_matrix<PrecisionType> MatrixType;
79 
80  void SetCovariance(const MatrixType& m)
81  {
82  m_Covariance = m;
83  }
84 
85  void SetCorrelation(const MatrixType& m)
86  {
87  m_Correlation = m;
88  }
89 
90  void SetNumberOfPixels(unsigned int n)
91  {
92  m_NumberOfPixels = n;
93  }
94 
95  void SetFAR(double falseAlarmRate)
96  {
97  if (falseAlarmRate < 0)
98  falseAlarmRate = 0;
99 
100  if (falseAlarmRate > 1)
101  falseAlarmRate = 1;
102 
103  m_FAR = falseAlarmRate;
104  }
105 
106  double GetFAR()
107  {
108  return m_FAR;
109  }
110 
111  void Compute();
112 
113  unsigned int GetNumberOfEndmembers()
114  {
115  return m_NumberOfEndmembers;
116  }
117 
118 protected:
121  {
122  }
123  void PrintSelf(std::ostream& os, itk::Indent indent) const override;
124 
125 private:
126  VirtualDimensionality(const Self&) = delete;
127  void operator=(const Self&) = delete;
128 
131 
132  unsigned int m_NumberOfPixels;
133  unsigned int m_NumberOfEndmembers;
134  double m_FAR;
135 };
136 }
137 
138 #ifndef OTB_MANUAL_INSTANTIATION
140 #endif
141 
142 #endif
Estimates the number of endmembers in a hyperspectral image with the Virtual Dimensionality (HFC) met...
VirtualDimensionality(const Self &)=delete
vnl_matrix< PrecisionType > MatrixType
void SetFAR(double falseAlarmRate)
itk::SmartPointer< const Self > ConstPointer
vnl_vector< PrecisionType > VectorType
void operator=(const Self &)=delete
void SetCorrelation(const MatrixType &m)
void SetNumberOfPixels(unsigned int n)
void SetCovariance(const MatrixType &m)
itk::SmartPointer< Self > Pointer
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.