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

itkLevelSetFunction.h

Go to the documentation of this file.
00001 /*=========================================================================
00002 
00003   Program:   Insight Segmentation & Registration Toolkit
00004   Module:    $RCSfile: itkLevelSetFunction.h,v $
00005   Language:  C++
00006   Date:      $Date: 2008-03-03 13:58:51 $
00007   Version:   $Revision: 1.24 $
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 __itkLevelSetFunction_h_
00018 #define __itkLevelSetFunction_h_
00019 
00020 #include "itkFiniteDifferenceFunction.h"
00021 #include "vnl/vnl_matrix_fixed.h"
00022 
00023 namespace itk {
00024 
00063 template <class TImageType>
00064 class ITK_EXPORT LevelSetFunction
00065   : public FiniteDifferenceFunction<TImageType>
00066 {
00067 public:
00069   typedef LevelSetFunction Self;
00070   typedef FiniteDifferenceFunction<TImageType> Superclass;
00071   typedef SmartPointer<Self> Pointer;
00072   typedef SmartPointer<const Self> ConstPointer;
00073 
00075   itkNewMacro(Self);
00076 
00078   itkTypeMacro( LevelSetFunction, FiniteDifferenceFunction );
00079 
00081   itkStaticConstMacro(ImageDimension, unsigned int,Superclass::ImageDimension);
00082 
00084   typedef double TimeStepType;
00085   typedef typename Superclass::ImageType  ImageType;
00086   typedef typename Superclass::PixelType  PixelType;
00087   typedef                      PixelType  ScalarValueType;
00088   typedef typename Superclass::PixelRealType  PixelRealType;
00089   typedef typename Superclass::RadiusType RadiusType;
00090   typedef typename Superclass::NeighborhoodType NeighborhoodType;
00091   typedef typename Superclass::NeighborhoodScalesType NeighborhoodScalesType;
00092   typedef typename Superclass::FloatOffsetType FloatOffsetType;
00093 
00095   //  typedef
00096   typedef FixedArray<ScalarValueType, itkGetStaticConstMacro(ImageDimension)> VectorType;
00097 
00103   struct GlobalDataStruct
00104   {
00105     ScalarValueType m_MaxAdvectionChange;
00106     ScalarValueType m_MaxPropagationChange;
00107     ScalarValueType m_MaxCurvatureChange;
00108 
00110     vnl_matrix_fixed<ScalarValueType,
00111                      itkGetStaticConstMacro(ImageDimension),
00112                      itkGetStaticConstMacro(ImageDimension)> m_dxy;
00113 
00115     ScalarValueType m_dx[itkGetStaticConstMacro(ImageDimension)];
00116 
00117     ScalarValueType m_dx_forward[itkGetStaticConstMacro(ImageDimension)];
00118     ScalarValueType m_dx_backward[itkGetStaticConstMacro(ImageDimension)];
00119     
00120     ScalarValueType m_GradMagSqr;
00121   };
00122 
00124   virtual VectorType AdvectionField(const NeighborhoodType &,
00125                                     const FloatOffsetType &, GlobalDataStruct * = 0)  const
00126     { return m_ZeroVectorConstant; }
00127 
00130   virtual ScalarValueType PropagationSpeed(
00131     const NeighborhoodType& ,
00132     const FloatOffsetType &, GlobalDataStruct * = 0 ) const
00133     { return NumericTraits<ScalarValueType>::Zero; }
00134 
00137   virtual ScalarValueType CurvatureSpeed(const NeighborhoodType &,
00138                                          const FloatOffsetType &, GlobalDataStruct * = 0
00139                                          ) const
00140     { return NumericTraits<ScalarValueType>::One; }
00141 
00144   virtual ScalarValueType LaplacianSmoothingSpeed(
00145     const NeighborhoodType &,
00146     const FloatOffsetType &, GlobalDataStruct * = 0) const
00147     { return NumericTraits<ScalarValueType>::One; }
00148 
00150   virtual void SetAdvectionWeight(const ScalarValueType a)
00151     { m_AdvectionWeight = a; }
00152   ScalarValueType GetAdvectionWeight() const
00153     { return m_AdvectionWeight; }
00155 
00157   virtual void SetPropagationWeight(const ScalarValueType p)
00158     { m_PropagationWeight = p; }
00159   ScalarValueType GetPropagationWeight() const
00160     { return m_PropagationWeight; }
00162 
00164   virtual void SetCurvatureWeight(const ScalarValueType c)
00165     { m_CurvatureWeight = c; }
00166   ScalarValueType GetCurvatureWeight() const
00167     { return m_CurvatureWeight; }
00169 
00171   void SetLaplacianSmoothingWeight(const ScalarValueType c)
00172     { m_LaplacianSmoothingWeight = c; }
00173   ScalarValueType GetLaplacianSmoothingWeight() const
00174     { return m_LaplacianSmoothingWeight; }
00176 
00178   void SetEpsilonMagnitude(const ScalarValueType e)
00179     { m_EpsilonMagnitude = e; }
00180   ScalarValueType GetEpsilonMagnitude() const
00181     { return m_EpsilonMagnitude; }
00183 
00185   virtual PixelType ComputeUpdate(const NeighborhoodType &neighborhood,
00186                                   void *globalData,
00187                                   const FloatOffsetType& = FloatOffsetType(0.0));
00188 
00195   virtual TimeStepType ComputeGlobalTimeStep(void *GlobalData) const;
00196 
00204   virtual void *GetGlobalDataPointer() const
00205     {
00206       GlobalDataStruct *ans = new GlobalDataStruct();
00207       ans->m_MaxAdvectionChange   = NumericTraits<ScalarValueType>::Zero;
00208       ans->m_MaxPropagationChange = NumericTraits<ScalarValueType>::Zero;
00209       ans->m_MaxCurvatureChange   = NumericTraits<ScalarValueType>::Zero;
00210       return ans; 
00211     }
00213 
00217   virtual void Initialize(const RadiusType &r);
00218 
00223   virtual void ReleaseGlobalDataPointer(void *GlobalData) const
00224     { delete (GlobalDataStruct *) GlobalData; }
00225 
00227   virtual ScalarValueType ComputeCurvatureTerm(const NeighborhoodType &,
00228                                                const FloatOffsetType &,
00229                                                GlobalDataStruct *gd = 0
00230                                                );
00231   virtual ScalarValueType ComputeMeanCurvature(const NeighborhoodType &,
00232                                                const FloatOffsetType &,
00233                                                GlobalDataStruct *gd = 0
00234                                                );
00236 
00237   virtual ScalarValueType ComputeMinimalCurvature(const NeighborhoodType &,
00238                                                   const FloatOffsetType &,
00239                                                   GlobalDataStruct *gd = 0
00240                                                   );
00241   
00242   virtual ScalarValueType Compute3DMinimalCurvature(const NeighborhoodType &,
00243                                                     const FloatOffsetType &,
00244                                                     GlobalDataStruct *gd = 0
00245                                                     );
00246   
00248   void SetUseMinimalCurvature( bool b )
00249   {
00250     m_UseMinimalCurvature = b;
00251   }
00252   bool GetUseMinimalCurvature() const
00253   {
00254     return m_UseMinimalCurvature;
00255   }
00256   void UseMinimalCurvatureOn()
00257   {
00258     this->SetUseMinimalCurvature(true);
00259   }
00260   void UseMinimalCurvatureOff()
00261   {
00262     this->SetUseMinimalCurvature(false);
00263   }
00265 
00270   static void SetMaximumCurvatureTimeStep(double n)
00271   {
00272     m_DT = n;
00273   }
00274   static double GetMaximumCurvatureTimeStep()
00275   {
00276     return m_DT;
00277   }
00279 
00284   static void SetMaximumPropagationTimeStep(double n)
00285   {
00286     m_WaveDT = n;
00287   }
00288   static double GetMaximumPropagationTimeStep()
00289   {
00290     return m_WaveDT;
00291   }
00293 
00294 protected:
00295   LevelSetFunction()
00296   {
00297     m_EpsilonMagnitude = 1.0e-5;
00298     m_AdvectionWeight = m_PropagationWeight 
00299       = m_CurvatureWeight = m_LaplacianSmoothingWeight 
00300       = NumericTraits<ScalarValueType>::Zero;
00301     m_UseMinimalCurvature = false;
00302   }
00303   virtual ~LevelSetFunction() {}
00304   void PrintSelf(std::ostream &s, Indent indent) const;
00305   
00307   static double m_WaveDT;
00308   static double m_DT;
00309 
00311   std::slice x_slice[itkGetStaticConstMacro(ImageDimension)];
00312 
00314   ::size_t m_Center;
00315 
00317   ::size_t m_xStride[itkGetStaticConstMacro(ImageDimension)];
00318 
00319   bool m_UseMinimalCurvature;
00320 
00323   static VectorType InitializeZeroVectorConstant();
00324 
00326   static VectorType m_ZeroVectorConstant;
00327 
00329   ScalarValueType m_EpsilonMagnitude;
00330 
00332   ScalarValueType m_AdvectionWeight;
00333 
00335   ScalarValueType m_PropagationWeight;
00336 
00338   ScalarValueType m_CurvatureWeight;
00339 
00341   ScalarValueType m_LaplacianSmoothingWeight;
00342 
00343 private:
00344   LevelSetFunction(const Self&); //purposely not implemented
00345   void operator=(const Self&);   //purposely not implemented
00346 };
00347 
00348 } // namespace itk
00349 
00350 // Define instantiation macro for this template.
00351 #define ITK_TEMPLATE_LevelSetFunction(_, EXPORT, x, y) namespace itk { \
00352   _(1(class EXPORT LevelSetFunction< ITK_TEMPLATE_1 x >)) \
00353   namespace Templates { typedef LevelSetFunction< ITK_TEMPLATE_1 x > \
00354                                      LevelSetFunction##y; } \
00355   }
00356 
00357 #if ITK_TEMPLATE_EXPLICIT
00358 # include "Templates/itkLevelSetFunction+-.h"
00359 #endif
00360 
00361 #if ITK_TEMPLATE_TXX
00362 # include "itkLevelSetFunction.txx"
00363 #endif
00364 
00365 #endif
00366 
00367 

Generated at Tue Jul 29 21:04:32 2008 for ITK by doxygen 1.5.1 written by Dimitri van Heesch, © 1997-2000