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

itkOptBSplineInterpolateImageFunction.h

Go to the documentation of this file.
00001 /*=========================================================================
00002 
00003   Program:   Insight Segmentation & Registration Toolkit
00004   Module:    $RCSfile: itkOptBSplineInterpolateImageFunction.h,v $
00005   Language:  C++
00006   Date:      $Date: 2008/02/07 15:07:57 $
00007   Version:   $Revision: 1.6 $
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   Portions of this code are covered under the VTK copyright.
00013   See VTKCopyright.txt or http://www.kitware.com/VTKCopyright.htm for details.
00014 
00015      This software is distributed WITHOUT ANY WARRANTY; without even
00016      the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
00017      PURPOSE.  See the above copyright notices for more information.
00018 
00019 =========================================================================*/
00020 #ifndef __itkOptBSplineInterpolateImageFunction_h
00021 #define __itkOptBSplineInterpolateImageFunction_h
00022 
00023 #include <vector>
00024 
00025 #include "itkImageLinearIteratorWithIndex.h"
00026 #include "itkInterpolateImageFunction.h"
00027 #include "vnl/vnl_matrix.h"
00028 
00029 #include "itkBSplineDecompositionImageFilter.h"
00030 #include "itkConceptChecking.h"
00031 #include "itkCovariantVector.h"
00032 
00033 namespace itk
00034 {
00068 template <
00069   class TImageType,
00070   class TCoordRep = double,
00071   class TCoefficientType = double >
00072 class ITK_EXPORT BSplineInterpolateImageFunction :
00073     public InterpolateImageFunction<TImageType,TCoordRep>
00074 {
00075 public:
00077   typedef BSplineInterpolateImageFunction                   Self;
00078   typedef InterpolateImageFunction<TImageType,TCoordRep>    Superclass;
00079   typedef SmartPointer<Self>                                Pointer;
00080   typedef SmartPointer<const Self>                          ConstPointer;
00081 
00083   itkTypeMacro(BSplineInterpolateImageFunction, InterpolateImageFunction);
00084 
00085 
00087   itkNewMacro( Self );
00088 
00090   typedef typename Superclass::OutputType OutputType;
00091 
00093   typedef typename Superclass::InputImageType InputImageType;
00094 
00096   itkStaticConstMacro(ImageDimension, unsigned int,Superclass::ImageDimension);
00097 
00099   typedef typename Superclass::IndexType IndexType;
00100 
00102   typedef typename Superclass::ContinuousIndexType ContinuousIndexType;
00103 
00105   typedef typename Superclass::PointType PointType;
00106 
00108   typedef ImageLinearIteratorWithIndex<TImageType> Iterator;
00109 
00111   typedef TCoefficientType CoefficientDataType;
00112   typedef Image<CoefficientDataType,
00113                      itkGetStaticConstMacro(ImageDimension) >
00114                                                            CoefficientImageType;
00115 
00117   typedef BSplineDecompositionImageFilter<TImageType,
00118                                                CoefficientImageType>
00119                                                               CoefficientFilter;
00120   typedef typename CoefficientFilter::Pointer CoefficientFilterPointer;
00121 
00123   typedef CovariantVector<OutputType,
00124                           itkGetStaticConstMacro(ImageDimension) >
00125                                                             CovariantVectorType;
00126 
00127 
00136   virtual OutputType Evaluate( const PointType & point ) const
00137     {
00138     ContinuousIndexType index;
00139     this->GetInputImage()->TransformPhysicalPointToContinuousIndex( point,
00140                                                                     index );
00141     // No thread info passed in, so call method that doesn't need thread ID.
00142     return ( this->EvaluateAtContinuousIndex( index ) );
00143     }
00145 
00146   virtual OutputType Evaluate( const PointType & point,
00147                                unsigned int threadID ) const
00148     {
00149     ContinuousIndexType index;
00150     this->GetInputImage()->TransformPhysicalPointToContinuousIndex( point,
00151                                                                     index );
00152     return ( this->EvaluateAtContinuousIndex( index, threadID ) );
00153     }
00154 
00155   virtual OutputType EvaluateAtContinuousIndex( const ContinuousIndexType &
00156                                                                  index ) const
00157     {
00158     // Don't know thread information, make evaluateIndex, weights on the stack.
00159     // Slower, but safer.
00160     vnl_matrix<long>        evaluateIndex(ImageDimension, ( m_SplineOrder + 1 ));
00161     vnl_matrix<double>      weights(ImageDimension, ( m_SplineOrder + 1 ));
00162 
00163     // Pass evaluateIndex, weights by reference. They're only good as long
00164     // as this method is in scope.
00165     return this->EvaluateAtContinuousIndexInternal( index,
00166                                                     evaluateIndex,
00167                                                     weights);
00168     }
00169 
00170   virtual OutputType EvaluateAtContinuousIndex( const ContinuousIndexType &
00171                                                                         index,
00172                                                 unsigned int threadID ) const;
00173 
00174   CovariantVectorType EvaluateDerivative( const PointType & point ) const
00175     {
00176     ContinuousIndexType index;
00177     this->GetInputImage()->TransformPhysicalPointToContinuousIndex( point,
00178                                                                     index );
00179     // No thread info passed in, so call method that doesn't need thread ID.
00180     return ( this->EvaluateDerivativeAtContinuousIndex( index ) );
00181     }
00182 
00183   CovariantVectorType EvaluateDerivative( const PointType & point,
00184                                           unsigned int threadID ) const
00185     {
00186     ContinuousIndexType index;
00187     this->GetInputImage()->TransformPhysicalPointToContinuousIndex( point,
00188                                                                     index );
00189     return ( this->EvaluateDerivativeAtContinuousIndex( index, threadID ) );
00190     }
00191 
00192   CovariantVectorType EvaluateDerivativeAtContinuousIndex(
00193                                          const ContinuousIndexType & x ) const
00194     {
00195     // Don't know thread information, make evaluateIndex, weights, weightsDerivative
00196     // on the stack.
00197     // Slower, but safer.
00198     vnl_matrix<long>          evaluateIndex(ImageDimension, ( m_SplineOrder + 1 ));
00199     vnl_matrix<double>        weights(ImageDimension, ( m_SplineOrder + 1 ));
00200     vnl_matrix<double>        weightsDerivative(ImageDimension, ( m_SplineOrder + 1));
00201 
00202     // Pass evaluateIndex, weights, weightsDerivative by reference. They're only good
00203     // as long as this method is in scope.
00204     return this->EvaluateDerivativeAtContinuousIndexInternal( x,
00205                                                               evaluateIndex,
00206                                                               weights,
00207                                                               weightsDerivative );
00208     }
00209 
00210   CovariantVectorType EvaluateDerivativeAtContinuousIndex(
00211                                          const ContinuousIndexType & x,
00212                                          unsigned int threadID ) const;
00213 
00214   void EvaluateValueAndDerivative( const PointType & point,
00215                                    OutputType & value,
00216                                    CovariantVectorType & deriv ) const
00217     {
00218     ContinuousIndexType index;
00219     this->GetInputImage()->TransformPhysicalPointToContinuousIndex( point,
00220                                                                     index );
00221 
00222     // No thread info passed in, so call method that doesn't need thread ID.
00223     this->EvaluateValueAndDerivativeAtContinuousIndex( index,
00224                                                        value,
00225                                                        deriv );
00226     }
00227 
00228   void EvaluateValueAndDerivative( const PointType & point,
00229                                    OutputType & value,
00230                                    CovariantVectorType & deriv,
00231                                    unsigned int threadID ) const
00232     {
00233     ContinuousIndexType index;
00234     this->GetInputImage()->TransformPhysicalPointToContinuousIndex( point,
00235                                                                     index );
00236     this->EvaluateValueAndDerivativeAtContinuousIndex( index,
00237                                                        value,
00238                                                        deriv,
00239                                                        threadID );
00240     }
00241 
00242   void EvaluateValueAndDerivativeAtContinuousIndex(
00243                                                 const ContinuousIndexType & x,
00244                                                 OutputType & value,
00245                                                 CovariantVectorType & deriv
00246                                                 ) const
00247     {
00248     // Don't know thread information, make evaluateIndex, weights, weightsDerivative
00249     // on the stack.
00250     // Slower, but safer.
00251     vnl_matrix<long>          evaluateIndex(ImageDimension, ( m_SplineOrder + 1 ));
00252     vnl_matrix<double>        weights(ImageDimension, ( m_SplineOrder + 1 ));
00253     vnl_matrix<double>        weightsDerivative(ImageDimension, ( m_SplineOrder + 1));
00254 
00255     // Pass evaluateIndex, weights, weightsDerivative by reference. They're only good
00256     // as long as this method is in scope.
00257     this->EvaluateValueAndDerivativeAtContinuousIndexInternal(x,
00258                                                               value,
00259                                                               deriv,
00260                                                               evaluateIndex,
00261                                                               weights,
00262                                                               weightsDerivative );
00263     }
00264 
00265   void EvaluateValueAndDerivativeAtContinuousIndex(
00266                                                 const ContinuousIndexType & x,
00267                                                 OutputType & value,
00268                                                 CovariantVectorType & deriv,
00269                                                 unsigned int threadID ) const;
00270 
00271 
00274   void SetSplineOrder(unsigned int SplineOrder);
00275   itkGetMacro(SplineOrder, int);
00277 
00278   void SetNumberOfThreads(unsigned int numThreads);
00279   itkGetMacro(NumberOfThreads, int);
00280 
00282   virtual void SetInputImage(const TImageType * inputData);
00283 
00284 
00294   itkSetMacro( UseImageDirection, bool );
00295   itkGetMacro( UseImageDirection, bool );
00296   itkBooleanMacro( UseImageDirection );
00298 
00299 
00300 protected:
00301 
00320   virtual OutputType EvaluateAtContinuousIndexInternal( const ContinuousIndexType & index,
00321                                                         vnl_matrix<long>& evaluateIndex,
00322                                                         vnl_matrix<double>& weights) const;
00323 
00324   virtual void EvaluateValueAndDerivativeAtContinuousIndexInternal( const ContinuousIndexType & x,
00325                                                        OutputType & value,
00326                                                        CovariantVectorType & derivativeValue,
00327                                                        vnl_matrix<long>& evaluateIndex,
00328                                                        vnl_matrix<double>& weights,
00329                                                        vnl_matrix<double>& weightsDerivative
00330                                                        ) const;
00331 
00332   virtual CovariantVectorType EvaluateDerivativeAtContinuousIndexInternal( const ContinuousIndexType & x,
00333                                                                            vnl_matrix<long>& evaluateIndex,
00334                                                                            vnl_matrix<double>& weights,
00335                                                                            vnl_matrix<double>& weightsDerivative
00336                                                                            ) const;
00337 
00338 
00339   BSplineInterpolateImageFunction();
00340   ~BSplineInterpolateImageFunction();
00341   void operator=( const Self& ); //purposely not implemented
00342   void PrintSelf(std::ostream& os, Indent indent) const;
00343 
00344   // These are needed by the smoothing spline routine.
00345   // temp storage for processing of Coefficients
00346   std::vector<CoefficientDataType>    m_Scratch;
00347   // Image size
00348   typename TImageType::SizeType       m_DataLength;
00349   // User specified spline order (3rd or cubic is the default)
00350   unsigned int                        m_SplineOrder;
00351 
00352   // Spline coefficients
00353   typename CoefficientImageType::ConstPointer       m_Coefficients;
00354 
00355 private:
00356   BSplineInterpolateImageFunction( const Self& ); //purposely not implemented
00358   void SetInterpolationWeights( const ContinuousIndexType & x,
00359                                 const vnl_matrix<long> & EvaluateIndex,
00360                                 vnl_matrix<double> & weights,
00361                                 unsigned int splineOrder ) const;
00362 
00364   void SetDerivativeWeights( const ContinuousIndexType & x,
00365                              const vnl_matrix<long> & EvaluateIndex,
00366                              vnl_matrix<double> & weights,
00367                              unsigned int splineOrder ) const;
00368 
00371   void GeneratePointsToIndex(  );
00372 
00374   void DetermineRegionOfSupport( vnl_matrix<long> & evaluateIndex,
00375                                  const ContinuousIndexType & x,
00376                                  unsigned int splineOrder ) const;
00377 
00380   void ApplyMirrorBoundaryConditions(vnl_matrix<long> & evaluateIndex,
00381                                      unsigned int splineOrder) const;
00382 
00383 
00384   Iterator                  m_CIterator;    // Iterator for traversing spline coefficients.
00385   unsigned long             m_MaxNumberInterpolationPoints; // number of neighborhood points used for interpolation
00386   std::vector<IndexType>    m_PointsToIndex;  // Preallocation of interpolation neighborhood indicies
00387 
00388   CoefficientFilterPointer     m_CoefficientFilter;
00389 
00390   // flag to take or not the image direction into account when computing the
00391   // derivatives.
00392   bool m_UseImageDirection;
00393 
00394   unsigned int         m_NumberOfThreads;
00395   vnl_matrix<long>   * m_ThreadedEvaluateIndex;
00396   vnl_matrix<double> * m_ThreadedWeights;
00397   vnl_matrix<double> * m_ThreadedWeightsDerivative;
00398 };
00399 
00400 } // namespace itk
00401 
00402 #ifndef ITK_MANUAL_INSTANTIATION
00403 #include "itkOptBSplineInterpolateImageFunction.txx"
00404 #endif
00405 
00406 #endif
00407 

Generated at Wed Nov 5 23:17:46 2008 for ITK by doxygen 1.5.1 written by Dimitri van Heesch, © 1997-2000