OTB  10.0.0
Orfeo Toolbox
otbMultivariateAlterationDetectorImageFilter.hxx
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 otbMultivariateAlterationDetectorImageFilter_hxx
22 #define otbMultivariateAlterationDetectorImageFilter_hxx
23 
25 #include "otbMath.h"
26 
27 #include "vnl/algo/vnl_matrix_inverse.h"
28 #include "vnl/algo/vnl_generalized_eigensystem.h"
29 
30 #include "itkImageRegionIterator.h"
31 #include "itkProgressReporter.h"
32 
33 namespace otb
34 {
35 template <class TInputImage, class TOutputImage>
37 {
38  this->SetNumberOfRequiredInputs(2);
39  m_CovarianceEstimator = CovarianceEstimatorType::New();
40  this->DynamicMultiThreadingOn();
41 }
42 
43 template <class TInputImage, class TOutputImage>
45 {
46  // Process object is not const-correct so the const casting is required.
47  this->SetNthInput(0, const_cast<TInputImage*>(image1));
48 }
49 
50 template <class TInputImage, class TOutputImage>
53 {
54  if (this->GetNumberOfInputs() < 1)
55  {
56  return nullptr;
57  }
58  return static_cast<const TInputImage*>(this->itk::ProcessObject::GetInput(0));
59 }
60 
61 template <class TInputImage, class TOutputImage>
63 {
64  // Process object is not const-correct so the const casting is required.
65  this->SetNthInput(1, const_cast<TInputImage*>(image2));
66 }
67 
68 template <class TInputImage, class TOutputImage>
71 {
72  if (this->GetNumberOfInputs() < 2)
73  {
74  return nullptr;
75  }
76  return static_cast<const TInputImage*>(this->itk::ProcessObject::GetInput(1));
77 }
78 
79 
80 template <class TInputImage, class TOutputImage>
82 {
83  // Call superclass implementation
84  Superclass::GenerateOutputInformation();
85 
86  // Retrieve input images pointers
87  const TInputImage* input1Ptr = this->GetInput1();
88  const TInputImage* input2Ptr = this->GetInput2();
89  TOutputImage* outputPtr = const_cast<TOutputImage*>(this->GetOutput());
90 
91  // Get the number of components for each image
92  unsigned int nbComp1 = input1Ptr->GetNumberOfComponentsPerPixel();
93  unsigned int nbComp2 = input2Ptr->GetNumberOfComponentsPerPixel();
94  unsigned int outNbComp = std::max(nbComp1, nbComp2);
95 
96  outputPtr->SetNumberOfComponentsPerPixel(outNbComp);
97 
98  if (input1Ptr->GetLargestPossibleRegion() != input2Ptr->GetLargestPossibleRegion())
99  {
100  itkExceptionMacro(<< "Input images does not have the same size!");
101  }
102 
103  // First concatenate both images
104  ConcatenateImageFilterPointer concatenateFilter = ConcatenateImageFilterType::New();
105  concatenateFilter->SetInput1(input1Ptr);
106  concatenateFilter->SetInput2(input2Ptr);
107 
108  // The compute covariance matrix
109  // TODO: implement switch off of this computation
110  m_CovarianceEstimator->SetInput(concatenateFilter->GetOutput());
111  m_CovarianceEstimator->Update();
112  m_CovarianceMatrix = m_CovarianceEstimator->GetCovariance();
113  m_MeanValues = m_CovarianceEstimator->GetMean();
114 
115  // Extract sub-matrices of the covariance matrix
116  VnlMatrixType s11 = m_CovarianceMatrix.GetVnlMatrix().extract(nbComp1, nbComp1);
117  VnlMatrixType s22 = m_CovarianceMatrix.GetVnlMatrix().extract(nbComp2, nbComp2, nbComp1, nbComp1);
118  VnlMatrixType s12 = m_CovarianceMatrix.GetVnlMatrix().extract(nbComp1, nbComp2, 0, nbComp1);
119  VnlMatrixType s21 = s12.transpose();
120 
121  // Extract means
122  m_Mean1 = VnlVectorType(nbComp1, 0);
123  m_Mean2 = VnlVectorType(nbComp2, 0);
124 
125  for (unsigned int i = 0; i < nbComp1; ++i)
126  {
127  m_Mean1[i] = m_MeanValues[i];
128  }
129 
130  for (unsigned int i = 0; i < nbComp2; ++i)
131  {
132  m_Mean2[i] = m_MeanValues[nbComp1 + i];
133  }
134 
135  if (nbComp1 == nbComp2)
136  {
137  // Case where nbbands1 == nbbands2
138 
139  VnlMatrixType invs22 = vnl_matrix_inverse<RealType>(s22).as_matrix();
140 
141  // Build the generalized eigensystem
142  VnlMatrixType s12s22is21 = s12 * invs22 * s21;
143 
144  vnl_generalized_eigensystem ges(s12s22is21, s11);
145 
146  m_V1 = ges.V;
147 
148  // Compute canonical correlation matrix
149  m_Rho = ges.D.get_diagonal();
150  m_Rho = m_Rho.apply(&std::sqrt);
151 
152  // We do not need to scale v1 since the
153  // vnl_generalized_eigensystem already gives unit variance
154 
155  VnlMatrixType invstderr1 = s11.apply(&std::sqrt);
156  invstderr1 = invstderr1.apply(&InverseValue);
157  VnlVectorType diag1 = invstderr1.get_diagonal();
158  invstderr1.fill(0);
159  invstderr1.set_diagonal(diag1);
160 
161  VnlMatrixType sign1 = VnlMatrixType(nbComp1, nbComp1, 0);
162 
163  VnlMatrixType aux4 = invstderr1 * s11 * m_V1;
164 
165  VnlVectorType aux5 = VnlVectorType(nbComp1, 0);
166 
167  for (unsigned int i = 0; i < nbComp1; ++i)
168  {
169  aux5 = aux5 + aux4.get_row(i);
170  }
171 
172  sign1.set_diagonal(aux5);
173  sign1 = sign1.apply(&SignOfValue);
174 
175  m_V1 = m_V1 * sign1;
176 
177  m_V2 = invs22 * s21 * m_V1;
178 
179  // Scale v2 for unit variance
180  VnlMatrixType aux1 = m_V2.transpose() * (s22 * m_V2);
181  VnlVectorType aux2 = aux1.get_diagonal();
182  aux2 = aux2.apply(&std::sqrt);
183  aux2 = aux2.apply(&InverseValue);
184  VnlMatrixType aux3 = VnlMatrixType(aux2.size(), aux2.size(), 0);
185  aux3.fill(0);
186  aux3.set_diagonal(aux2);
187  m_V2 = m_V2 * aux3;
188  }
189  else
190  {
191  VnlMatrixType sl(nbComp1 + nbComp2, nbComp1 + nbComp2, 0);
192  VnlMatrixType sr(nbComp1 + nbComp2, nbComp1 + nbComp2, 0);
193 
194  sl.update(s12, 0, nbComp1);
195  sl.update(s21, nbComp1, 0);
196  sr.update(s11, 0, 0);
197  sr.update(s22, nbComp1, nbComp1);
198 
199  vnl_generalized_eigensystem ges(sl, sr);
200 
201  VnlMatrixType V = ges.V;
202 
203  V.fliplr();
204 
205  m_V1 = V.extract(nbComp1, nbComp1);
206  m_V2 = V.extract(nbComp2, nbComp2, nbComp1, 0);
207 
208  m_Rho = ges.D.get_diagonal().flip().extract(std::max(nbComp1, nbComp2), 0);
209 
210  // Scale v1 to get a unit variance
211  VnlMatrixType aux1 = m_V1.transpose() * (s11 * m_V1);
212 
213  VnlVectorType aux2 = aux1.get_diagonal();
214  aux2 = aux2.apply(&std::sqrt);
215  aux2 = aux2.apply(&InverseValue);
216 
217  VnlMatrixType aux3 = VnlMatrixType(aux2.size(), aux2.size(), 0);
218  aux3.set_diagonal(aux2);
219  m_V1 = m_V1 * aux3;
220 
221  VnlMatrixType invstderr1 = s11.apply(&std::sqrt);
222  invstderr1 = invstderr1.apply(&InverseValue);
223  VnlVectorType diag1 = invstderr1.get_diagonal();
224  invstderr1.fill(0);
225  invstderr1.set_diagonal(diag1);
226 
227  VnlMatrixType sign1 = VnlMatrixType(nbComp1, nbComp1, 0);
228 
229  VnlMatrixType aux4 = invstderr1 * s11 * m_V1;
230 
231  VnlVectorType aux5 = VnlVectorType(nbComp1, 0);
232 
233  for (unsigned int i = 0; i < nbComp1; ++i)
234  {
235  aux5 = aux5 + aux4.get_row(i);
236  }
237 
238  sign1.set_diagonal(aux5);
239  sign1 = sign1.apply(&SignOfValue);
240 
241  m_V1 = m_V1 * sign1;
242 
243  // Scale v2 for unit variance
244  aux1 = m_V2.transpose() * (s22 * m_V2);
245  aux2 = aux1.get_diagonal();
246  aux2 = aux2.apply(&std::sqrt);
247  aux2 = aux2.apply(&InverseValue);
248  aux3 = VnlMatrixType(aux2.size(), aux2.size(), 0);
249  aux3.fill(0);
250  aux3.set_diagonal(aux2);
251  m_V2 = m_V2 * aux3;
252 
253  VnlMatrixType sign2 = VnlMatrixType(nbComp2, nbComp2, 0);
254 
255  aux5 = (m_V1.transpose() * s12 * m_V2).transpose().get_diagonal();
256  sign2.set_diagonal(aux5);
257  sign2 = sign2.apply(&SignOfValue);
258  m_V2 = m_V2 * sign2;
259 
260  m_Rho.flip();
261  }
262 }
263 
264 template <class TInputImage, class TOutputImage>
266 {
267  // Retrieve input images pointers
268  const TInputImage* input1Ptr = this->GetInput1();
269  const TInputImage* input2Ptr = this->GetInput2();
270  TOutputImage* outputPtr = this->GetOutput();
271 
272 
273  typedef itk::ImageRegionConstIterator<InputImageType> ConstIteratorType;
274  typedef itk::ImageRegionIterator<OutputImageType> IteratorType;
275 
276  IteratorType outIt(outputPtr, outputRegionForThread);
277  ConstIteratorType inIt1(input1Ptr, outputRegionForThread);
278  ConstIteratorType inIt2(input2Ptr, outputRegionForThread);
279 
280  inIt1.GoToBegin();
281  inIt2.GoToBegin();
282  outIt.GoToBegin();
283 
284  // Get the number of components for each image
285  unsigned int nbComp1 = input1Ptr->GetNumberOfComponentsPerPixel();
286  unsigned int nbComp2 = input2Ptr->GetNumberOfComponentsPerPixel();
287  unsigned int outNbComp = outputPtr->GetNumberOfComponentsPerPixel();
288 
289 
290  while (!inIt1.IsAtEnd() && !inIt2.IsAtEnd() && !outIt.IsAtEnd())
291  {
292  VnlVectorType x1(nbComp1, 0);
293  VnlVectorType x1bis(outNbComp, 0);
294  VnlVectorType x2(nbComp2, 0);
295  VnlVectorType x2bis(outNbComp, 0);
296  VnlVectorType mad(outNbComp, 0);
297 
298  for (unsigned int i = 0; i < nbComp1; ++i)
299  {
300  x1[i] = inIt1.Get()[i];
301  }
302 
303  for (unsigned int i = 0; i < nbComp2; ++i)
304  {
305  x2[i] = inIt2.Get()[i];
306  }
307 
308  VnlVectorType first = (x1 - m_Mean1) * m_V1;
309  VnlVectorType second = (x2 - m_Mean2) * m_V2;
310 
311  for (unsigned int i = 0; i < nbComp1; ++i)
312  {
313  x1bis[i] = first[i];
314  }
315 
316  for (unsigned int i = 0; i < nbComp2; ++i)
317  {
318  x2bis[i] = second[i];
319  }
320 
321  mad = x1bis - x2bis;
322 
323  typename OutputImageType::PixelType outPixel(outNbComp);
324 
325  if (nbComp1 == nbComp2)
326  {
327  for (unsigned int i = 0; i < outNbComp; ++i)
328  {
329  outPixel[i] = mad[i];
330  }
331  }
332  else
333  {
334  for (unsigned int i = 0; i < outNbComp; ++i)
335  {
336  outPixel[i] = mad[outNbComp - i - 1];
337 
338  if (i < outNbComp - std::min(nbComp1, nbComp2))
339  {
340  outPixel[i] *= std::sqrt(2.);
341  }
342  }
343  }
344 
345  outIt.Set(outPixel);
346 
347  ++inIt1;
348  ++inIt2;
349  ++outIt;
350  }
351 }
352 }
353 
354 #endif
void DynamicThreadedGenerateData(const OutputImageRegionType &outputRegionForThread) override
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.
T InverseValue(const T &value)
Definition: otbMath.h:96
T SignOfValue(const T &value)
Definition: otbMath.h:102