OTB  9.0.0
Orfeo Toolbox
otbSpectralResponse.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 otbSpectralResponse_hxx
22 #define otbSpectralResponse_hxx
23 
24 #include "itkNumericTraits.h"
25 
26 #include "otbSpectralResponse.h"
27 
28 #include <algorithm>
29 #include <fstream>
30 
31 namespace otb
32 {
33 
34 template <class TPrecision, class TValuePrecision>
36 {
37  m_SensitivityThreshold = 0.01;
38  m_IntervalComputed = false;
39  m_PosGuess = 0;
40  m_UsePosGuess = false;
41 }
42 
43 template <class TPrecision, class TValuePrecision>
44 void SpectralResponse<TPrecision, TValuePrecision>::Load(const std::string& filename, ValuePrecisionType coefNormalization)
45 {
46  // Parse JPL file spectral response (ASCII file)
47  // Begin line 27
48  std::ifstream fin(filename);
49  if (fin.fail())
50  {
51  itkExceptionMacro(<< "Error opening file" << filename);
52  }
53 
54  int NumLigne = 26; // Go to the line 27
55  // Ignore first 26th lines which are metadatas information
56  for (int i = 0; i < NumLigne; ++i)
57  fin.ignore(std::numeric_limits<std::streamsize>::max(), '\n');
58 
59  while (!fin.eof())
60  {
61  // For each
62  std::pair<TPrecision, TValuePrecision> currentPair;
63 
64  fin >> currentPair.first;
65  fin >> currentPair.second;
66  currentPair.second = currentPair.second / coefNormalization;
67  if (currentPair.first != itk::NumericTraits<TPrecision>::ZeroValue() && currentPair.second != itk::NumericTraits<TValuePrecision>::ZeroValue())
68  // Add not null pair of values to the vector
69  m_Response.push_back(currentPair);
70  }
71  fin.close();
72  // Sort the vector using the specific functor sort_pair
73  std::sort(m_Response.begin(), m_Response.end(), sort_pair());
74 
75  m_IntervalComputed = false;
76 }
77 
78 template <class TPrecision, class TValuePrecision>
80 {
81  m_Response.clear();
82  m_IntervalComputed = false;
83  return true;
84 }
85 
86 template <class TPrecision, class TValuePrecision>
88 {
89  return m_Response.size();
90 }
91 
92 
93 template <class TPrecision, class TValuePrecision>
95 {
96  m_PosGuess = 0;
97  if (m_Response.size() <= 1)
98  {
99  itkExceptionMacro(<< "ERROR spectral response need at least 2 value to perform interpolation.");
100  }
101 
102  TPrecision lambdaMax = this->GetInterval().second;
103  if (lambda > lambdaMax)
104  return;
105  typename VectorPairType::const_iterator it = m_Response.begin();
106 
107  while (((*it).first < lambda))
108  {
109  m_PosGuess++;
110  ++it;
111  if (it == (m_Response.end()))
112  return;
113  }
114 
115  if (m_PosGuess > 0)
116  m_PosGuess--;
117  return;
118 }
119 
120 
121 template <class TPrecision, class TValuePrecision>
123 operator()(const PrecisionType& lambda)
124 {
125 
126  // Suppose that the vector is sorted
127 
128  // Guess a starting lambda
129  if (m_Response.size() <= 1)
130  {
131  itkExceptionMacro(<< "ERROR spectral response need at least 2 value to perform interpolation.");
132  }
133 
134  typename VectorPairType::const_iterator beg = m_Response.begin();
135  typename VectorPairType::const_iterator last = m_Response.end();
136  --last;
137 
138  TPrecision lambdaMin = this->GetInterval().first;
139  TPrecision lambdaMax = this->GetInterval().second;
140 
141  if (lambda < lambdaMin)
142  return static_cast<TValuePrecision>(0.0);
143  if (lambda > lambdaMax)
144  return static_cast<TValuePrecision>(0.0);
145 
146  typename VectorPairType::const_iterator it;
147 
148  if (m_UsePosGuess)
149  it = beg + m_PosGuess;
150  else
151  it = beg;
152 
153  TPrecision lambda1 = (*beg).first;
154  TValuePrecision SR1 = static_cast<TValuePrecision>(0.0);
155  while (((*it).first < lambda))
156  {
157 
158  lambda1 = (*it).first;
159  SR1 = (*it).second;
160  ++it;
161  if (it == (m_Response.end()))
162  {
163  return static_cast<TValuePrecision>(0.0);
164  }
165  }
166  TPrecision lambda2 = (*it).first;
167 
168  // if the guess is just right
169  if (lambda2 == lambda)
170  {
171  return (*it).second;
172  }
173  else
174  {
175 
176  TPrecision lambdaDist = lambda - lambda1;
177  TPrecision ratio = lambdaDist / (lambda2 - lambda1);
178 
179  TValuePrecision SR2 = (*it).second;
180 
181  return static_cast<TValuePrecision>(ratio * SR1 + (1 - ratio) * SR2);
182  }
183 }
184 
185 template <class TPrecision, class TValuePrecision>
187 {
188  typename ImageType::IndexType start;
189  start[0] = 0;
190  start[1] = 0;
191 
192  typename ImageType::SizeType size;
193  size[0] = 1;
194  size[1] = 1;
195 
196  typename ImageType::PointType origin;
197  origin[0] = 0;
198  origin[1] = 0;
199 
200  typename ImageType::SpacingType spacing;
201  spacing[0] = 1;
202  spacing[1] = 1;
203 
204  typename ImageType::RegionType region;
205  region.SetSize(size);
206  region.SetIndex(start);
207 
208  image->SetRegions(region);
209  image->SetNumberOfComponentsPerPixel(this->Size());
210  image->Allocate();
211 
212  typename ImageType::IndexType idx;
213  typename ImageType::PixelType pixel;
214  pixel.SetSize(this->Size());
215 
216  for (unsigned int j = 0; j < this->Size(); ++j)
217  {
218  pixel[j] = m_Response[j].second;
219  }
220  idx[0] = 0;
221  idx[1] = 0;
222  image->SetPixel(idx, pixel);
223  return image;
224 }
225 
226 template <class TPrecision, class TValuePrecision>
228 {
229 
230  typename ImageType::IndexType idx;
231  idx[0] = 0;
232  idx[1] = 0;
233 
234  for (unsigned int j = 0; j < this->Size(); ++j)
235  {
236  m_Response[j].second = image->GetPixel(idx)[j];
237  }
238  m_IntervalComputed = false;
239 }
240 
241 template <class TPrecision, class TValuePrecision>
244 {
245 
246  // Assume that the SR is sorted
247  typename FilterFunctionValuesType::ValuesVectorType valuesVector;
248  Self& responseCalculator = *this;
249  for (double i = m_Response.front().first; i <= m_Response.back().first; i += step)
250  {
251  valuesVector.push_back(responseCalculator(i));
252  }
253  FilterFunctionValuesPointerType functionValues = FilterFunctionValuesType::New();
254 
255  functionValues->SetFilterFunctionValues(valuesVector);
256  functionValues->SetMinSpectralValue(m_Response.front().first);
257  functionValues->SetMaxSpectralValue(m_Response.back().first);
258  functionValues->SetUserStep(step);
259 
260  return functionValues;
261 }
262 
263 template <class TPrecision, class TValuePrecision>
265 {
266  typename VectorPairType::const_iterator it = m_Response.begin();
267 
268  while ((*it).second <= m_SensitivityThreshold)
269  {
270  ++it;
271  if (it == (m_Response.end() - 1))
272  {
273  m_Interval.first = static_cast<TPrecision>(0.0);
274  m_Interval.second = static_cast<TPrecision>(0.0);
275  m_IntervalComputed = true;
276  return;
277  }
278  }
279  m_Interval.first = (*it).first;
280  it = (m_Response.end() - 1);
281  while ((*it).second <= m_SensitivityThreshold)
282  {
283  if (it == (m_Response.begin()))
284  {
285  m_Interval.second = (*it).first;
286  m_IntervalComputed = true;
287  return;
288  }
289  --it;
290  }
291 
292  m_Interval.second = (*it).first;
293  m_IntervalComputed = true;
294 }
295 
296 template <class TPrecision, class TValuePrecision>
297 void SpectralResponse<TPrecision, TValuePrecision>::PrintSelf(std::ostream& os, itk::Indent indent) const
298 {
299  Superclass::PrintSelf(os, indent);
300  os << std::endl;
301  os << indent << "[Wavelength (micrometers), Reflectance (percent)]" << std::endl;
302  for (typename VectorPairType::const_iterator it = m_Response.begin(); it != m_Response.end(); ++it)
303  {
304  os << indent << "Num " << it - m_Response.begin() << ": [" << (*it).first << "," << (*it).second << "]" << std::endl;
305  }
306 }
307 
308 } // end namespace otb
309 
310 #endif
otb::VectorImage::RegionType
Superclass::RegionType RegionType
Definition: otbVectorImage.h:113
otb::SpectralResponse::SetFromImage
void SetFromImage(ImagePointerType image)
Definition: otbSpectralResponse.hxx:227
otb::SpectralResponse::Size
virtual unsigned int Size() const
Definition: otbSpectralResponse.hxx:87
otb::SpectralResponse::PrintSelf
void PrintSelf(std::ostream &os, itk::Indent indent) const override
Definition: otbSpectralResponse.hxx:297
otb::SpectralResponse::FilterFunctionValuesPointerType
itk::SmartPointer< FilterFunctionValuesType > FilterFunctionValuesPointerType
Definition: otbSpectralResponse.h:79
otb::SpectralResponse::ValuePrecisionType
TValuePrecision ValuePrecisionType
Definition: otbSpectralResponse.h:66
otb::VectorImage::SpacingType
Superclass::SpacingType SpacingType
Definition: otbVectorImage.h:117
otb::SpectralResponse::GetFilterFunctionValues
FilterFunctionValuesPointerType GetFilterFunctionValues(double step=0.0025)
Definition: otbSpectralResponse.hxx:243
otb::VectorImage::SizeType
Superclass::SizeType SizeType
Definition: otbVectorImage.h:107
otb::SpectralResponse::GetImage
ImagePointerType GetImage(ImagePointerType image) const
Definition: otbSpectralResponse.hxx:186
otb::SpectralResponse::ImagePointerType
ImageType::Pointer ImagePointerType
Definition: otbSpectralResponse.h:74
otb::SpectralResponse::PrecisionType
TPrecision PrecisionType
Definition: otbSpectralResponse.h:65
otb
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.
Definition: otbJoinContainer.h:32
otb::SpectralResponse::Load
void Load(const std::string &filename, ValuePrecisionType coefNormalization=1.0)
Definition: otbSpectralResponse.hxx:44
otb::VectorImage::PointType
Superclass::PointType PointType
Definition: otbVectorImage.h:121
otb::SpectralResponse::SpectralResponse
SpectralResponse()
Definition: otbSpectralResponse.hxx:35
otb::SpectralResponse::Clear
virtual bool Clear()
Definition: otbSpectralResponse.hxx:79
otb::SpectralResponse::operator()
ValuePrecisionType operator()(const PrecisionType &lambda)
Definition: otbSpectralResponse.hxx:123
otbSpectralResponse.h
otb::VectorImage::PixelType
Superclass::PixelType PixelType
Definition: otbVectorImage.h:63
otb::VectorImage::IndexType
Superclass::IndexType IndexType
Definition: otbVectorImage.h:101
otb::SpectralResponse::ComputeInterval
void ComputeInterval()
Definition: otbSpectralResponse.hxx:264
otb::SpectralResponse
This class represents the spectral response of an object (or a satellite band).
Definition: otbSpectralResponse.h:55
otb::SpectralResponse::SetPosGuessMin
void SetPosGuessMin(const PrecisionType &lambda)
Definition: otbSpectralResponse.hxx:94
otb::SpectralResponse::sort_pair
Definition: otbSpectralResponse.h:128
otb::FilterFunctionValues::ValuesVectorType
std::vector< WavelengthSpectralBandType > ValuesVectorType
Definition: otbFilterFunctionValues.h:60