OTB  10.0.0
Orfeo Toolbox
otbMultiChannelsPolarimetricSynthesisFilter.h
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 otbMultiChannelsPolarimetricSynthesisFilter_h
22 #define otbMultiChannelsPolarimetricSynthesisFilter_h
23 
24 #include "itkInPlaceImageFilter.h"
26 #include "otbPolarimetricData.h"
27 #include <complex>
28 
29 namespace otb
30 {
31 
46 template <class TInputImage, class TOutputImage,
47  class TFunction = Functor::PolarimetricSynthesisFunctor<typename TInputImage::InternalPixelType, typename TInputImage::InternalPixelType,
48  typename TInputImage::InternalPixelType, typename TInputImage::InternalPixelType,
49  typename TOutputImage::PixelType>>
50 class ITK_EXPORT MultiChannelsPolarimetricSynthesisFilter : public itk::InPlaceImageFilter<TInputImage, TOutputImage>
51 {
52 public:
55  typedef itk::InPlaceImageFilter<TInputImage, TOutputImage> Superclass;
56  typedef itk::SmartPointer<Self> Pointer;
57  typedef itk::SmartPointer<const Self> ConstPointer;
58 
60  itkNewMacro(Self);
61 
63  itkTypeMacro(MultiChannelsPolarimetricSynthesisFilter, InPlaceImageFilter);
64 
66  typedef std::complex<double> InputPixelType;
67  typedef TFunction FunctorType;
68  typedef TInputImage InputImageType;
69  typedef typename InputImageType::ConstPointer InputImagePointer;
70  typedef typename InputImageType::RegionType InputImageRegionType;
71  typedef typename InputImageType::PixelType InputImagePixelType;
72  typedef TOutputImage OutputImageType;
73  typedef typename OutputImageType::Pointer OutputImagePointer;
74  typedef typename OutputImageType::RegionType OutputImageRegionType;
75  typedef typename OutputImageType::PixelType OutputImagePixelType;
76  typedef typename std::complex<double> ComplexType;
77  typedef typename itk::FixedArray<ComplexType, 2> ComplexArrayType;
78  typedef typename itk::FixedArray<int, 4> IndexArrayType;
79 
85  {
86  return m_Functor;
87  }
88  const FunctorType& GetFunctor() const
89  {
90  return m_Functor;
91  }
93 
100  void SetFunctor(const FunctorType& functor)
101  {
102  if (m_Functor != functor)
103  {
104  m_Functor = functor;
105  this->Modified();
106  }
107  }
108 
111  {
112  m_Ei = ei;
113  this->GetFunctor().SetEi(ei);
114  this->Modified();
115  }
116 
119  {
120  m_Er = er;
121  this->GetFunctor().SetEr(er);
122  this->Modified();
123  }
125 
127  itkSetMacro(PsiI, double);
128  itkGetMacro(PsiI, double);
129 
131  itkSetMacro(KhiI, double);
132  itkGetMacro(KhiI, double);
133 
135  itkSetMacro(PsiR, double);
136  itkGetMacro(PsiR, double);
137 
139  itkSetMacro(KhiR, double);
140  itkGetMacro(KhiR, double);
141 
143  itkSetMacro(EmissionH, bool);
144  itkGetMacro(EmissionH, bool);
145 
147  itkSetMacro(EmissionV, bool);
148  itkGetMacro(EmissionV, bool);
149 
151  itkSetMacro(Mode, int);
152  itkGetMacro(Mode, int);
153 
155  itkSetMacro(Gain, double);
156  itkGetMacro(Gain, double);
157 
159  void ForceCoPolar();
160 
162  void ForceCrossPolar();
163 
164 protected:
167 
170  {
171  }
172 
183  void GenerateOutputInformation() override;
184 
185  void BeforeThreadedGenerateData() override;
186 
197  void DynamicThreadedGenerateData(const OutputImageRegionType& outputRegionForThread) override;
198 
200  void ComputeElectromagneticFields();
201 
203  void VerifyAndForceInputs();
204 
205  void PrintSelf(std::ostream& os, itk::Indent indent) const override;
206 
207 private:
209 
211  double m_PsiI;
212 
214  double m_KhiI;
215 
217  double m_PsiR;
218 
220  double m_KhiR;
221 
223  double m_Gain;
224 
226  int m_Mode;
227 
230 
233 
236 
238 
242 };
243 
244 } // end namespace otb
245 
246 #ifndef OTB_MANUAL_INSTANTIATION
248 #endif
249 
250 #endif
This class computes the polarimetric synthesis from two to four radar images, depending on the polari...
itk::InPlaceImageFilter< TInputImage, TOutputImage > Superclass
MultiChannelsPolarimetricSynthesisFilter(const Self &)=delete
itk::SmartPointer< Self > Pointer
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.