ITK  5.2.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  * http://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 
25 namespace itk
26 {
27 namespace Function
28 {
36 template <unsigned int VRadius, typename TInput = double, typename TOutput = double>
37 class ITK_TEMPLATE_EXPORT CosineWindowFunction
38 {
39 public:
40  inline TOutput
41  operator()(const TInput & A) const
42  {
43  return static_cast<TOutput>(std::cos(A * m_Factor));
44  }
46 
47 private:
49  static const double m_Factor;
50 };
51 
59 template <unsigned int VRadius, typename TInput = double, typename TOutput = double>
60 class ITK_TEMPLATE_EXPORT HammingWindowFunction
61 {
62 public:
63  inline TOutput
64  operator()(const TInput & A) const
65  {
66  return static_cast<TOutput>(0.54 + 0.46 * std::cos(A * m_Factor));
67  }
69 
70 private:
72  static const double m_Factor;
73 };
74 
82 template <unsigned int VRadius, typename TInput = double, typename TOutput = double>
83 class ITK_TEMPLATE_EXPORT WelchWindowFunction
84 {
85 public:
86  inline TOutput
87  operator()(const TInput & A) const
88  {
89  return static_cast<TOutput>(1.0 - A * m_Factor * A);
90  }
92 
93 private:
95  static const double m_Factor;
96 };
97 
107 template <unsigned int VRadius, typename TInput = double, typename TOutput = double>
108 class ITK_TEMPLATE_EXPORT LanczosWindowFunction
109 {
110 public:
111  inline TOutput
112  operator()(const TInput & A) const
113  {
114  if (A == 0.0)
115  {
116  return static_cast<TOutput>(1.0);
117  }
118  double z = m_Factor * A;
119  return static_cast<TOutput>(std::sin(z) / z);
120  }
122 
123 private:
125  static const double m_Factor;
126 };
127 
135 template <unsigned int VRadius, typename TInput = double, typename TOutput = double>
136 class ITK_TEMPLATE_EXPORT BlackmanWindowFunction
137 {
138 public:
139  inline TOutput
140  operator()(const TInput & A) const
141  {
142  return static_cast<TOutput>(0.42 + 0.5 * std::cos(A * m_Factor1) + 0.08 * std::cos(A * m_Factor2));
143  }
145 
146 private:
148  static const double m_Factor1;
149 
151  static const double m_Factor2;
152 };
153 } // namespace Function
154 
265 template <typename TInputImage,
266  unsigned int VRadius,
267  typename TWindowFunction = Function::HammingWindowFunction<VRadius>,
269  class TCoordRep = double>
271 {
272 public:
273  ITK_DISALLOW_COPY_AND_MOVE(WindowedSincInterpolateImageFunction);
275 
279 
282 
285 
287  itkNewMacro(Self);
288 
291 
294 
296  using RealType = typename Superclass::RealType;
297 
299  static constexpr unsigned int ImageDimension = Superclass::ImageDimension;
300 
304 
306  using SizeType = typename Superclass::SizeType;
307 
309  using ImageType = TInputImage;
310 
313 
314  void
315  SetInputImage(const ImageType * image) override;
316 
323  OutputType
324  EvaluateAtContinuousIndex(const ContinuousIndexType & index) const override;
325 
326  SizeType
327  GetRadius() const override
328  {
329  SizeType radius;
330  radius.Fill(VRadius);
331  return radius;
332  }
333 
334 protected:
337  void
338  PrintSelf(std::ostream & os, Indent indent) const override;
339 
340 private:
341  // Internal type alias
343 
344  // Constant to store twice the radius
345  static const unsigned int m_WindowSize;
346 
348  TWindowFunction m_WindowFunction;
349 
352  unsigned int * m_OffsetTable;
353 
355  unsigned int m_OffsetTableSize;
356 
358  unsigned int ** m_WeightOffsetTable;
359 
361  inline double
362  Sinc(double x) const
363  {
364  const double px = itk::Math::pi * x;
365 
366  return (x == 0.0) ? 1.0 : std::sin(px) / px;
367  }
368 };
369 } // namespace itk
370 
371 #ifndef ITK_MANUAL_INSTANTIATION
372 # include "itkWindowedSincInterpolateImageFunction.hxx"
373 #endif
374 
375 #endif // _itkWindowedSincInterpolateImageFunction_h
itk::WindowedSincInterpolateImageFunction::m_WindowFunction
TWindowFunction m_WindowFunction
Definition: itkWindowedSincInterpolateImageFunction.h:348
itk::Function::HammingWindowFunction::operator()
TOutput operator()(const TInput &A) const
Definition: itkWindowedSincInterpolateImageFunction.h:64
itk::InterpolateImageFunction::SizeType
typename InputImageType::SizeType SizeType
Definition: itkInterpolateImageFunction.h:79
itk::WindowedSincInterpolateImageFunction::EvaluateAtContinuousIndex
OutputType EvaluateAtContinuousIndex(const ContinuousIndexType &index) const override
itkConstNeighborhoodIterator.h
itk::WindowedSincInterpolateImageFunction::ImageDimension
static constexpr unsigned int ImageDimension
Definition: itkWindowedSincInterpolateImageFunction.h:299
itk::WindowedSincInterpolateImageFunction::SetInputImage
void SetInputImage(const ImageType *image) override
itk::Function::CosineWindowFunction::operator()
TOutput operator()(const TInput &A) const
Definition: itkWindowedSincInterpolateImageFunction.h:41
itk::Function::WelchWindowFunction::operator()
TOutput operator()(const TInput &A) const
Definition: itkWindowedSincInterpolateImageFunction.h:87
itk::WindowedSincInterpolateImageFunction::m_OffsetTableSize
unsigned int m_OffsetTableSize
Definition: itkWindowedSincInterpolateImageFunction.h:355
itk::WindowedSincInterpolateImageFunction::ImageType
TInputImage ImageType
Definition: itkWindowedSincInterpolateImageFunction.h:309
itk::Function::BlackmanWindowFunction::m_Factor2
static const double m_Factor2
Definition: itkWindowedSincInterpolateImageFunction.h:151
itk::WindowedSincInterpolateImageFunction::~WindowedSincInterpolateImageFunction
~WindowedSincInterpolateImageFunction() override
itk::Function::CosineWindowFunction
Window function for sinc interpolation.
Definition: itkWindowedSincInterpolateImageFunction.h:37
itk::SmartPointer< Self >
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::WindowedSincInterpolateImageFunction::m_WeightOffsetTable
unsigned int ** m_WeightOffsetTable
Definition: itkWindowedSincInterpolateImageFunction.h:358
itk::Function::CosineWindowFunction::m_Factor
static const double m_Factor
Definition: itkWindowedSincInterpolateImageFunction.h:49
itk::LightObject
Light weight base class for most itk classes.
Definition: itkLightObject.h:59
itk::WindowedSincInterpolateImageFunction::PrintSelf
void PrintSelf(std::ostream &os, Indent indent) const override
itk::Function::WelchWindowFunction::m_Factor
static const double m_Factor
Definition: itkWindowedSincInterpolateImageFunction.h:95
itk::WindowedSincInterpolateImageFunction::m_OffsetTable
unsigned int * m_OffsetTable
Definition: itkWindowedSincInterpolateImageFunction.h:352
itk::InterpolateImageFunction::IndexValueType
typename Superclass::IndexValueType IndexValueType
Definition: itkInterpolateImageFunction.h:76
itk::InterpolateImageFunction::ContinuousIndexType
typename Superclass::ContinuousIndexType ContinuousIndexType
Definition: itkInterpolateImageFunction.h:82
itk::Function::BlackmanWindowFunction::operator()
TOutput operator()(const TInput &A) const
Definition: itkWindowedSincInterpolateImageFunction.h:140
itk::WindowedSincInterpolateImageFunction::Sinc
double Sinc(double x) const
Definition: itkWindowedSincInterpolateImageFunction.h:362
itk::Function::HammingWindowFunction::m_Factor
static const double m_Factor
Definition: itkWindowedSincInterpolateImageFunction.h:72
itk::InterpolateImageFunction::InputImageType
typename Superclass::InputImageType InputImageType
Definition: itkInterpolateImageFunction.h:66
itk::Function::LanczosWindowFunction
Window function for sinc interpolation.
Definition: itkWindowedSincInterpolateImageFunction.h:108
itk::WindowedSincInterpolateImageFunction
Use the windowed sinc function to interpolate.
Definition: itkWindowedSincInterpolateImageFunction.h:270
itk::InterpolateImageFunction::ImageDimension
static constexpr unsigned int ImageDimension
Definition: itkInterpolateImageFunction.h:69
itk::Function::LanczosWindowFunction::operator()
TOutput operator()(const TInput &A) const
Definition: itkWindowedSincInterpolateImageFunction.h:112
itk::InterpolateImageFunction::OutputType
typename Superclass::OutputType OutputType
Definition: itkInterpolateImageFunction.h:63
itk::WindowedSincInterpolateImageFunction::WindowedSincInterpolateImageFunction
WindowedSincInterpolateImageFunction()
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnnulusOperator.h:24
itk::ConstNeighborhoodIterator
Const version of NeighborhoodIterator, defining iteration of a local N-dimensional neighborhood of pi...
Definition: itkConstNeighborhoodIterator.h:51
itk::WindowedSincInterpolateImageFunction::GetRadius
SizeType GetRadius() const override
Definition: itkWindowedSincInterpolateImageFunction.h:327
itk::WindowedSincInterpolateImageFunction::m_WindowSize
static const unsigned int m_WindowSize
Definition: itkWindowedSincInterpolateImageFunction.h:345
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:64
itk::Function::BlackmanWindowFunction::m_Factor1
static const double m_Factor1
Definition: itkWindowedSincInterpolateImageFunction.h:148
itk::Function::BlackmanWindowFunction
Window function for sinc interpolation.
Definition: itkWindowedSincInterpolateImageFunction.h:136
itkZeroFluxNeumannBoundaryCondition.h
itk::Function::LanczosWindowFunction::m_Factor
static const double m_Factor
Definition: itkWindowedSincInterpolateImageFunction.h:125
itk::InterpolateImageFunction
Base class for all image interpolators.
Definition: itkInterpolateImageFunction.h:45
itk::InterpolateImageFunction::IndexType
typename Superclass::IndexType IndexType
Definition: itkInterpolateImageFunction.h:75
itk::Function::HammingWindowFunction
Window function for sinc interpolation.
Definition: itkWindowedSincInterpolateImageFunction.h:60
itk::Function::WelchWindowFunction
Window function for sinc interpolation.
Definition: itkWindowedSincInterpolateImageFunction.h:83
itk::InterpolateImageFunction::RealType
typename NumericTraits< typename TInputImage::PixelType >::RealType RealType
Definition: itkInterpolateImageFunction.h:85