Main Page   Groups   Namespace List   Class Hierarchy   Alphabetical List   Compound List   File List   Namespace Members   Compound Members   File Members   Concepts

itkWindowedSincInterpolateImageFunction.h

Go to the documentation of this file.
00001 #ifndef _itkWindowedSincInterpolateImageFunction_h
00002 #define _itkWindowedSincInterpolateImageFunction_h
00003 
00004 #include "itkConstNeighborhoodIterator.h"
00005 #include "itkConstantBoundaryCondition.h"
00006 #include "itkInterpolateImageFunction.h"
00007 
00008 namespace itk
00009 {
00010 
00011 namespace Function {
00012 
00019 template< unsigned int VRadius, 
00020   class TInput=double, class TOutput=double>
00021 class CosineWindowFunction
00022 {
00023 public:
00024   inline TOutput operator()( const TInput & A ) const
00025     { return (TOutput) cos( A * m_Factor ); }
00026 private:
00028   static const double m_Factor;
00029 }; 
00030 
00037 template< unsigned int VRadius, 
00038   class TInput=double, class TOutput=double>
00039 class HammingWindowFunction
00040 {
00041 public:
00042   inline TOutput operator()( const TInput & A ) const
00043    { return (TOutput) 0.54 + 0.46 * cos( A * m_Factor ); }
00044 private:
00046   static const double m_Factor;
00047 }; 
00048 
00055 template< unsigned int VRadius, 
00056   class TInput=double, class TOutput=double>
00057 class WelchWindowFunction
00058 {
00059 public:
00060   inline TOutput operator()( const TInput & A ) const
00061    { return (TOutput) (1.0 - A * m_Factor * A); }
00062 private:
00064   static const double m_Factor;
00065 }; 
00066 
00075 template< unsigned int VRadius, 
00076   class TInput=double, class TOutput=double>
00077 class LancosWindowFunction
00078 {
00079 public:
00080   inline TOutput operator()( const TInput & A ) const
00081     {
00082     if(A == 0.0) return (TOutput) 1.0;
00083     double z = m_Factor * A;
00084     return (TOutput) ( sin(z) / z ); 
00085     }
00086 private:
00088   static const double m_Factor;
00089 }; 
00090 
00097 template< unsigned int VRadius, 
00098   class TInput=double, class TOutput=double>
00099 class BlackmanWindowFunction
00100 {
00101 public:
00102   inline TOutput operator()( const TInput & A ) const
00103     {
00104     return (TOutput)
00105       (0.42 + 0.5 * cos(A * m_Factor1) + 0.08 * cos(A * m_Factor2));
00106     }
00107 private:
00109   static const double m_Factor1;
00110   
00112   static const double m_Factor2;
00113 }; 
00114 
00115 } // namespace Function
00116 
00226 template <
00227   class TInputImage, 
00228   unsigned int VRadius, 
00229   class TWindowFunction = Function::HammingWindowFunction<VRadius>,
00230   class TBoundaryCondition = ConstantBoundaryCondition<TInputImage>,
00231   class TCoordRep=double >
00232 class ITK_EXPORT WindowedSincInterpolateImageFunction : 
00233   public InterpolateImageFunction<TInputImage, TCoordRep>
00234 {
00235 public:
00237   typedef WindowedSincInterpolateImageFunction Self;
00238   typedef InterpolateImageFunction<TInputImage,TCoordRep> Superclass;
00239   typedef SmartPointer<Self> Pointer;
00240   typedef SmartPointer<const Self>  ConstPointer;
00241   
00243   itkTypeMacro(WindowedSincInterpolateImageFunction, 
00244     InterpolateImageFunction);
00245 
00247   itkNewMacro(Self);  
00248 
00250   typedef typename Superclass::OutputType OutputType;
00251 
00253   typedef typename Superclass::InputImageType InputImageType;
00254 
00256   typedef typename Superclass::RealType RealType;
00257 
00259   itkStaticConstMacro(ImageDimension, unsigned int,Superclass::ImageDimension);
00260 
00262   typedef typename Superclass::IndexType IndexType;
00263   
00265   typedef TInputImage ImageType;
00266 
00268   typedef typename Superclass::ContinuousIndexType ContinuousIndexType;
00269 
00270   virtual void SetInputImage(const ImageType *image);
00271 
00278   virtual OutputType EvaluateAtContinuousIndex( 
00279     const ContinuousIndexType & index ) const;
00280 
00281 protected:
00282   WindowedSincInterpolateImageFunction();
00283   virtual ~WindowedSincInterpolateImageFunction();
00284   void PrintSelf(std::ostream& os, Indent indent) const;
00285 
00286 private:
00287   WindowedSincInterpolateImageFunction( const Self& ); //purposely not implemented
00288   void operator=( const Self& ); //purposely not implemented
00289 
00290   // Internal typedefs
00291   typedef ConstNeighborhoodIterator<
00292     ImageType, TBoundaryCondition> IteratorType;
00293 
00294   // Constant to store twice the radius
00295   static const unsigned int m_WindowSize;
00296 
00298   TWindowFunction m_WindowFunction;
00299   
00302   unsigned int *m_OffsetTable;
00303 
00305   unsigned int m_OffsetTableSize;
00306 
00308   unsigned int **m_WeightOffsetTable;
00309 
00311   inline double Sinc(double x) const
00312     { 
00313     double px = vnl_math::pi * x;
00314     return (x == 0.0) ? 1.0 : sin(px) / px;
00315     }
00316 };
00317 
00318 } // namespace itk
00319 
00320 #ifndef ITK_MANUAL_INSTANTIATION
00321 #include "itkWindowedSincInterpolateImageFunction.txx"
00322 #endif
00323 
00324 #endif // _itkWindowedSincInterpolateImageFunction_h

Generated at Thu May 25 00:18:47 2006 for ITK by doxygen 1.3.5 written by Dimitri van Heesch, © 1997-2000