21 #ifndef otbMultivariateAlterationDetectorImageFilter_hxx
22 #define otbMultivariateAlterationDetectorImageFilter_hxx
27 #include "vnl/algo/vnl_matrix_inverse.h"
28 #include "vnl/algo/vnl_generalized_eigensystem.h"
30 #include "itkImageRegionIterator.h"
31 #include "itkProgressReporter.h"
35 template <
class TInputImage,
class TOutputImage>
38 this->SetNumberOfRequiredInputs(2);
39 m_CovarianceEstimator = CovarianceEstimatorType::New();
40 this->DynamicMultiThreadingOn();
43 template <
class TInputImage,
class TOutputImage>
47 this->SetNthInput(0,
const_cast<TInputImage*
>(image1));
50 template <
class TInputImage,
class TOutputImage>
54 if (this->GetNumberOfInputs() < 1)
58 return static_cast<const TInputImage*
>(this->itk::ProcessObject::GetInput(0));
61 template <
class TInputImage,
class TOutputImage>
65 this->SetNthInput(1,
const_cast<TInputImage*
>(image2));
68 template <
class TInputImage,
class TOutputImage>
72 if (this->GetNumberOfInputs() < 2)
76 return static_cast<const TInputImage*
>(this->itk::ProcessObject::GetInput(1));
80 template <
class TInputImage,
class TOutputImage>
84 Superclass::GenerateOutputInformation();
87 const TInputImage* input1Ptr = this->GetInput1();
88 const TInputImage* input2Ptr = this->GetInput2();
89 TOutputImage* outputPtr =
const_cast<TOutputImage*
>(this->GetOutput());
92 unsigned int nbComp1 = input1Ptr->GetNumberOfComponentsPerPixel();
93 unsigned int nbComp2 = input2Ptr->GetNumberOfComponentsPerPixel();
94 unsigned int outNbComp = std::max(nbComp1, nbComp2);
96 outputPtr->SetNumberOfComponentsPerPixel(outNbComp);
98 if (input1Ptr->GetLargestPossibleRegion() != input2Ptr->GetLargestPossibleRegion())
100 itkExceptionMacro(<<
"Input images does not have the same size!");
105 concatenateFilter->SetInput1(input1Ptr);
106 concatenateFilter->SetInput2(input2Ptr);
110 m_CovarianceEstimator->SetInput(concatenateFilter->GetOutput());
111 m_CovarianceEstimator->Update();
112 m_CovarianceMatrix = m_CovarianceEstimator->GetCovariance();
113 m_MeanValues = m_CovarianceEstimator->GetMean();
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);
125 for (
unsigned int i = 0; i < nbComp1; ++i)
127 m_Mean1[i] = m_MeanValues[i];
130 for (
unsigned int i = 0; i < nbComp2; ++i)
132 m_Mean2[i] = m_MeanValues[nbComp1 + i];
135 if (nbComp1 == nbComp2)
139 VnlMatrixType invs22 = vnl_matrix_inverse<RealType>(s22).as_matrix();
144 vnl_generalized_eigensystem ges(s12s22is21, s11);
149 m_Rho = ges.D.get_diagonal();
150 m_Rho = m_Rho.apply(&std::sqrt);
159 invstderr1.set_diagonal(diag1);
167 for (
unsigned int i = 0; i < nbComp1; ++i)
169 aux5 = aux5 + aux4.get_row(i);
172 sign1.set_diagonal(aux5);
177 m_V2 = invs22 * s21 * m_V1;
182 aux2 = aux2.apply(&std::sqrt);
186 aux3.set_diagonal(aux2);
194 sl.update(s12, 0, nbComp1);
195 sl.update(s21, nbComp1, 0);
196 sr.update(s11, 0, 0);
197 sr.update(s22, nbComp1, nbComp1);
199 vnl_generalized_eigensystem ges(sl, sr);
205 m_V1 = V.extract(nbComp1, nbComp1);
206 m_V2 = V.extract(nbComp2, nbComp2, nbComp1, 0);
208 m_Rho = ges.D.get_diagonal().flip().extract(std::max(nbComp1, nbComp2), 0);
214 aux2 = aux2.apply(&std::sqrt);
218 aux3.set_diagonal(aux2);
225 invstderr1.set_diagonal(diag1);
233 for (
unsigned int i = 0; i < nbComp1; ++i)
235 aux5 = aux5 + aux4.get_row(i);
238 sign1.set_diagonal(aux5);
244 aux1 = m_V2.transpose() * (s22 * m_V2);
245 aux2 = aux1.get_diagonal();
246 aux2 = aux2.apply(&std::sqrt);
250 aux3.set_diagonal(aux2);
255 aux5 = (m_V1.transpose() * s12 * m_V2).transpose().get_diagonal();
256 sign2.set_diagonal(aux5);
264 template <
class TInputImage,
class TOutputImage>
268 const TInputImage* input1Ptr = this->GetInput1();
269 const TInputImage* input2Ptr = this->GetInput2();
270 TOutputImage* outputPtr = this->GetOutput();
273 typedef itk::ImageRegionConstIterator<InputImageType> ConstIteratorType;
274 typedef itk::ImageRegionIterator<OutputImageType> IteratorType;
276 IteratorType outIt(outputPtr, outputRegionForThread);
277 ConstIteratorType inIt1(input1Ptr, outputRegionForThread);
278 ConstIteratorType inIt2(input2Ptr, outputRegionForThread);
285 unsigned int nbComp1 = input1Ptr->GetNumberOfComponentsPerPixel();
286 unsigned int nbComp2 = input2Ptr->GetNumberOfComponentsPerPixel();
287 unsigned int outNbComp = outputPtr->GetNumberOfComponentsPerPixel();
290 while (!inIt1.IsAtEnd() && !inIt2.IsAtEnd() && !outIt.IsAtEnd())
298 for (
unsigned int i = 0; i < nbComp1; ++i)
300 x1[i] = inIt1.Get()[i];
303 for (
unsigned int i = 0; i < nbComp2; ++i)
305 x2[i] = inIt2.Get()[i];
311 for (
unsigned int i = 0; i < nbComp1; ++i)
316 for (
unsigned int i = 0; i < nbComp2; ++i)
318 x2bis[i] = second[i];
323 typename OutputImageType::PixelType outPixel(outNbComp);
325 if (nbComp1 == nbComp2)
327 for (
unsigned int i = 0; i < outNbComp; ++i)
329 outPixel[i] = mad[i];
334 for (
unsigned int i = 0; i < outNbComp; ++i)
336 outPixel[i] = mad[outNbComp - i - 1];
338 if (i < outNbComp - std::min(nbComp1, nbComp2))
340 outPixel[i] *= std::sqrt(2.);
vnl_vector< RealType > VnlVectorType
void SetInput2(const TInputImage *image2)
vnl_matrix< RealType > VnlMatrixType
MultivariateAlterationDetectorImageFilter()
void SetInput1(const TInputImage *image1)
TInputImage InputImageType
ConcatenateImageFilterType::Pointer ConcatenateImageFilterPointer
const TInputImage * GetInput2()
OutputImageType::RegionType OutputImageRegionType
void GenerateOutputInformation() override
void DynamicThreadedGenerateData(const OutputImageRegionType &outputRegionForThread) override
const TInputImage * GetInput1()
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.
T InverseValue(const T &value)
T SignOfValue(const T &value)