OTB  10.0.0
Orfeo Toolbox
otbGenericRSResampleImageFilter.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 otbGenericRSResampleImageFilter_h
22 #define otbGenericRSResampleImageFilter_h
23 
26 #include <string>
27 
28 namespace otb
29 {
30 
55 template <class TInputImage, class TOutputImage>
56 class ITK_EXPORT GenericRSResampleImageFilter : public itk::ImageToImageFilter<TInputImage, TOutputImage>
57 {
58 public:
61  typedef itk::ImageToImageFilter<TInputImage, TOutputImage> Superclass;
62  typedef itk::SmartPointer<Self> Pointer;
63  typedef itk::SmartPointer<const Self> ConstPointer;
64 
66  itkNewMacro(Self);
67 
69  itkTypeMacro(GenericRSResampleImageFilter, itk::ImageToImageFilter);
70 
72  typedef TInputImage InputImageType;
73  typedef TOutputImage OutputImageType;
74  typedef typename OutputImageType::InternalPixelType OutputInternalPixelType;
75  typedef typename OutputImageType::PointType OutputPointType;
76 
87 
91 
94 
95 
101 
102  typedef itk::ImageBase<OutputImageType::ImageDimension> ImageBaseType;
103 
105  otbSetObjectMemberMacro(Resampler, DisplacementFieldSpacing, SpacingType);
106 
107  otbGetObjectMemberConstReferenceMacro(Resampler, DisplacementFieldSpacing, SpacingType);
108 
112  void SetOutputOrigin(const OriginType& origin)
113  {
114  m_Resampler->SetOutputOrigin(origin);
115  this->Modified();
116  }
119 
121  otbSetObjectMemberMacro(Resampler, OutputStartIndex, IndexType);
122  otbGetObjectMemberConstReferenceMacro(Resampler, OutputStartIndex, IndexType);
124 
126  otbSetObjectMemberMacro(Resampler, OutputSize, SizeType);
129 
131  otbSetObjectMemberMacro(Resampler, OutputSpacing, SpacingType);
134 
136  void SetInterpolator(InterpolatorType* interpolator)
137  {
138  m_Resampler->SetInterpolator(interpolator);
139  this->Modified();
140  }
141  otbGetObjectMemberConstMacro(Resampler, Interpolator, const InterpolatorType*);
143 
145  otbSetObjectMemberMacro(Resampler, EdgePaddingValue, typename OutputImageType::PixelType);
146  otbGetObjectMemberMacro(Resampler, EdgePaddingValue, typename OutputImageType::PixelType);
148 
154  void SetInputProjectionRef(const std::string& ref)
155  {
156  m_Transform->SetOutputProjectionRef(ref);
157  this->Modified();
158  }
160 
161  std::string GetInputProjectionRef() const
162  {
163  return m_Transform->GetOutputProjectionRef();
164  }
165 
166  void SetOutputProjectionRef(const std::string& ref)
167  {
168  m_Transform->SetInputProjectionRef(ref);
169  this->Modified();
170  }
171 
172  std::string GetOutputProjectionRef() const
173  {
174  return m_Transform->GetInputProjectionRef();
175  }
176 
180  {
181  m_Transform->SetOutputImageMetadata(imd);
182  this->Modified();
183  }
185 
187  {
188  return m_Transform->GetOutputImageMetadata();
189  }
190 
194  {
195  m_Transform->SetInputImageMetadata(imd);
196  this->Modified();
197  }
199 
201  {
202  return m_Transform->GetInputImageMetadata();
203  }
204 
206  void SetOutputParametersFromImage(const ImageBaseType* image);
207 
211  template <class TImageType>
212  void SetOutputParametersFromImage(const TImageType* image);
213 
217  void SetOutputParametersFromMap(const std::string& map, const SpacingType& spacing);
218 
222  void SetOutputParametersFromMap(const std::string& projectionRef);
223 
225  void SetInputRpcGridSize(const SizeType& gridSize)
226  {
227  m_InputRpcEstimator->SetGridSize(gridSize);
228  this->Modified();
229  }
230 
232  void SetInputRpcGridSize(unsigned int gridSize)
233  {
234  m_InputRpcEstimator->SetGridSize(gridSize);
235  this->Modified();
236  }
237 
240  {
241  return m_InputRpcEstimator->GetGridSize();
242  }
243 
245  itkSetMacro(EstimateInputRpcModel, bool);
246  itkGetMacro(EstimateInputRpcModel, bool);
247  itkBooleanMacro(EstimateInputRpcModel);
249 
251  void SetOutputRpcGridSize(const SizeType& gridSize)
252  {
253  m_OutputRpcEstimator->SetGridSize(gridSize);
254  this->Modified();
255  }
256 
258  void SetOutputRpcGridSize(unsigned int gridSize)
259  {
260  m_OutputRpcEstimator->SetGridSize(gridSize);
261  this->Modified();
262  }
263 
266  {
267  return m_OutputRpcEstimator->GetGridSize();
268  }
269 
271  itkSetMacro(EstimateOutputRpcModel, bool);
272  itkGetMacro(EstimateOutputRpcModel, bool);
273  itkBooleanMacro(EstimateOutputRpcModel);
275 
277  void SetDisplacementFilterNumberOfThreads(unsigned int nbThread)
278  {
279  m_Resampler->SetDisplacementFilterNumberOfThreads(nbThread);
280  }
281 
283  void PropagateRequestedRegion(itk::DataObject* output) override;
284 
285 protected:
289 
290  void GenerateData() override;
291 
292  void GenerateOutputInformation() override;
293 
294  virtual void UpdateTransform();
295 
296  void PrintSelf(std::ostream& os, itk::Indent indent) const override;
297 
298 private:
300  void operator=(const Self&) = delete;
301 
302  // Method to estimate the input & output rpc model
303  void EstimateOutputRpcModel();
304  void EstimateInputRpcModel();
305 
306  // boolean that allow the estimation of the input rpc model
310 
311  // Filters pointers
316 };
317 
318 } // namespace otb
319 
320 #ifndef OTB_MANUAL_INSTANTIATION
322 #endif
323 
324 #endif
This class is a composite filter that allows you to project an input image from any coordinate system...
void SetOutputProjectionRef(const std::string &ref)
itk::ImageToImageFilter< TInputImage, TOutputImage > Superclass
otbSetObjectMemberMacro(Resampler, DisplacementFieldSpacing, SpacingType)
otbGetObjectMemberConstReferenceMacro(Resampler, OutputSpacing, SpacingType)
InputRpcModelEstimatorPointerType m_InputRpcEstimator
otbGetObjectMemberMacro(Resampler, EdgePaddingValue, typename OutputImageType::PixelType)
OutputRpcModelEstimatorPointerType m_OutputRpcEstimator
itk::ImageBase< OutputImageType::ImageDimension > ImageBaseType
otbGetObjectMemberConstReferenceMacro(Resampler, OutputStartIndex, IndexType)
void SetInputImageMetadata(const ImageMetadata *imd)
OutputImageType::InternalPixelType OutputInternalPixelType
otbGetObjectMemberConstMacro(Resampler, Interpolator, const InterpolatorType *)
otbSetObjectMemberMacro(Resampler, EdgePaddingValue, typename OutputImageType::PixelType)
ResamplerType::InterpolatorType InterpolatorType
void SetInterpolator(InterpolatorType *interpolator)
OutputRpcModelEstimatorType::Pointer OutputRpcModelEstimatorPointerType
void operator=(const Self &)=delete
otbGetObjectMemberConstReferenceMacro(Resampler, DisplacementFieldSpacing, SpacingType)
StreamingResampleImageFilter< InputImageType, OutputImageType > ResamplerType
void SetInputRpcGridSize(const SizeType &gridSize)
otbGetObjectMemberConstReferenceMacro(Resampler, OutputSize, SizeType)
void SetOutputImageMetadata(const ImageMetadata *imd)
void SetOutputRpcGridSize(const SizeType &gridSize)
otbSetObjectMemberMacro(Resampler, OutputSpacing, SpacingType)
GenericRSResampleImageFilter(const Self &)=delete
otbSetObjectMemberMacro(Resampler, OutputStartIndex, IndexType)
InputRpcModelEstimatorType::Pointer InputRpcModelEstimatorPointerType
PhysicalToRPCSensorModelImageFilter< InputImageType > InputRpcModelEstimatorType
otbGetObjectMemberConstReferenceMacro(Resampler, OutputOrigin, OriginType)
PhysicalToRPCSensorModelImageFilter< OutputImageType > OutputRpcModelEstimatorType
void SetInputProjectionRef(const std::string &ref)
GenericRSTransformType::Pointer GenericRSTransformPointerType
itk::SmartPointer< const Self > ConstPointer
otbSetObjectMemberMacro(Resampler, OutputSize, SizeType)
void SetOutputOrigin(const OriginType &origin)
void SetDisplacementFilterNumberOfThreads(unsigned int nbThread)
This is the class to handle generic remote sensing transform.
itk::SmartPointer< Self > Pointer
Generic class containing image metadata used in OTB.
This filter estimates a RPC sensor models from a physical model.
This class is a composite filter resampling an input image by setting a transform....
DisplacementFieldGeneratorType::SizeType SizeType
DisplacementFieldGeneratorType::RegionType RegionType
DisplacementFieldGeneratorType::OriginType OriginType
DisplacementFieldGeneratorType::TransformType TransformType
itk::InterpolateImageFunction< InputImageType, TInterpolatorPrecisionType > InterpolatorType
DisplacementFieldGeneratorType::SpacingType SpacingType
DisplacementFieldGeneratorType::IndexType IndexType
itk::ImageBase< 2 > ImageBaseType
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.