ITK  4.0.0
Insight Segmentation and Registration Toolkit
itkLevelSetFunction.h
Go to the documentation of this file.
00001 /*=========================================================================
00002  *
00003  *  Copyright Insight Software Consortium
00004  *
00005  *  Licensed under the Apache License, Version 2.0 (the "License");
00006  *  you may not use this file except in compliance with the License.
00007  *  You may obtain a copy of the License at
00008  *
00009  *         http://www.apache.org/licenses/LICENSE-2.0.txt
00010  *
00011  *  Unless required by applicable law or agreed to in writing, software
00012  *  distributed under the License is distributed on an "AS IS" BASIS,
00013  *  WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
00014  *  See the License for the specific language governing permissions and
00015  *  limitations under the License.
00016  *
00017  *=========================================================================*/
00018 #ifndef __itkLevelSetFunction_h
00019 #define __itkLevelSetFunction_h
00020 
00021 #include "itkFiniteDifferenceFunction.h"
00022 #include "vnl/vnl_matrix_fixed.h"
00023 
00024 namespace itk
00025 {
00065 template< class TImageType >
00066 class ITK_EXPORT LevelSetFunction:
00067   public FiniteDifferenceFunction< TImageType >
00068 {
00069 public:
00071   typedef LevelSetFunction                       Self;
00072   typedef FiniteDifferenceFunction< TImageType > Superclass;
00073   typedef SmartPointer< Self >                   Pointer;
00074   typedef SmartPointer< const Self >             ConstPointer;
00075 
00077   itkNewMacro(Self);
00078 
00080   itkTypeMacro(LevelSetFunction, FiniteDifferenceFunction);
00081 
00083   itkStaticConstMacro(ImageDimension, unsigned int, Superclass::ImageDimension);
00084 
00086   typedef double                                      TimeStepType;
00087   typedef typename Superclass::ImageType              ImageType;
00088   typedef typename Superclass::PixelType              PixelType;
00089   typedef                      PixelType              ScalarValueType;
00090   typedef typename Superclass::PixelRealType          PixelRealType;
00091   typedef typename Superclass::RadiusType             RadiusType;
00092   typedef typename Superclass::NeighborhoodType       NeighborhoodType;
00093   typedef typename Superclass::NeighborhoodScalesType NeighborhoodScalesType;
00094   typedef typename Superclass::FloatOffsetType        FloatOffsetType;
00095 
00097   //  typedef
00098   typedef FixedArray< ScalarValueType, itkGetStaticConstMacro(ImageDimension) > VectorType;
00099 
00105   struct GlobalDataStruct {
00106     ScalarValueType m_MaxAdvectionChange;
00107     ScalarValueType m_MaxPropagationChange;
00108     ScalarValueType m_MaxCurvatureChange;
00109 
00111     vnl_matrix_fixed< ScalarValueType,
00112                       itkGetStaticConstMacro(ImageDimension),
00113                       itkGetStaticConstMacro(ImageDimension) > m_dxy;
00114 
00116     ScalarValueType m_dx[itkGetStaticConstMacro(ImageDimension)];
00117 
00118     ScalarValueType m_dx_forward[itkGetStaticConstMacro(ImageDimension)];
00119     ScalarValueType m_dx_backward[itkGetStaticConstMacro(ImageDimension)];
00120 
00121     ScalarValueType m_GradMagSqr;
00122   };
00123 
00125   virtual VectorType AdvectionField(const NeighborhoodType &,
00126                                     const FloatOffsetType &, GlobalDataStruct * = 0)  const
00127   { return m_ZeroVectorConstant; }
00128 
00131   virtual ScalarValueType PropagationSpeed(
00132     const NeighborhoodType &,
00133     const FloatOffsetType &, GlobalDataStruct * = 0) const
00134   { return NumericTraits< ScalarValueType >::Zero; }
00135 
00138   virtual ScalarValueType CurvatureSpeed(const NeighborhoodType &,
00139                                          const FloatOffsetType &, GlobalDataStruct * = 0
00140                                          ) const
00141   { return NumericTraits< ScalarValueType >::One; }
00142 
00145   virtual ScalarValueType LaplacianSmoothingSpeed(
00146     const NeighborhoodType &,
00147     const FloatOffsetType &, GlobalDataStruct * = 0) const
00148   { return NumericTraits< ScalarValueType >::One; }
00149 
00151   virtual void SetAdvectionWeight(const ScalarValueType a)
00152   { m_AdvectionWeight = a; }
00153   ScalarValueType GetAdvectionWeight() const
00154   { return m_AdvectionWeight; }
00156 
00158   virtual void SetPropagationWeight(const ScalarValueType p)
00159   { m_PropagationWeight = p; }
00160   ScalarValueType GetPropagationWeight() const
00161   { return m_PropagationWeight; }
00163 
00165   virtual void SetCurvatureWeight(const ScalarValueType c)
00166   { m_CurvatureWeight = c; }
00167   ScalarValueType GetCurvatureWeight() const
00168   { return m_CurvatureWeight; }
00170 
00172   void SetLaplacianSmoothingWeight(const ScalarValueType c)
00173   { m_LaplacianSmoothingWeight = c; }
00174   ScalarValueType GetLaplacianSmoothingWeight() const
00175   { return m_LaplacianSmoothingWeight; }
00177 
00179   void SetEpsilonMagnitude(const ScalarValueType e)
00180   { m_EpsilonMagnitude = e; }
00181   ScalarValueType GetEpsilonMagnitude() const
00182   { return m_EpsilonMagnitude; }
00184 
00186   virtual PixelType ComputeUpdate( const NeighborhoodType & neighborhood,
00187                                    void *globalData,
00188                                    const FloatOffsetType & = FloatOffsetType(0.0) );
00189 
00196   virtual TimeStepType ComputeGlobalTimeStep(void *GlobalData) const;
00197 
00205   virtual void * GetGlobalDataPointer() const
00206   {
00207     GlobalDataStruct *ans = new GlobalDataStruct();
00208 
00209     ans->m_MaxAdvectionChange   = NumericTraits< ScalarValueType >::Zero;
00210     ans->m_MaxPropagationChange = NumericTraits< ScalarValueType >::Zero;
00211     ans->m_MaxCurvatureChange   = NumericTraits< ScalarValueType >::Zero;
00212     return ans;
00213   }
00214 
00218   virtual void Initialize(const RadiusType & r);
00219 
00224   virtual void ReleaseGlobalDataPointer(void *GlobalData) const
00225   { delete (GlobalDataStruct *)GlobalData; }
00226 
00228   virtual ScalarValueType ComputeCurvatureTerm(const NeighborhoodType &,
00229                                                const FloatOffsetType &,
00230                                                GlobalDataStruct *gd = 0
00231                                                );
00232 
00233   virtual ScalarValueType ComputeMeanCurvature(const NeighborhoodType &,
00234                                                const FloatOffsetType &,
00235                                                GlobalDataStruct *gd = 0
00236                                                );
00237 
00238   virtual ScalarValueType ComputeMinimalCurvature(const NeighborhoodType &,
00239                                                   const FloatOffsetType &,
00240                                                   GlobalDataStruct *gd = 0
00241                                                   );
00242 
00243   virtual ScalarValueType Compute3DMinimalCurvature(const NeighborhoodType &,
00244                                                     const FloatOffsetType &,
00245                                                     GlobalDataStruct *gd = 0
00246                                                     );
00247 
00249   void SetUseMinimalCurvature(bool b)
00250   {
00251     m_UseMinimalCurvature = b;
00252   }
00253 
00254   bool GetUseMinimalCurvature() const
00255   {
00256     return m_UseMinimalCurvature;
00257   }
00258 
00259   void UseMinimalCurvatureOn()
00260   {
00261     this->SetUseMinimalCurvature(true);
00262   }
00263 
00264   void UseMinimalCurvatureOff()
00265   {
00266     this->SetUseMinimalCurvature(false);
00267   }
00268 
00273   static void SetMaximumCurvatureTimeStep(double n)
00274   {
00275     m_DT = n;
00276   }
00277 
00278   static double GetMaximumCurvatureTimeStep()
00279   {
00280     return m_DT;
00281   }
00282 
00287   static void SetMaximumPropagationTimeStep(double n)
00288   {
00289     m_WaveDT = n;
00290   }
00291 
00292   static double GetMaximumPropagationTimeStep()
00293   {
00294     return m_WaveDT;
00295   }
00296 
00297 protected:
00298   LevelSetFunction()
00299   {
00300     m_EpsilonMagnitude = static_cast< ScalarValueType >( 1.0e-5 );
00301     m_AdvectionWeight = m_PropagationWeight =
00302                           m_CurvatureWeight = m_LaplacianSmoothingWeight =
00303                                                 NumericTraits< ScalarValueType >::Zero;
00304     m_UseMinimalCurvature = false;
00305   }
00306 
00307   virtual ~LevelSetFunction() {}
00308   void PrintSelf(std::ostream & s, Indent indent) const;
00309 
00311   static double m_WaveDT;
00312   static double m_DT;
00313 
00315   std::slice x_slice[itkGetStaticConstMacro(ImageDimension)];
00316 
00318   OffsetValueType m_Center;
00319 
00321   OffsetValueType m_xStride[itkGetStaticConstMacro(ImageDimension)];
00322 
00323   bool m_UseMinimalCurvature;
00324 
00327   static VectorType InitializeZeroVectorConstant();
00328 
00330   static VectorType m_ZeroVectorConstant;
00331 
00333   ScalarValueType m_EpsilonMagnitude;
00334 
00336   ScalarValueType m_AdvectionWeight;
00337 
00339   ScalarValueType m_PropagationWeight;
00340 
00342   ScalarValueType m_CurvatureWeight;
00343 
00345   ScalarValueType m_LaplacianSmoothingWeight;
00346 private:
00347   LevelSetFunction(const Self &); //purposely not implemented
00348   void operator=(const Self &);   //purposely not implemented
00349 };
00350 } // namespace itk
00352 
00353 // Define instantiation macro for this template.
00354 #define ITK_TEMPLATE_LevelSetFunction(_, EXPORT, TypeX, TypeY)     \
00355   namespace itk                                                    \
00356   {                                                                \
00357   _( 1 ( class EXPORT LevelSetFunction< ITK_TEMPLATE_1 TypeX > ) ) \
00358   namespace Templates                                              \
00359   {                                                                \
00360   typedef LevelSetFunction< ITK_TEMPLATE_1 TypeX >                 \
00361   LevelSetFunction##TypeY;                                       \
00362   }                                                                \
00363   }
00364 
00365 #if ITK_TEMPLATE_EXPLICIT
00366 #include "Templates/itkLevelSetFunction+-.h"
00367 #endif
00368 
00369 #if ITK_TEMPLATE_TXX
00370 #include "itkLevelSetFunction.hxx"
00371 #endif
00372 
00373 #endif
00374