OTB  7.4.0
Orfeo Toolbox
otbSARStreamingDEMCheckFilter.h
Go to the documentation of this file.
1 /*
2  * Copyright (C) 2005-2018 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 otbSARStreamingDEMCheckFilter_h
22 #define otbSARStreamingDEMCheckFilter_h
23 
25 #include "itkNumericTraits.h"
26 #include "itkArray.h"
27 #include "itkSimpleDataObjectDecorator.h"
28 #include "itkPoint.h"
29 
31 
33 #include "otbImageKeywordlist.h"
34 
35 #include <complex>
36 #include <cmath>
37 
38 namespace otb
39 {
40 
53 template<class TInputImage>
54 class ITK_EXPORT PersistentDEMCheckFilter :
55  public PersistentImageFilter<TInputImage, TInputImage>
56 {
57 public:
61  typedef itk::SmartPointer<Self> Pointer;
62  typedef itk::SmartPointer<const Self> ConstPointer;
63 
65  itkNewMacro(Self);
66 
69 
71  typedef TInputImage ImageType;
72  typedef typename TInputImage::Pointer InputImagePointer;
73 
74  typedef typename TInputImage::RegionType RegionType;
75  typedef typename TInputImage::SizeType SizeType;
77  typedef typename TInputImage::PixelType PixelType;
78  typedef typename ImageType::IndexValueType IndexValueType;
79  typedef typename ImageType::SizeValueType SizeValueType;
80 
81  itkStaticConstMacro(InputImageDimension, unsigned int,
82  TInputImage::ImageDimension);
83 
85  itkStaticConstMacro(ImageDimension, unsigned int,
86  TInputImage::ImageDimension);
87 
89  // typedef typename itk::NumericTraits<PixelType>::RealType RealType;
91 
93  typedef typename itk::DataObject::Pointer DataObjectPointer;
94  typedef itk::ProcessObject::DataObjectPointerArraySizeType DataObjectPointerArraySizeType;
95 
97  typedef itk::SimpleDataObjectDecorator<RealType> RealObjectType;
98  typedef itk::SimpleDataObjectDecorator<long int> LongObjectType;
99  typedef itk::SimpleDataObjectDecorator<PixelType> PixelObjectType;
100 
101  // Set SAR Image Size and Origin
102  void setSARSizeAndOrigin(int sizeC, int sizeL, int originC, int originL)
103  {
104  m_SizeSARC = sizeC;
105  m_SizeSARL = sizeL;
106  m_OriginSARC = originC;
107  m_OriginSARL = originL;
108  }
109 
111  long int GetCounter() const
112  {
113  return this->GetCounter_Output()->Get();
114  }
115  LongObjectType* GetCounter_Output();
116  const LongObjectType* GetCounter_Output() const;
118 
121  DataObjectPointer MakeOutput(DataObjectPointerArraySizeType idx) ITK_OVERRIDE;
122  using Superclass::MakeOutput;
123 
127  void Synthetize(void) ITK_OVERRIDE;
128  void Reset(void) ITK_OVERRIDE;
130 
131  itkSetMacro(UserIgnoredValue, RealType);
132  itkGetMacro(UserIgnoredValue, RealType);
133 
134 protected:
136  ~PersistentDEMCheckFilter() ITK_OVERRIDE {}
137  void PrintSelf(std::ostream& os, itk::Indent indent) const ITK_OVERRIDE;
138 
139 
143  virtual void GenerateInputRequestedRegion() ITK_OVERRIDE;
144 
146  void ThreadedGenerateData(const RegionType&
147  outputRegionForThread,
148  itk::ThreadIdType threadId) ITK_OVERRIDE;
149 
153  void AllocateOutputs() ITK_OVERRIDE;
154  void GenerateOutputInformation() ITK_OVERRIDE;
156 
157 private:
158  PersistentDEMCheckFilter(const Self &); //purposely not implemented
159  void operator =(const Self&); //purposely not implemented
160 
161  itk::Array<long int> m_Thread_Counter;
162 
163  /* Ignored values */
165 
166  // Main characterics for SAR Image (Origin and Size)
171 }; // end of class PersistentDEMCheckFilter
172 
173 /*===========================================================================*/
174 
209 template<class TInputImage>
210 class ITK_EXPORT SARStreamingDEMCheckFilter :
211  public PersistentFilterStreamingDecorator<PersistentDEMCheckFilter<TInputImage> >
212 {
213 public:
218  typedef itk::SmartPointer<Self> Pointer;
219  typedef itk::SmartPointer<const Self> ConstPointer;
220 
222  itkNewMacro(Self);
223 
225  itkTypeMacro(StreamingDEMCheckFilter, PersistentFilterStreamingDecorator);
226 
227  typedef typename Superclass::FilterType StatFilterType;
230  typedef TInputImage InputImageType;
231 
233  typedef itk::SimpleDataObjectDecorator<RealType> RealObjectType;
234  typedef itk::SimpleDataObjectDecorator<long int> LongObjectType;
235  typedef itk::SimpleDataObjectDecorator<PixelType> PixelObjectType;
236 
238  typedef TInputImage ImageType;
239  typedef typename TInputImage::Pointer InputImagePointer;
240 
241  typedef typename TInputImage::RegionType RegionType;
242  typedef typename TInputImage::SizeType SizeType;
244 
245  // Define Point2DType and Point3DType
246  using Point2DType = itk::Point<double,2>;
247  using Point3DType = itk::Point<double,3>;
248 
249  using Superclass::SetInput;
250  void SetInput(InputImageType * input)
251  {
252  this->GetFilter()->SetInput(input);
253  }
254  const InputImageType * GetInput()
255  {
256  return this->GetFilter()->GetInput();
257  }
258 
259  void setSARSizeAndOrigin(int sizeC, int sizeL, int originC, int originL)
260  {
261  this->GetFilter()->setSARSizeAndOrigin(sizeC, sizeL, originC, originL);
262  }
263 
265  long int GetCounter() const
266  {
267  return this->GetFilter()->GetCounter_Output()->Get();
268  }
269  LongObjectType* GetCounter_Output()
270  {
271  return this->GetFilter()->GetCounter_Output();
272  }
273  const LongObjectType* GetCounter_Output() const
274  {
275  return this->GetFilter()->GetCounter_Output();
276  }
278 
279  // Get some information abour DEM (number of DEM points with a projection into SAR geometry)
280  void GetDEMCheck(long int & nbDEMPts_into_SAR)
281  {
282  nbDEMPts_into_SAR = this->GetCounter();
283  }
284 
285  otbSetObjectMemberMacro(Filter, UserIgnoredValue, RealType);
286  otbGetObjectMemberMacro(Filter, UserIgnoredValue, RealType);
287 
288 protected:
291 
293  ~SARStreamingDEMCheckFilter() ITK_OVERRIDE {}
294 
295 private:
296  SARStreamingDEMCheckFilter(const Self &); //purposely not implemented
297  void operator =(const Self&); //purposely not implemented
298 
299 };
300 
301 } // end namespace otb
302 
303 #ifndef OTB_MANUAL_INSTANTIATION
304 #include "otbSARStreamingDEMCheckFilter.txx"
305 #endif
306 
307 #endif
itk::Size< Monteverdi_DIMENSION > SizeType
Definition: mvdTypes.h:137
itk::Index< Monteverdi_DIMENSION > IndexType
Definition: mvdTypes.h:133
itk::SmartPointer< const Self > ConstPointer
itk::SimpleDataObjectDecorator< long int > LongObjectType
itk::NumericTraits< double >::RealType RealType
itk::SimpleDataObjectDecorator< RealType > RealObjectType
void GetDEMCheck(long int &nbDEMPts_into_SAR)
PersistentImageFilter< TInputImage, TInputImage > Superclass
#define otbGetObjectMemberMacro(object, name, type)
Definition: otbMacro.h:89
PersistentFilterStreamingDecorator< PersistentDEMCheckFilter< TInputImage > > Superclass
const LongObjectType * GetCounter_Output() const
static const std::string Filter
itk::SimpleDataObjectDecorator< PixelType > PixelObjectType
itk::ProcessObject::DataObjectPointerArraySizeType DataObjectPointerArraySizeType
itk::SimpleDataObjectDecorator< PixelType > PixelObjectType
itk::SmartPointer< const Self > ConstPointer
itk::SimpleDataObjectDecorator< RealType > RealObjectType
Monteverdi_FLOATING_TYPE RealType
Definition: mvdTypes.h:84
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.
void setSARSizeAndOrigin(int sizeC, int sizeL, int originC, int originL)
itk::SimpleDataObjectDecorator< long int > LongObjectType
#define otbSetObjectMemberMacro(object, name, type)
Definition: otbMacro.h:79
Retrive some pixels from the input DEM.
This class streams the whole input image through the PersistentDEMCheckFilter.
void setSARSizeAndOrigin(int sizeC, int sizeL, int originC, int originL)
This filter is the base class for all filter persisting data through multiple update. For instance, a filter computing global statistics on an image with streaming capabilities will have to keep the temporary results for each streamed piece of the image in order to synthesize the global statistics at the end. This filter is an itk::ImageToImageFilter, providing two additional methods. The first one, Synthetize(), allows the user to synthesize temporary data produced by the multiple updates on different pieces of the image to the global result. The second one, Reset(), allows the user to reset the temporary data for a new input image for instance.
This filter link a persistent filter with a StreamingImageVirtualWriter.