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

itkBSplineScatteredDataPointSetToImageFilter.h

Go to the documentation of this file.
00001 /*=========================================================================
00002 
00003   Program:   Insight Segmentation & Registration Toolkit
00004   Module:    $RCSfile: itkBSplineScatteredDataPointSetToImageFilter.h,v $
00005   Language:  C++
00006   Date:      $Date: 2009-04-20 14:53:42 $
00007   Version:   $Revision: 1.8 $
00008 
00009   Copyright (c) Insight Software Consortium. All rights reserved.
00010   See ITKCopyright.txt or http://www.itk.org/HTML/Copyright.htm for details.
00011 
00012      This software is distributed WITHOUT ANY WARRANTY; without even
00013      the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
00014      PURPOSE.  See the above copyright notices for more information.
00015 
00016 =========================================================================*/
00017 #ifndef __itkBSplineScatteredDataPointSetToImageFilter_h
00018 #define __itkBSplineScatteredDataPointSetToImageFilter_h
00019 
00020 #include "itkPointSetToImageFilter.h"
00021 #include "itkBSplineKernelFunction.h"
00022 #include "itkCoxDeBoorBSplineKernelFunction.h"
00023 #include "itkFixedArray.h"
00024 #include "itkVariableSizeMatrix.h"
00025 #include "itkVector.h"
00026 #include "itkVectorContainer.h"
00027 
00028 #include "vnl/vnl_matrix.h"
00029 
00030 namespace itk
00031 {
00032 
00072 template <class TInputPointSet, class TOutputImage>
00073 class BSplineScatteredDataPointSetToImageFilter
00074 : public PointSetToImageFilter<TInputPointSet, TOutputImage>
00075 {
00076 public:
00077   typedef BSplineScatteredDataPointSetToImageFilter  Self;
00078   typedef PointSetToImageFilter
00079     <TInputPointSet, TOutputImage>                   Superclass;
00080   typedef SmartPointer<Self>                         Pointer;
00081   typedef SmartPointer<const Self>                   ConstPointer;
00082 
00084   itkNewMacro(Self);
00085 
00087   itkStaticConstMacro( ImageDimension, unsigned int,
00088                        TOutputImage::ImageDimension );
00089 
00090   typedef TOutputImage                               ImageType;
00091   typedef TInputPointSet                             PointSetType;
00092 
00094   typedef typename ImageType::PixelType              PixelType;
00095   typedef typename ImageType::RegionType             RegionType;
00096   typedef typename ImageType::SizeType               SizeType;
00097   typedef typename ImageType::IndexType              IndexType;
00098   typedef typename ImageType::PointType              ContinuousIndexType;
00099 
00101   typedef typename PointSetType::PointType           PointType;
00102   typedef typename PointSetType::PixelType           PointDataType;
00103   typedef typename PointSetType::PointDataContainer  PointDataContainerType;
00104 
00106   typedef float                                      RealType;
00107   typedef VectorContainer<unsigned, RealType>        WeightsContainerType;
00108   typedef Image<PointDataType,
00109     itkGetStaticConstMacro( ImageDimension )>        PointDataImageType;
00110   typedef Image<RealType,
00111     itkGetStaticConstMacro( ImageDimension )>        RealImageType;
00112   typedef FixedArray<unsigned,
00113     itkGetStaticConstMacro( ImageDimension )>        ArrayType;
00114   typedef VariableSizeMatrix<RealType>               GradientType;
00116 
00120   typedef CoxDeBoorBSplineKernelFunction<3>          KernelType;
00121   typedef BSplineKernelFunction<0>                   KernelOrder0Type;
00122   typedef BSplineKernelFunction<1>                   KernelOrder1Type;
00123   typedef BSplineKernelFunction<2>                   KernelOrder2Type;
00124   typedef BSplineKernelFunction<3>                   KernelOrder3Type;
00125 
00128   void SetNumberOfLevels( unsigned int );
00129   void SetNumberOfLevels( ArrayType );
00130   itkGetConstReferenceMacro( NumberOfLevels, ArrayType );
00131 
00132   void SetSplineOrder( unsigned int );
00133   void SetSplineOrder( ArrayType );
00134   itkGetConstReferenceMacro( SplineOrder, ArrayType );
00135 
00136   itkSetMacro( NumberOfControlPoints, ArrayType );
00137   itkGetConstReferenceMacro( NumberOfControlPoints, ArrayType );
00138   itkGetConstReferenceMacro( CurrentNumberOfControlPoints, ArrayType );
00139 
00140   itkSetMacro( CloseDimension, ArrayType );
00141   itkGetConstReferenceMacro( CloseDimension, ArrayType );
00142 
00143   itkSetMacro( GenerateOutputImage, bool );
00144   itkGetConstReferenceMacro( GenerateOutputImage, bool );
00145   itkBooleanMacro( GenerateOutputImage );
00146 
00147   void SetPointWeights( WeightsContainerType * weights );
00148 
00152   itkSetMacro( PhiLattice, typename PointDataImageType::Pointer );
00153   itkGetConstMacro( PhiLattice, typename PointDataImageType::Pointer );
00155 
00160   void EvaluateAtPoint( PointType, PointDataType & );
00161   void EvaluateAtIndex( IndexType, PointDataType & );
00162   void EvaluateAtContinuousIndex( ContinuousIndexType, PointDataType & );
00164 
00170   void Evaluate( PointType, PointDataType & );
00171 
00176   void EvaluateGradientAtPoint( PointType, GradientType & );
00177   void EvaluateGradientAtIndex( IndexType, GradientType & );
00178   void EvaluateGradientAtContinuousIndex( ContinuousIndexType, GradientType & );
00180 
00187   void EvaluateGradient( PointType, GradientType & );
00188 
00189 protected:
00190   BSplineScatteredDataPointSetToImageFilter();
00191   virtual ~BSplineScatteredDataPointSetToImageFilter();
00192   void PrintSelf( std::ostream& os, Indent indent ) const;
00193 
00194   void GenerateData();
00195 
00196 private:
00197 
00198   //purposely not implemented
00199   BSplineScatteredDataPointSetToImageFilter(const Self&);
00200   void operator=(const Self&);
00201 
00202   void GenerateControlLattice();
00203   void RefineControlLattice();
00204   void UpdatePointSet();
00205   void GenerateOutputImage();
00206   void GenerateOutputImageFast();
00207   void CollapsePhiLattice( PointDataImageType *,
00208                            PointDataImageType *,
00209                            RealType, unsigned int );
00210 
00211 
00212   bool                                           m_DoMultilevel;
00213   bool                                           m_GenerateOutputImage;
00214   bool                                           m_UsePointWeights;
00215   unsigned int                                   m_MaximumNumberOfLevels;
00216   unsigned int                                   m_CurrentLevel;
00217   ArrayType                                      m_NumberOfControlPoints;
00218   ArrayType                                      m_CurrentNumberOfControlPoints;
00219   ArrayType                                      m_CloseDimension;
00220   ArrayType                                      m_SplineOrder;
00221   ArrayType                                      m_NumberOfLevels;
00222   typename WeightsContainerType::Pointer         m_PointWeights;
00223   typename PointDataImageType::Pointer           m_PhiLattice;
00224   typename PointDataImageType::Pointer           m_PsiLattice;
00225   typename PointDataContainerType::Pointer       m_InputPointData;
00226   typename PointDataContainerType::Pointer       m_OutputPointData;
00227 
00228   typename KernelType::Pointer                   m_Kernel[ImageDimension];
00229   typename KernelOrder0Type::Pointer             m_KernelOrder0;
00230   typename KernelOrder1Type::Pointer             m_KernelOrder1;
00231   typename KernelOrder2Type::Pointer             m_KernelOrder2;
00232   typename KernelOrder3Type::Pointer             m_KernelOrder3;
00233 
00234   RealType                                       m_BSplineEpsilon;
00235 
00236   vnl_matrix<RealType>
00237     m_RefinedLatticeCoefficients[ImageDimension];
00238 
00239   inline typename RealImageType::IndexType
00240   NumberToIndex( unsigned int number, typename RealImageType::SizeType size )
00241     {
00242     typename RealImageType::IndexType k;
00243     k[0] = 1;
00244 
00245     for ( unsigned int i = 1; i < ImageDimension; i++ )
00246       {
00247       k[i] = size[ImageDimension-i-1]*k[i-1];
00248       }
00249     typename RealImageType::IndexType index;
00250     for ( unsigned int i = 0; i < ImageDimension; i++ )
00251       {
00252       index[ImageDimension-i-1]
00253         = static_cast<unsigned int>( number/k[ImageDimension-i-1] );
00254       number %= k[ImageDimension-i-1];
00255       }
00256     return index;
00257     }
00258 };
00259 
00260 } // end namespace itk
00261 
00262 #ifndef ITK_MANUAL_INSTANTIATION
00263 #include "itkBSplineScatteredDataPointSetToImageFilter.txx"
00264 #endif
00265 
00266 #endif
00267 

Generated at Tue Sep 15 02:28:55 2009 for ITK by doxygen 1.5.8 written by Dimitri van Heesch, © 1997-2000