OTB  10.0.0
Orfeo Toolbox
otbReduceSpectralResponse.h
Go to the documentation of this file.
1 /*
2  * Copyright (C) 2005-2024 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 otbReduceSpectralResponse_h
22 #define otbReduceSpectralResponse_h
23 
24 
25 #include "itkDataObject.h"
26 #include <itkObjectFactory.h>
27 #include <vector>
28 #include <utility>
29 #include <limits>
30 
31 namespace otb
32 {
63 template <class TSpectralResponse, class TRSR>
64 class ReduceSpectralResponse : public itk::DataObject
65 {
66 
67 public:
70  typedef itk::DataObject Superclass;
71  typedef itk::SmartPointer<Self> Pointer;
72  typedef itk::SmartPointer<const Self> ConstPointer;
73 
75  typedef TSpectralResponse InputSpectralResponseType;
76  typedef TRSR InputRSRType;
77 
78  typedef typename InputSpectralResponseType::PairType PairType;
79  typedef typename InputRSRType::Pointer InputRSRPointerType;
80  typedef typename InputSpectralResponseType::Pointer InputSpectralResponsePointerType;
81 
83  typedef typename InputRSRType::ValuePrecisionType ValuePrecisionType;
84 
85  typedef typename InputSpectralResponseType::VectorPairType VectorPairType;
86 
87  // TODO need a specific class for the integration of stectral responses (now it is in the functor)
88  typedef typename InputRSRType::RSRVectorType RSRVectorType;
89  typedef typename std::vector<ValuePrecisionType> ReduceSpectralResponseVectorType;
91  itkNewMacro(Self);
92  ;
93  itkTypeMacro(ReduceSpectralResponse, DataObject);
94  ;
96 
97  itkGetConstObjectMacro(InputSatRSR, InputRSRType);
98  ;
99  itkSetObjectMacro(InputSatRSR, InputRSRType);
100  ;
101 
102  itkGetConstObjectMacro(InputSpectralResponse, InputSpectralResponseType);
103  ;
104  itkSetObjectMacro(InputSpectralResponse, InputSpectralResponseType);
105  ;
106 
109  ;
111 
112  itkSetMacro(ReflectanceMode, bool);
113  itkGetConstMacro(ReflectanceMode, bool);
114 
116  virtual bool Clear();
117 
119  void PrintSelf(std::ostream& os, itk::Indent indent) const override;
120 
125  inline ValuePrecisionType operator()(const unsigned int numBand);
126 
128  void CalculateResponse();
129 
131  void LoadInputsFromFiles(const std::string& spectralResponseFile, const std::string& RSRFile, const unsigned int nbRSRBands,
132  ValuePrecisionType coefNormSpectre = 1.0, ValuePrecisionType coefNormRSR = 1.0);
133 
134 protected:
137 
139  // ReduceSpectralResponse( const std::string & filename );
140 
143 
145  // void PrintSelf(std::ostream& os, itk::Indent indent) const;
146 
149 
152 
153 private:
154  ReduceSpectralResponse(const Self&) = delete;
155  void operator=(const Self&) = delete;
156 
159 
162 };
163 
164 } // end namespace otb
165 
166 
167 #ifndef OTB_MANUAL_INSTANTIATION
169 #endif
170 
171 #endif
This class computes the reduced spectral response of each band of a sensor.
InputSpectralResponsePointerType m_ReduceResponse
InputRSRType::PrecisionType PrecisionType
void LoadInputsFromFiles(const std::string &spectralResponseFile, const std::string &RSRFile, const unsigned int nbRSRBands, ValuePrecisionType coefNormSpectre=1.0, ValuePrecisionType coefNormRSR=1.0)
InputRSRType::ValuePrecisionType ValuePrecisionType
InputSpectralResponseType::Pointer InputSpectralResponsePointerType
itk::SmartPointer< const Self > ConstPointer
void PrintSelf(std::ostream &os, itk::Indent indent) const override
InputRSRType::Pointer InputRSRPointerType
InputSpectralResponsePointerType m_InputSpectralResponse
ValuePrecisionType operator()(const unsigned int numBand)
ReduceSpectralResponse(const Self &)=delete
InputRSRType::RSRVectorType RSRVectorType
InputSpectralResponseType::PairType PairType
std::vector< ValuePrecisionType > ReduceSpectralResponseVectorType
itk::SmartPointer< Self > Pointer
InputSpectralResponseType::VectorPairType VectorPairType
itkGetObjectMacro(ReduceResponse, InputSpectralResponseType)
void operator=(const Self &)=delete
double PrecisionType
Definition: otbDateUtils.h:29
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.