ITK  5.0.0
Insight Segmentation and Registration Toolkit
itkWindowedSincInterpolateImageFunction.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright Insight Software Consortium
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,
37  typename TInput = double, typename TOutput = double >
38 class ITK_TEMPLATE_EXPORT CosineWindowFunction
39 {
40 public:
41  inline TOutput operator()(const TInput & A) const
42  { return static_cast<TOutput>(std::cos(A * m_Factor)); }
44 
45 private:
47  static const double m_Factor;
48 };
49 
57 template< unsigned int VRadius,
58  typename TInput = double, typename TOutput = double >
59 class ITK_TEMPLATE_EXPORT HammingWindowFunction
60 {
61 public:
62  inline TOutput operator()(const TInput & A) const
63  { return static_cast<TOutput>(0.54 + 0.46 * std::cos(A * m_Factor) ); }
65 
66 private:
68  static const double m_Factor;
69 };
70 
78 template< unsigned int VRadius,
79  typename TInput = double, typename TOutput = double >
80 class ITK_TEMPLATE_EXPORT WelchWindowFunction
81 {
82 public:
83  inline TOutput operator()(const TInput & A) const
84  { return static_cast<TOutput>( 1.0 - A * m_Factor * A ); }
86 
87 private:
89  static const double m_Factor;
90 };
91 
101 template< unsigned int VRadius,
102  typename TInput = double, typename TOutput = double >
103 class ITK_TEMPLATE_EXPORT LanczosWindowFunction
104 {
105 public:
106  inline TOutput operator()(const TInput & A) const
107  {
108  if ( A == 0.0 ) { return static_cast<TOutput>(1.0); }
109  double z = m_Factor * A;
110  return static_cast<TOutput>( std::sin(z) / z );
111  }
113 
114 private:
116  static const double m_Factor;
117 };
118 
126 template< unsigned int VRadius,
127  typename TInput = double, typename TOutput = double >
128 class ITK_TEMPLATE_EXPORT BlackmanWindowFunction
129 {
130 public:
131  inline TOutput operator()(const TInput & A) const
132  {
133  return static_cast<TOutput>
134  ( 0.42 + 0.5 * std::cos(A * m_Factor1) + 0.08 * std::cos(A * m_Factor2) );
135  }
137 
138 private:
140  static const double m_Factor1;
141 
143  static const double m_Factor2;
144 };
145 } // namespace Function
146 
257 template<
258  typename TInputImage,
259  unsigned int VRadius,
260  typename TWindowFunction = Function::HammingWindowFunction< VRadius >,
262  class TCoordRep = double >
264  public InterpolateImageFunction< TInputImage, TCoordRep >
265 {
266 public:
267  ITK_DISALLOW_COPY_AND_ASSIGN(WindowedSincInterpolateImageFunction);
269 
273 
276 
280 
282  itkNewMacro(Self);
283 
286 
289 
291  using RealType = typename Superclass::RealType;
292 
294  static constexpr unsigned int ImageDimension = Superclass::ImageDimension;
295 
299 
301  using ImageType = TInputImage;
302 
305 
306  void SetInputImage(const ImageType *image) override;
307 
315  const ContinuousIndexType & index) const override;
316 
317 protected:
320  void PrintSelf(std::ostream & os, Indent indent) const override;
321 
322 private:
323  // Internal type alias
325  ImageType, TBoundaryCondition >;
326 
327  // Constant to store twice the radius
328  static const unsigned int m_WindowSize;
329 
331  TWindowFunction m_WindowFunction;
332 
335  unsigned int *m_OffsetTable;
336 
338  unsigned int m_OffsetTableSize;
339 
341  unsigned int **m_WeightOffsetTable;
342 
344  inline double Sinc(double x) const
345  {
346  double px = itk::Math::pi * x;
347 
348  return ( x == 0.0 ) ? 1.0 : std::sin(px) / px;
349  }
350 };
351 } // namespace itk
352 
353 #ifndef ITK_MANUAL_INSTANTIATION
354 #include "itkWindowedSincInterpolateImageFunction.hxx"
355 #endif
356 
357 #endif // _itkWindowedSincInterpolateImageFunction_h
A function object that determines a neighborhood of values at an image boundary according to a Neuman...
Light weight base class for most itk classes.
typename Superclass::ContinuousIndexType ContinuousIndexType
Use the windowed sinc function to interpolate.
typename Superclass::OutputType OutputType
static constexpr double pi
Definition: itkMath.h:63
OutputType EvaluateAtContinuousIndex(const ContinuousIndexType &index) const override
Const version of NeighborhoodIterator, defining iteration of a local N-dimensional neighborhood of pi...
typename NumericTraits< typename TInputImage::PixelType >::RealType RealType
static constexpr unsigned int ImageDimension
typename Superclass::IndexType IndexType
Window function for sinc interpolation. .
Window function for sinc interpolation. .
void SetInputImage(const ImageType *image) override
Base class for all image interpolaters.
typename Superclass::InputImageType InputImageType
Window function for sinc interpolation. Note: Paper referenced in WindowedSincInterpolateImageFuncti...
Control indentation during Print() invocation.
Definition: itkIndent.h:49
Window function for sinc interpolation. .
typename Superclass::IndexValueType IndexValueType
void PrintSelf(std::ostream &os, Indent indent) const override