OTB  10.0.0
Orfeo Toolbox
otbStreamingResampleImageFilter.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 otbStreamingResampleImageFilter_h
22 #define otbStreamingResampleImageFilter_h
23 
24 #include "itkImageToImageFilter.h"
26 #include "itkTransformToDisplacementFieldSource.h"
27 #include "itkLinearInterpolateImageFunction.h"
28 #include "otbImage.h"
29 #include "itkVector.h"
30 
31 #include "otbMacro.h"
32 
33 namespace otb
34 {
35 
58 template <class TInputImage, class TOutputImage, class TInterpolatorPrecisionType = double>
59 class ITK_EXPORT StreamingResampleImageFilter : public itk::ImageToImageFilter<TInputImage, TOutputImage>
60 {
61 public:
64  typedef itk::ImageToImageFilter<TInputImage, TOutputImage> Superclass;
65  typedef itk::SmartPointer<Self> Pointer;
66  typedef itk::SmartPointer<const Self> ConstPointer;
67 
69  itkNewMacro(Self);
70 
72  itkTypeMacro(StreamingResampleImageFilter, itk::ImageToImageFilter);
73 
75  typedef TInputImage InputImageType;
76  typedef TOutputImage OutputImageType;
77 
79  typedef itk::Vector<double, TOutputImage::ImageDimension> DisplacementType;
81 
84 
86  typedef itk::TransformToDisplacementFieldSource<DisplacementFieldType, double> DisplacementFieldGeneratorType;
87  typedef typename DisplacementFieldGeneratorType::TransformType TransformType;
88  typedef typename DisplacementFieldGeneratorType::SizeType SizeType;
89  typedef typename DisplacementFieldGeneratorType::SpacingType SpacingType;
90  typedef typename DisplacementFieldGeneratorType::OriginType OriginType;
91  typedef typename DisplacementFieldGeneratorType::IndexType IndexType;
92  typedef typename DisplacementFieldGeneratorType::RegionType RegionType;
93 
95  typedef itk::InterpolateImageFunction<InputImageType, TInterpolatorPrecisionType> InterpolatorType;
96  typedef typename InterpolatorType::Pointer InterpolatorPointerType;
97  typedef itk::LinearInterpolateImageFunction<InputImageType, TInterpolatorPrecisionType> DefaultInterpolatorType;
98 
100  typedef itk::ImageBase<OutputImageType::ImageDimension> ImageBaseType;
101 
102 
104  void SetTransform(TransformType* transform)
105  {
106  m_DisplacementFilter->SetTransform(transform);
107  this->Modified();
108  }
111 
113  void SetDisplacementFieldSpacing(SpacingType spacing);
114 
116  {
117  return m_SignedOutputSpacing;
118  };
119 
121  // Output Origin
122  void SetOutputOrigin(const OriginType& origin)
123  {
124  m_DisplacementFilter->SetOutputOrigin(origin);
125  m_WarpFilter->SetOutputOrigin(origin);
126  this->Modified();
127  }
130 
131  // Output Start index
132  otbSetObjectMemberMacro(WarpFilter, OutputStartIndex, IndexType);
133  otbGetObjectMemberConstReferenceMacro(WarpFilter, OutputStartIndex, IndexType);
134 
135  // Output Size
136  otbSetObjectMemberMacro(WarpFilter, OutputSize, SizeType);
138 
139  // Output Spacing
140  otbSetObjectMemberMacro(WarpFilter, OutputSpacing, SpacingType);
142 
144  void SetInterpolator(InterpolatorType* interpolator)
145  {
146  m_WarpFilter->SetInterpolator(interpolator);
147  this->Modified();
148  }
149  otbGetObjectMemberConstMacro(WarpFilter, Interpolator, const InterpolatorType*);
151 
153  otbSetObjectMemberMacro(WarpFilter, EdgePaddingValue, typename OutputImageType::PixelType);
154  otbGetObjectMemberMacro(WarpFilter, EdgePaddingValue, typename OutputImageType::PixelType);
156 
158  void SetOutputParametersFromImage(const ImageBaseType* image);
159 
160  /* Set number of threads for Deformation field generator*/
161  void SetDisplacementFilterNumberOfThreads(unsigned int nbThread)
162  {
163  m_DisplacementFilter->SetNumberOfWorkUnits(nbThread);
164  }
165 
167  void PropagateRequestedRegion(itk::DataObject* output) override;
168 
169 protected:
171 
174 
175  void GenerateData() override;
176 
177  void GenerateOutputInformation() override;
178 
179  void PrintSelf(std::ostream& os, itk::Indent indent) const override;
180 
181 private:
183  void operator=(const Self&) = delete;
184 
185  // We need this to respect ConstRef macro and to be compliant with itk positive
186  // spacing
188 
189  typename DisplacementFieldGeneratorType::Pointer m_DisplacementFilter;
191 };
192 
193 } // namespace otb
194 
195 #ifndef OTB_MANUAL_INSTANTIATION
197 #endif
198 
199 #endif
Creation of an "otb" image which contains metadata.
Definition: otbImage.h:92
This class is a composite filter resampling an input image by setting a transform....
itk::TransformToDisplacementFieldSource< DisplacementFieldType, double > DisplacementFieldGeneratorType
itk::ImageToImageFilter< TInputImage, TOutputImage > Superclass
const SpacingType & GetDisplacementFieldSpacing() const
otbGetObjectMemberMacro(WarpFilter, EdgePaddingValue, typename OutputImageType::PixelType)
otbSetObjectMemberMacro(WarpFilter, OutputSpacing, SpacingType)
otbSetObjectMemberMacro(WarpFilter, OutputSize, SizeType)
otbGetObjectMemberConstReferenceMacro(WarpFilter, OutputOrigin, OriginType)
itk::SmartPointer< const Self > ConstPointer
itk::LinearInterpolateImageFunction< InputImageType, TInterpolatorPrecisionType > DefaultInterpolatorType
void SetDisplacementFilterNumberOfThreads(unsigned int nbThread)
otbGetObjectMemberConstMacro(WarpFilter, Interpolator, const InterpolatorType *)
otb::Image< DisplacementType > DisplacementFieldType
otbGetObjectMemberConstMacro(DisplacementFilter, Transform, const TransformType *)
DisplacementFieldGeneratorType::SizeType SizeType
DisplacementFieldGeneratorType::RegionType RegionType
DisplacementFieldGeneratorType::Pointer m_DisplacementFilter
itk::ImageBase< OutputImageType::ImageDimension > ImageBaseType
otbGetObjectMemberConstReferenceMacro(WarpFilter, OutputStartIndex, IndexType)
DisplacementFieldGeneratorType::OriginType OriginType
void SetOutputOrigin(const OriginType &origin)
DisplacementFieldGeneratorType::TransformType TransformType
otbGetObjectMemberConstReferenceMacro(WarpFilter, OutputSpacing, SpacingType)
void SetInterpolator(InterpolatorType *interpolator)
itk::InterpolateImageFunction< InputImageType, TInterpolatorPrecisionType > InterpolatorType
DisplacementFieldGeneratorType::SpacingType SpacingType
itk::Vector< double, TOutputImage::ImageDimension > DisplacementType
DisplacementFieldGeneratorType::IndexType IndexType
otbSetObjectMemberMacro(WarpFilter, EdgePaddingValue, typename OutputImageType::PixelType)
otbSetObjectMemberMacro(WarpFilter, OutputStartIndex, IndexType)
otbGetObjectMemberConstReferenceMacro(WarpFilter, OutputSize, SizeType)
void operator=(const Self &)=delete
StreamingResampleImageFilter(const Self &)=delete
StreamingWarpImageFilter< InputImageType, OutputImageType, DisplacementFieldType > WarpImageFilterType
This class acts like the itk::WarpImageFilter, but it does not request the largest possible region of...
Class to overload method passed to virtual pure in ITK V4.
Definition: otbTransform.h:42
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.