OTB  10.0.0
Orfeo Toolbox
otbWindowedSincInterpolateImageLanczosFunction.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 otbWindowedSincInterpolateImageLanczosFunction_h
22 #define otbWindowedSincInterpolateImageLanczosFunction_h
23 
25 #include "vnl/vnl_math.h"
26 
27 namespace otb
28 {
29 
30 namespace Function
31 {
42 template <class TInput = double, class TOutput = double>
44 {
45 public:
47  {
48  } // default factor is 1 at construction
49  void SetRadius(unsigned int radius)
50  {
51  m_Radius = radius;
52  m_Factor = CONST_PI / static_cast<double>(radius);
53  }
54  unsigned int GetRadius() const
55  {
56  return m_Radius;
57  }
58  double GetFactor()
59  {
60  return m_Factor;
61  }
63 
64  inline TOutput operator()(const TInput& A) const
65  {
66  double x = static_cast<double>(A);
67  double px = CONST_PI * x;
68  double temp;
69  if (x == 0.0)
70  {
71  temp = 1.0;
72  }
73  else
74  {
75  double z = m_Factor * x;
76  temp = std::sin(z) / z;
77  }
78  return (x == 0.0) ? static_cast<TOutput>(temp) : static_cast<TOutput>(temp * std::sin(px) / px);
79  }
80 
81 private:
82  unsigned int m_Radius;
83  // Equal to \f$ \frac{\pi}{m} \f$
84  double m_Factor;
85 };
86 } // namespace Function
87 
103 template <class TInputImage, class TBoundaryCondition = itk::ConstantBoundaryCondition<TInputImage>, class TCoordRep = double,
104  class TInputInterpolator = double, class TOutputInterpolator = double>
106  : public WindowedSincInterpolateImageFunctionBase<TInputImage, typename Function::LanczosWindowFunction<TInputInterpolator, TOutputInterpolator>,
107  TBoundaryCondition, TCoordRep>
108 {
109 public:
110 
114  TBoundaryCondition, TCoordRep>
116  typedef itk::SmartPointer<Self> Pointer;
117  typedef itk::SmartPointer<const Self> ConstPointer;
118 
121 
123  itkNewMacro(Self);
124 
126  typedef typename Superclass::InputImageType InputImageType;
127  typedef typename Superclass::OutputType OutputType;
128 
130  itkStaticConstMacro(ImageDimension, unsigned int, Superclass::ImageDimension);
131 
133  typedef typename Superclass::IndexType IndexType;
134  typedef typename Superclass::SizeType SizeType;
135  typedef typename Superclass::RealType RealType;
136  typedef typename Superclass::IteratorType IteratorType;
137  typedef typename Superclass::ContinuousIndexType ContinuousIndexType;
138 
139 protected:
142  {
143  }
144  void PrintSelf(std::ostream& os, itk::Indent indent) const override
145  {
146  Superclass::PrintSelf(os, indent);
147  }
148 
149 private:
151  void operator=(const Self&) = delete;
152 };
153 
154 } // end namespace otb
155 
156 #endif
Generic interpolation of an otb::Image.
Use the WindowedSincInterpolateImageFunctionBase with a Lanczos Function.
WindowedSincInterpolateImageLanczosFunction(const Self &)=delete
WindowedSincInterpolateImageFunctionBase< TInputImage, typename Function::LanczosWindowFunction< TInputInterpolator, TOutputInterpolator >, TBoundaryCondition, TCoordRep > Superclass
void PrintSelf(std::ostream &os, itk::Indent indent) const override
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.
constexpr double CONST_PI
Definition: otbMath.h:49