ITK  4.9.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 >
39 {
40 public:
41  inline TOutput operator()(const TInput & A) const
42  { return (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 >
60 {
61 public:
62  inline TOutput operator()(const TInput & A) const
63  { return (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 >
81 {
82 public:
83  inline TOutput operator()(const TInput & A) const
84  { return (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 >
104 {
105 public:
106  inline TOutput operator()(const TInput & A) const
107  {
108  if ( A == 0.0 ) { return (TOutput)1.0; }
109  double z = m_Factor * A;
110  return (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 >
129 {
130 public:
131  inline TOutput operator()(const TInput & A) const
132  {
133  return (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 
271 
274 
278 
280  itkNewMacro(Self);
281 
283  typedef typename Superclass::OutputType OutputType;
284 
287 
289  typedef typename Superclass::RealType RealType;
290 
292  itkStaticConstMacro(ImageDimension, unsigned int, Superclass::ImageDimension);
293 
297 
299  typedef TInputImage ImageType;
300 
303 
304  virtual void SetInputImage(const ImageType *image) ITK_OVERRIDE;
305 
313  const ContinuousIndexType & index) const ITK_OVERRIDE;
314 
315 protected:
318  void PrintSelf(std::ostream & os, Indent indent) const ITK_OVERRIDE;
319 
320 private:
321  WindowedSincInterpolateImageFunction(const Self &); //not implemented
322  void operator=(const Self &) ITK_DELETE_FUNCTION;
323 
324  // Internal typedefs
326  ImageType, TBoundaryCondition > IteratorType;
327 
328  // Constant to store twice the radius
329  static const unsigned int m_WindowSize;
330 
332  TWindowFunction m_WindowFunction;
333 
336  unsigned int *m_OffsetTable;
337 
339  unsigned int m_OffsetTableSize;
340 
342  unsigned int **m_WeightOffsetTable;
343 
345  inline double Sinc(double x) const
346  {
347  double px = vnl_math::pi * x;
348 
349  return ( x == 0.0 ) ? 1.0 : std::sin(px) / px;
350  }
351 };
352 } // namespace itk
353 
354 #ifndef ITK_MANUAL_INSTANTIATION
355 #include "itkWindowedSincInterpolateImageFunction.hxx"
356 #endif
357 
358 #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.
Use the windowed sinc function to interpolate.
InterpolateImageFunction< TInputImage, TCoordRep > Superclass
static const double pi
Definition: itkMath.h:57
Const version of NeighborhoodIterator, defining iteration of a local N-dimensional neighborhood of pi...
virtual OutputType EvaluateAtContinuousIndex(const ContinuousIndexType &index) const override
static const unsigned int ImageDimension
Window function for sinc interpolation. .
void operator=(const Self &) ITK_DELETE_FUNCTION
Window function for sinc interpolation. .
Superclass::InputImageType InputImageType
Superclass::ContinuousIndexType ContinuousIndexType
Base class for all image interpolaters.
virtual void SetInputImage(const ImageType *image) override
Window function for sinc interpolation. Note: Paper referenced in WindowedSincInterpolateImageFuncti...
Control indentation during Print() invocation.
Definition: itkIndent.h:49
Window function for sinc interpolation. .
void PrintSelf(std::ostream &os, Indent indent) const override
Superclass::IndexValueType IndexValueType
NumericTraits< typename TInputImage::PixelType >::RealType RealType