ITK  6.0.0
Insight Toolkit
itkWindowedSincInterpolateImageFunction.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright NumFOCUS
4  *
5  * Licensed under the Apache License, Version 2.0 (the "License");
6  * you may not use this file except in compliance with the License.
7  * You may obtain a copy of the License at
8  *
9  * https://www.apache.org/licenses/LICENSE-2.0.txt
10  *
11  * Unless required by applicable law or agreed to in writing, software
12  * distributed under the License is distributed on an "AS IS" BASIS,
13  * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
14  * See the License for the specific language governing permissions and
15  * limitations under the License.
16  *
17  *=========================================================================*/
18 #ifndef itkWindowedSincInterpolateImageFunction_h
19 #define itkWindowedSincInterpolateImageFunction_h
20 
24 #include "itkMath.h"
25 
26 namespace itk
27 {
28 namespace Function
29 {
37 template <unsigned int VRadius, typename TInput = double, typename TOutput = double>
38 class ITK_TEMPLATE_EXPORT CosineWindowFunction
39 {
40 public:
41  inline TOutput
42  operator()(const TInput & A) const
43  {
44  return static_cast<TOutput>(std::cos(A * m_Factor));
45  }
48 private:
50  static const double m_Factor;
51 };
52 
60 template <unsigned int VRadius, typename TInput = double, typename TOutput = double>
61 class ITK_TEMPLATE_EXPORT HammingWindowFunction
62 {
63 public:
64  inline TOutput
65  operator()(const TInput & A) const
66  {
67  return static_cast<TOutput>(0.54 + 0.46 * std::cos(A * m_Factor));
68  }
71 private:
73  static const double m_Factor;
74 };
75 
83 template <unsigned int VRadius, typename TInput = double, typename TOutput = double>
84 class ITK_TEMPLATE_EXPORT WelchWindowFunction
85 {
86 public:
87  inline TOutput
88  operator()(const TInput & A) const
89  {
90  return static_cast<TOutput>(1.0 - A * m_Factor * A);
91  }
94 private:
96  static const double m_Factor;
97 };
98 
108 template <unsigned int VRadius, typename TInput = double, typename TOutput = double>
109 class ITK_TEMPLATE_EXPORT LanczosWindowFunction
110 {
111 public:
112  inline TOutput
113  operator()(const TInput & A) const
114  {
115  if (A == 0.0)
116  {
117  return static_cast<TOutput>(1.0);
118  }
119  const double z = m_Factor * A;
120  return static_cast<TOutput>(std::sin(z) / z);
121  }
124 private:
126  static const double m_Factor;
127 };
128 
136 template <unsigned int VRadius, typename TInput = double, typename TOutput = double>
137 class ITK_TEMPLATE_EXPORT BlackmanWindowFunction
138 {
139 public:
140  inline TOutput
141  operator()(const TInput & A) const
142  {
143  return static_cast<TOutput>(0.42 + 0.5 * std::cos(A * m_Factor1) + 0.08 * std::cos(A * m_Factor2));
144  }
147 private:
149  static const double m_Factor1;
150 
152  static const double m_Factor2;
153 };
154 } // namespace Function
155 
266 template <typename TInputImage,
267  unsigned int VRadius,
268  typename TWindowFunction = Function::HammingWindowFunction<VRadius>,
270  class TCoordinate = double>
271 class ITK_TEMPLATE_EXPORT WindowedSincInterpolateImageFunction
272  : public InterpolateImageFunction<TInputImage, TCoordinate>
273 {
274 public:
275  ITK_DISALLOW_COPY_AND_MOVE(WindowedSincInterpolateImageFunction);
281 
284 
286  itkOverrideGetNameOfClassMacro(WindowedSincInterpolateImageFunction);
287 
289  itkNewMacro(Self);
290 
292  using typename Superclass::OutputType;
293 
295  using typename Superclass::InputImageType;
296 
298  using typename Superclass::RealType;
299 
301  static constexpr unsigned int ImageDimension = Superclass::ImageDimension;
302 
304  using typename Superclass::IndexType;
305  using typename Superclass::IndexValueType;
306 
308  using typename Superclass::SizeType;
309 
311  using ImageType = TInputImage;
312 
314  using typename Superclass::ContinuousIndexType;
315 
316  void
317  SetInputImage(const ImageType * image) override;
318 
325  OutputType
326  EvaluateAtContinuousIndex(const ContinuousIndexType & index) const override;
327 
328  SizeType
329  GetRadius() const override
330  {
331  constexpr auto radius = SizeType::Filled(VRadius);
332  return radius;
333  }
334 
335 protected:
337  ~WindowedSincInterpolateImageFunction() override = default;
338  void
339  PrintSelf(std::ostream & os, Indent indent) const override;
340 
341 private:
342  // Internal type alias
344 
345  // Constant to store twice the radius
346  static constexpr unsigned int m_WindowSize{ 2 * VRadius };
347 
349  TWindowFunction m_WindowFunction{};
350 
352  static constexpr unsigned int m_OffsetTableSize = Math::UnsignedPower(m_WindowSize, ImageDimension);
353 
356  unsigned int m_OffsetTable[m_OffsetTableSize]{};
357 
359  unsigned int m_WeightOffsetTable[m_OffsetTableSize][ImageDimension]{};
360 
362  inline double
363  Sinc(double x) const
364  {
365  const double px = itk::Math::pi * x;
366 
367  return (x == 0.0) ? 1.0 : std::sin(px) / px;
368  }
369 };
370 } // namespace itk
371 
372 #ifndef ITK_MANUAL_INSTANTIATION
373 # include "itkWindowedSincInterpolateImageFunction.hxx"
374 #endif
375 
376 #endif // _itkWindowedSincInterpolateImageFunction_h
itk::Function::HammingWindowFunction::operator()
TOutput operator()(const TInput &A) const
Definition: itkWindowedSincInterpolateImageFunction.h:65
itk::WindowedSincInterpolateImageFunction::Sinc
double Sinc(double x) const
Definition: itkWindowedSincInterpolateImageFunction.h:363
itkConstNeighborhoodIterator.h
itk::Function::CosineWindowFunction::operator()
TOutput operator()(const TInput &A) const
Definition: itkWindowedSincInterpolateImageFunction.h:42
itk::Function::WelchWindowFunction::operator()
TOutput operator()(const TInput &A) const
Definition: itkWindowedSincInterpolateImageFunction.h:88
itk::Function::BlackmanWindowFunction::m_Factor2
static const double m_Factor2
Definition: itkWindowedSincInterpolateImageFunction.h:152
itk::GTest::TypedefsAndConstructors::Dimension2::SizeType
ImageBaseType::SizeType SizeType
Definition: itkGTestTypedefsAndConstructors.h:49
itk::Function::CosineWindowFunction
Window function for sinc interpolation.
Definition: itkWindowedSincInterpolateImageFunction.h:38
itk::SmartPointer< Self >
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::WindowedSincInterpolateImageFunction::ImageType
TInputImage ImageType
Definition: itkWindowedSincInterpolateImageFunction.h:311
itk::FunctionBase< Point< TCoordinate, TInputImage::ImageDimension >, NumericTraits< TInputImage::PixelType >::RealType >::OutputType
NumericTraits< TInputImage::PixelType >::RealType OutputType
Definition: itkFunctionBase.h:62
itk::IndexValueType
long IndexValueType
Definition: itkIntTypes.h:93
itk::Function::CosineWindowFunction::m_Factor
static const double m_Factor
Definition: itkWindowedSincInterpolateImageFunction.h:50
itk::GTest::TypedefsAndConstructors::Dimension2::IndexType
ImageBaseType::IndexType IndexType
Definition: itkGTestTypedefsAndConstructors.h:50
itk::LightObject
Light weight base class for most itk classes.
Definition: itkLightObject.h:55
itk::Function::WelchWindowFunction::m_Factor
static const double m_Factor
Definition: itkWindowedSincInterpolateImageFunction.h:96
itk::Math::UnsignedPower
constexpr TReturnType UnsignedPower(const uintmax_t base, const uintmax_t exponent) noexcept
Definition: itkMath.h:793
itk::InterpolateImageFunction::SizeType
typename InputImageType::SizeType SizeType
Definition: itkInterpolateImageFunction.h:79
itk::Function::BlackmanWindowFunction::operator()
TOutput operator()(const TInput &A) const
Definition: itkWindowedSincInterpolateImageFunction.h:141
itk::Function::HammingWindowFunction::m_Factor
static const double m_Factor
Definition: itkWindowedSincInterpolateImageFunction.h:73
itk::Function::LanczosWindowFunction
Window function for sinc interpolation.
Definition: itkWindowedSincInterpolateImageFunction.h:109
itk::WindowedSincInterpolateImageFunction
Use the windowed sinc function to interpolate.
Definition: itkWindowedSincInterpolateImageFunction.h:271
itk::Function::LanczosWindowFunction::operator()
TOutput operator()(const TInput &A) const
Definition: itkWindowedSincInterpolateImageFunction.h:113
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnatomicalOrientation.h:29
itk::ConstNeighborhoodIterator
Const version of NeighborhoodIterator, defining iteration of a local N-dimensional neighborhood of pi...
Definition: itkConstNeighborhoodIterator.h:51
itk::ContinuousIndex< TCoordinate, Self::ImageDimension >
itk::ZeroFluxNeumannBoundaryCondition
A function object that determines a neighborhood of values at an image boundary according to a Neuman...
Definition: itkZeroFluxNeumannBoundaryCondition.h:58
itkInterpolateImageFunction.h
itk::Math::pi
static constexpr double pi
Definition: itkMath.h:66
itk::Function::BlackmanWindowFunction::m_Factor1
static const double m_Factor1
Definition: itkWindowedSincInterpolateImageFunction.h:149
itk::WindowedSincInterpolateImageFunction::GetRadius
SizeType GetRadius() const override
Definition: itkWindowedSincInterpolateImageFunction.h:329
itk::Function::BlackmanWindowFunction
Window function for sinc interpolation.
Definition: itkWindowedSincInterpolateImageFunction.h:137
itkZeroFluxNeumannBoundaryCondition.h
itk::Function::LanczosWindowFunction::m_Factor
static const double m_Factor
Definition: itkWindowedSincInterpolateImageFunction.h:126
itkMath.h
itk::InterpolateImageFunction
Base class for all image interpolators.
Definition: itkInterpolateImageFunction.h:45
itk::Function::HammingWindowFunction
Window function for sinc interpolation.
Definition: itkWindowedSincInterpolateImageFunction.h:61
itk::Function::WelchWindowFunction
Window function for sinc interpolation.
Definition: itkWindowedSincInterpolateImageFunction.h:84