OTB  9.0.0
Orfeo Toolbox
otbSpectralAngleDistanceImageFilter.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 otbSpectralAngleDistanceImageFilter_hxx
22 #define otbSpectralAngleDistanceImageFilter_hxx
23 
25 #include "itkImageRegionIterator.h"
26 #include "itkProgressReporter.h"
27 #include "otbMacro.h"
28 #include "otbMath.h"
29 
30 namespace otb
31 {
35 template <class TInputImage, class TOutputImage>
37 {
38  m_ReferencePixel = 0;
39 }
40 
41 template <class TInputImage, class TOutputImage>
43 {
44  if (this->GetInput()->GetNumberOfComponentsPerPixel() == 1)
45  {
46  itkExceptionMacro(<< "Not valid input image : mono channel image gives a nul output image.");
47  }
48 }
49 
50 template <class TInputImage, class TOutputImage>
51 void SpectralAngleDistanceImageFilter<TInputImage, TOutputImage>::ThreadedGenerateData(const OutputImageRegionType& outputRegionForThread,
52  itk::ThreadIdType threadId)
53 {
54 
55  if (m_ReferencePixel.Size() == 0)
56  {
57  itkExceptionMacro(<< "Reference pixel is not set!");
58  }
59 
60  InputImageConstPointerType inputPtr = this->GetInput();
61  OutputImagePointerType outputPtr = this->GetOutput();
62 
63  // inputPtr->UpdateOutputInformation();
64  // Check if the reference pixel size matches the input image number of components.
65  if (m_ReferencePixel.GetSize() != inputPtr->GetNumberOfComponentsPerPixel())
66  {
67  itkExceptionMacro(<< "Reference pixel size (" << m_ReferencePixel.GetSize() << " and input image pixel size (" << inputPtr->GetNumberOfComponentsPerPixel()
68  << ") don't match!");
69  }
70 
71  // Define the portion of the input to walk for this thread, using
72  // the CallCopyOutputRegionToInputRegion method allows for the input
73  // and output images to be different dimensions
74  InputImageRegionType inputRegionForThread;
75  this->CallCopyOutputRegionToInputRegion(inputRegionForThread, outputRegionForThread);
76 
77  // Define the iterators
78  itk::ImageRegionConstIterator<InputImageType> inputIt(inputPtr, inputRegionForThread);
79  itk::ImageRegionIterator<OutputImageType> outputIt(outputPtr, outputRegionForThread);
80  itk::ProgressReporter progress(this, threadId, outputRegionForThread.GetNumberOfPixels());
81 
82  inputIt.GoToBegin();
83  outputIt.GoToBegin();
84 
85  while (!inputIt.IsAtEnd() && !outputIt.IsAtEnd())
86  {
87  double dist = 0.0;
88  double scalarProd = 0.0;
89  double normProd = 0.0;
90  double normProd1 = 0.0;
91  double normProd2 = 0.0;
92  InputPixelType pixel = inputIt.Get();
93  for (unsigned int i = 0; i < pixel.Size(); ++i)
94  {
95  scalarProd += pixel[i] * m_ReferencePixel[i];
96  normProd1 += pixel[i] * pixel[i];
97  normProd2 += m_ReferencePixel[i] * m_ReferencePixel[i];
98  }
99  normProd = normProd1 * normProd2;
100 
101  if (normProd == 0.0)
102  {
103  dist = 0.0;
104  }
105  else
106  {
107  dist = std::acos(scalarProd / std::sqrt(normProd));
108  }
109  //------ This part was suppressed since the filter must perform only the spectral angle computation ---
110  // Spectral angle normalization
111  // dist = dist/(CONST_PI_2);
112  // square ponderation
113  // dist = std::sqrt(dist);
114  outputIt.Set(static_cast<OutputPixelType>(dist));
115  ++inputIt;
116  ++outputIt;
117  progress.CompletedPixel(); // potential exception thrown here
118  }
119 }
otbMath.h
otb
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.
Definition: otbJoinContainer.h:32
otbMacro.h
otbSpectralAngleDistanceImageFilter.h
otb::SpectralAngleDistanceImageFilter::ThreadedGenerateData
void ThreadedGenerateData(const OutputImageRegionType &outputRegionForThread, itk::ThreadIdType threadId) override
otb::SpectralAngleDistanceImageFilter::SpectralAngleDistanceImageFilter
SpectralAngleDistanceImageFilter()
otb::SpectralAngleDistanceImageFilter::BeforeThreadedGenerateData
void BeforeThreadedGenerateData() override