OTB  10.0.0
Orfeo Toolbox
otbStreamingMatrixTransposeMatrixImageFilter.h
Go to the documentation of this file.
1 /*
2  * Copyright (C) 1999-2011 Insight Software Consortium
3  * Copyright (C) 2005-2024 Centre National d'Etudes Spatiales (CNES)
4  *
5  * This file is part of Orfeo Toolbox
6  *
7  * https://www.orfeo-toolbox.org/
8  *
9  * Licensed under the Apache License, Version 2.0 (the "License");
10  * you may not use this file except in compliance with the License.
11  * You may obtain a copy of the License at
12  *
13  * http://www.apache.org/licenses/LICENSE-2.0
14  *
15  * Unless required by applicable law or agreed to in writing, software
16  * distributed under the License is distributed on an "AS IS" BASIS,
17  * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
18  * See the License for the specific language governing permissions and
19  * limitations under the License.
20  */
21 
22 #ifndef otbStreamingMatrixTransposeMatrixImageFilter_h
23 #define otbStreamingMatrixTransposeMatrixImageFilter_h
24 
26 #include "itkSimpleDataObjectDecorator.h"
27 #include "otbStreamingTraits.h"
28 #include "itkVariableSizeMatrix.h"
29 #include "itkVariableLengthVector.h"
31 
32 namespace otb
33 {
34 
54 template <class TInputImage, class TInputImage2>
55 class ITK_EXPORT PersistentMatrixTransposeMatrixImageFilter : 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  // First Input
72  typedef TInputImage ImageType;
73  typedef typename TInputImage::Pointer InputImagePointer;
74  typedef typename TInputImage::RegionType RegionType;
75  typedef typename TInputImage::SizeType SizeType;
76  typedef typename TInputImage::IndexType IndexType;
77  typedef typename TInputImage::PixelType PixelType;
78  typedef typename TInputImage::InternalPixelType InternalPixelType;
79 
80  typedef typename TInputImage2::IndexType IndexType2;
81  typedef typename TInputImage2::PixelType PixelType2;
82  typedef typename TInputImage2::InternalPixelType InternalPixelType2;
83 
84  itkStaticConstMacro(InputImageDimension, unsigned int, TInputImage::ImageDimension);
85 
86  itkSetMacro(UsePadFirstInput, bool);
87  itkGetMacro(UsePadFirstInput, bool);
88  itkSetMacro(UsePadSecondInput, bool);
89  itkGetMacro(UsePadSecondInput, bool);
90 
92  itkStaticConstMacro(ImageDimension, unsigned int, TInputImage::ImageDimension);
93 
95  // First Input
96  typedef double RealType;
97  typedef itk::VariableLengthVector<RealType> RealPixelType;
98 
100  typedef typename itk::DataObject::Pointer DataObjectPointer;
101  typedef itk::ProcessObject::DataObjectPointerArraySizeType DataObjectPointerArraySizeType;
102 
104  typedef typename itk::Array<long> ArrayLongPixelType;
105  typedef typename itk::VariableSizeMatrix<RealType> MatrixType;
106  typedef typename std::vector<MatrixType> ArrayMatrixType;
107  typedef typename std::vector<RealPixelType> ArrayRealPixelType;
108  typedef typename std::vector<PixelType> ArrayPixelType;
109  typedef itk::SimpleDataObjectDecorator<RealPixelType> RealPixelObjectType;
110  typedef itk::SimpleDataObjectDecorator<PixelType> PixelObjectType;
111  typedef itk::SimpleDataObjectDecorator<MatrixType> MatrixObjectType;
112 
115  {
116  return this->GetResultOutput()->Get();
117  }
118  MatrixObjectType* GetResultOutput();
119  const MatrixObjectType* GetResultOutput() const;
121 
125  DataObjectPointer MakeOutput(DataObjectPointerArraySizeType idx) override;
126  using Superclass::MakeOutput;
127 
131  void AllocateOutputs() override;
132  void GenerateOutputInformation() override;
133  // Override since the filter needs all the data for the algorithm
134  void GenerateInputRequestedRegion() override;
135  void Reset(void) override;
136  void Synthetize(void) override;
138 
140  void SetFirstInput(const TInputImage* input1)
141  {
142  this->SetInput(0, input1);
143  }
144  void SetSecondInput(const TInputImage2* input2)
145  {
146  this->SetInput(1, input2);
147  }
149 
150  const TInputImage* GetFirstInput()
151  {
152  if (this->GetNumberOfInputs() < 1)
153  {
154  return nullptr;
155  }
156  else
157  return (static_cast<const TInputImage*>(this->itk::ProcessObject::GetInput(0)));
158  }
159 
160  const TInputImage2* GetSecondInput()
161  {
162  if (this->GetNumberOfInputs() < 2)
163  {
164  return nullptr;
165  }
166  else
167  return (static_cast<const TInputImage2*>(this->itk::ProcessObject::GetInput(1)));
168  }
169 
170 protected:
173  {
174  }
175  void PrintSelf(std::ostream& os, itk::Indent indent) const override;
177  void ThreadedGenerateData(const RegionType& outputRegionForThread, itk::ThreadIdType threadId) override;
178 
179 private:
181  void operator=(const Self&) = delete;
182 
186 
188  unsigned int m_NumberOfComponents1;
189  unsigned int m_NumberOfComponents2;
190 }; // end of class
191 
215 template <class TInputImage1, class TInputImage2>
217  : public PersistentFilterStreamingDecorator<PersistentMatrixTransposeMatrixImageFilter<TInputImage1, TInputImage2>>
218 {
219 public:
223  typedef itk::SmartPointer<Self> Pointer;
224  typedef itk::SmartPointer<const Self> ConstPointer;
225 
227  itkNewMacro(Self);
228 
231 
234  typedef itk::SimpleDataObjectDecorator<MatrixType> MatrixObjectType;
235 
236  typedef TInputImage1 InputImage1Type;
237  typedef TInputImage2 InputImage2Type;
238 
241  {
242  this->GetFilter()->SetFirstInput(input1);
243  }
245  {
246  this->GetFilter()->SetSecondInput(input2);
247  }
248  void SetUsePadFirstInput(bool pad)
249  {
250  this->GetFilter()->SetUsePadFirstInput(pad);
251  }
252  void SetUsePadSecondInput(bool pad)
253  {
254  this->GetFilter()->SetUsePadSecondInput(pad);
255  }
257  {
258  return this->GetFilter()->GetUsePadFirstInput();
259  }
261  {
262  return this->GetFilter()->GetUsePadSecondInput();
263  }
264 
266  MatrixType GetResult(void) const
267  {
268  return this->GetResultOutput()->Get();
269  }
271  {
272  return this->GetFilter()->GetResultOutput();
273  }
275  {
276  return this->GetFilter()->GetResultOutput();
277  }
279 
280 protected:
283 
286  {
287  }
288 
289 private:
291  void operator=(const Self&) = delete;
292 };
293 
294 } // end namespace otb
295 
296 #ifndef OTB_MANUAL_INSTANTIATION
298 #endif
299 
300 #endif
This filter link a persistent filter with a StreamingImageVirtualWriter.
This filter is the base class for all filter persisting data through multiple update....
PersistentMatrixTransposeMatrixImageFilter(const Self &)=delete
itk::SimpleDataObjectDecorator< RealPixelType > RealPixelObjectType
itk::ProcessObject::DataObjectPointerArraySizeType DataObjectPointerArraySizeType
PersistentImageFilter< TInputImage, TInputImage > Superclass
This class streams the whole input image through the PersistentMatrixTransposeMatrixImageFilter.
StreamingMatrixTransposeMatrixImageFilter(const Self &)=delete
PersistentFilterStreamingDecorator< PersistentMatrixTransposeMatrixImageFilter< TInputImage1, TInputImage2 > > Superclass
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.