OTB  6.7.0
Orfeo Toolbox
otbMultiChannelsPolarimetricSynthesisFilter.h
Go to the documentation of this file.
1 /*
2  * Copyright (C) 2005-2019 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<
48  typename TInputImage::InternalPixelType,
49  typename TInputImage::InternalPixelType,
50  typename TInputImage::InternalPixelType,
51  typename TInputImage::InternalPixelType,
52  typename TOutputImage::PixelType> >
53 class ITK_EXPORT MultiChannelsPolarimetricSynthesisFilter : public itk::InPlaceImageFilter<TInputImage, TOutputImage>
54 {
55 public:
61 
63  itkNewMacro(Self);
64 
66  itkTypeMacro(MultiChannelsPolarimetricSynthesisFilter, InPlaceImageFilter);
67 
69  typedef std::complex <double> InputPixelType;
70  typedef TFunction FunctorType;
71  typedef TInputImage InputImageType;
72  typedef typename InputImageType::ConstPointer InputImagePointer;
73  typedef typename InputImageType::RegionType InputImageRegionType;
74  typedef typename InputImageType::PixelType InputImagePixelType;
75  typedef TOutputImage OutputImageType;
76  typedef typename OutputImageType::Pointer OutputImagePointer;
77  typedef typename OutputImageType::RegionType OutputImageRegionType;
78  typedef typename OutputImageType::PixelType OutputImagePixelType;
79  typedef typename std::complex <double> ComplexType;
82 
88  {
89  return m_Functor;
90  }
91  const FunctorType& GetFunctor() const
92  {
93  return m_Functor;
94  }
96 
103  void SetFunctor(const FunctorType& functor)
104  {
105  if (m_Functor != functor)
106  {
107  m_Functor = functor;
108  this->Modified();
109  }
110  }
111 
114  {
115  m_Ei = ei;
116  this->GetFunctor().SetEi(ei);
117  this->Modified();
118  }
119 
122  {
123  m_Er = er;
124  this->GetFunctor().SetEr(er);
125  this->Modified();
126  }
128 
130  itkSetMacro(PsiI, double);
131  itkGetMacro(PsiI, double);
132 
134  itkSetMacro(KhiI, double);
135  itkGetMacro(KhiI, double);
136 
138  itkSetMacro(PsiR, double);
139  itkGetMacro(PsiR, double);
140 
142  itkSetMacro(KhiR, double);
143  itkGetMacro(KhiR, double);
144 
146  itkSetMacro(EmissionH, bool);
147  itkGetMacro(EmissionH, bool);
148 
150  itkSetMacro(EmissionV, bool);
151  itkGetMacro(EmissionV, bool);
152 
154  itkSetMacro(Mode, int);
155  itkGetMacro(Mode, int);
156 
158  itkSetMacro(Gain, double);
159  itkGetMacro(Gain, double);
160 
162  void ForceCoPolar();
163 
165  void ForceCrossPolar();
166 
167 protected:
170 
173 
184  void GenerateOutputInformation() override;
185 
186  void BeforeThreadedGenerateData() override;
187 
198  void ThreadedGenerateData(const OutputImageRegionType& outputRegionForThread,
199  itk::ThreadIdType threadId) override;
200 
202  void ComputeElectromagneticFields();
203 
205  void VerifyAndForceInputs();
206 
207  void PrintSelf(std::ostream& os, itk::Indent indent) const override;
208 
209 private:
210  MultiChannelsPolarimetricSynthesisFilter(const Self &) = delete;
211 
213  double m_PsiI;
214 
216  double m_KhiI;
217 
219  double m_PsiR;
220 
222  double m_KhiR;
223 
225  double m_Gain;
226 
228  int m_Mode;
229 
232 
235 
238 
240 
244 
245 };
246 
247 } // end namespace otb
248 
249 #ifndef OTB_MANUAL_INSTANTIATION
251 #endif
252 
253 #endif
This class computes the polarimetric synthesis from two to four radar images, depening on the polarim...
unsigned int ThreadIdType
itk::InPlaceImageFilter< TInputImage, TOutputImage > Superclass