21 #ifndef otbMultiChannelsPolarimetricSynthesisFilter_hxx
22 #define otbMultiChannelsPolarimetricSynthesisFilter_hxx
27 #include "itkImageRegionIterator.h"
28 #include "itkProgressReporter.h"
37 template <
class TInputImage,
class TOutputImage,
class TFunction>
39 : m_PsiI(0.0), m_KhiI(0.0), m_PsiR(0.0), m_KhiR(0.0), m_Gain(1.0), m_Mode(0), m_EmissionH(false), m_EmissionV(false)
41 this->SetNumberOfRequiredInputs(1);
50 template <
class TInputImage,
class TOutputImage,
class TFunction>
57 typename Superclass::OutputImagePointer outputPtr = this->GetOutput();
58 typename Superclass::InputImageConstPointer inputPtr = this->GetInput();
60 if (!outputPtr || !inputPtr)
68 this->CallCopyInputRegionToOutputRegion(outputLargestPossibleRegion, inputPtr->GetLargestPossibleRegion());
69 outputPtr->SetLargestPossibleRegion(outputLargestPossibleRegion);
72 const itk::ImageBase<Superclass::InputImageDimension>* phyData;
74 phyData =
dynamic_cast<const itk::ImageBase<Superclass::InputImageDimension>*
>(this->GetInput());
82 const typename InputImageType::SpacingType& inputSpacing = inputPtr->GetSignedSpacing();
83 const typename InputImageType::PointType& inputOrigin = inputPtr->GetOrigin();
84 const typename InputImageType::DirectionType& inputDirection = inputPtr->GetDirection();
86 typename OutputImageType::SpacingType outputSpacing;
87 typename OutputImageType::PointType outputOrigin;
88 typename OutputImageType::DirectionType outputDirection;
92 for (i = 0; i < Superclass::InputImageDimension; ++i)
94 outputSpacing[i] = inputSpacing[i];
95 outputOrigin[i] = inputOrigin[i];
96 for (j = 0; j < Superclass::OutputImageDimension; ++j)
98 if (j < Superclass::InputImageDimension)
100 outputDirection[j][i] = inputDirection[j][i];
104 outputDirection[j][i] = 0.0;
108 for (; i < Superclass::OutputImageDimension; ++i)
110 outputSpacing[i] = 1.0;
111 outputOrigin[i] = 0.0;
112 for (j = 0; j < Superclass::OutputImageDimension; ++j)
116 outputDirection[j][i] = 1.0;
120 outputDirection[j][i] = 0.0;
126 outputPtr->SetSignedSpacing(outputSpacing);
127 outputPtr->SetOrigin(outputOrigin);
128 outputPtr->SetDirection(outputDirection);
133 itkExceptionMacro(<<
"otb::MultiChannelsPolarimetricSynthesisFilter::GenerateOutputInformation "
134 <<
"cannot cast input to " <<
typeid(itk::ImageBase<Superclass::InputImageDimension>*).name());
141 template <
class TInputImage,
class TOutputImage,
class TFunction>
143 itk::ThreadIdType threadId)
153 this->CallCopyOutputRegionToInputRegion(inputRegionForThread, outputRegionForThread);
156 itk::ImageRegionConstIterator<TInputImage> inputIt(inputPtr, inputRegionForThread);
157 itk::ImageRegionIterator<TOutputImage> outputIt(outputPtr, outputRegionForThread);
159 itk::ProgressReporter progress(
this, threadId, outputRegionForThread.GetNumberOfPixels());
162 outputIt.GoToBegin();
170 while (!inputIt.IsAtEnd())
172 outputIt.Set(m_Gain * GetFunctor()(inputIt.Get()[0], inputIt.Get()[1], inputIt.Get()[2], inputIt.Get()[3]));
175 progress.CompletedPixel();
181 while (!inputIt.IsAtEnd())
183 outputIt.Set(m_Gain * GetFunctor()(inputIt.Get()[0], inputIt.Get()[1], inputIt.Get()[1], inputIt.Get()[2]));
186 progress.CompletedPixel();
192 while (!inputIt.IsAtEnd())
194 outputIt.Set(m_Gain * GetFunctor()(inputIt.Get()[0], inputIt.Get()[1], 0, 0));
197 progress.CompletedPixel();
203 while (!inputIt.IsAtEnd())
205 outputIt.Set(m_Gain * GetFunctor()(0, 0, inputIt.Get()[2], inputIt.Get()[3]));
208 progress.CompletedPixel();
213 itkExceptionMacro(
"Unknown architecture : Polarimetric synthesis is impossible !");
221 template <
class TInputImage,
class TOutputImage,
class TFunction>
230 real = std::cos(DTOR * m_PsiI) * std::cos(DTOR * m_KhiI);
231 imag = -std::sin(DTOR * m_PsiI) * std::sin(DTOR * m_KhiI);
234 real = std::sin(DTOR * m_PsiI) * std::cos(DTOR * m_KhiI);
235 imag = std::cos(DTOR * m_PsiI) * std::sin(DTOR * m_KhiI);
238 real = std::cos(DTOR * m_PsiR) * std::cos(DTOR * m_KhiR);
239 imag = -std::sin(DTOR * m_PsiR) * std::sin(DTOR * m_KhiR);
242 real = std::sin(DTOR * m_PsiR) * std::cos(DTOR * m_KhiR);
243 imag = std::cos(DTOR * m_PsiR) * std::sin(DTOR * m_KhiR);
258 template <
class TInputImage,
class TOutputImage,
class TFunction>
290 itkExceptionMacro(
"Unknown architecture : Polarimetric synthesis is impossible !!");
296 else if (GetMode() == 2)
303 template <
class TInputImage,
class TOutputImage,
class TFunction>
307 int NumberOfImages = this->GetInput()->GetNumberOfComponentsPerPixel();
310 m_ArchitectureType->DetermineArchitecture(NumberOfImages, GetEmissionH(), GetEmissionV());
313 VerifyAndForceInputs();
316 ComputeElectromagneticFields();
322 template <
class TInputImage,
class TOutputImage,
class TFunction>
333 template <
class TInputImage,
class TOutputImage,
class TFunction>
336 SetPsiR(m_PsiI + 90);
345 template <
class TInputImage,
class TOutputImage,
class TFunction>
348 this->Superclass::PrintSelf(os, indent);
349 os << indent <<
"PsiI: " << m_PsiI << std::endl;
350 os << indent <<
"KhiI: " << m_KhiI << std::endl;
351 os << indent <<
"PsiR: " << m_PsiR << std::endl;
352 os << indent <<
"KhiR: " << m_KhiR << std::endl;
355 os << indent <<
"Ei0 im: " << m_Ei[0].imag() << std::endl;
356 os << indent <<
"Ei0 re: " << m_Ei[0].real() << std::endl;
357 os << indent <<
"Ei1 im: " << m_Ei[1].imag() << std::endl;
358 os << indent <<
"Ei1 re: " << m_Ei[1].real() << std::endl;
360 os << indent <<
"Er0 im: " << m_Er[0].imag() << std::endl;
361 os << indent <<
"Er0 re: " << m_Er[0].real() << std::endl;
362 os << indent <<
"Er1 im: " << m_Er[1].imag() << std::endl;
363 os << indent <<
"Er1 re: " << m_Er[1].real() << std::endl;