00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017 #ifndef __itkSymmetricForcesDemonsRegistrationFunction_h
00018 #define __itkSymmetricForcesDemonsRegistrationFunction_h
00019
00020 #include "itkPDEDeformableRegistrationFunction.h"
00021 #include "itkPoint.h"
00022 #include "itkCovariantVector.h"
00023 #include "itkInterpolateImageFunction.h"
00024 #include "itkLinearInterpolateImageFunction.h"
00025 #include "itkCentralDifferenceImageFunction.h"
00026
00027 namespace itk {
00028
00060 template<class TFixedImage, class TMovingImage, class TDeformationField>
00061 class ITK_EXPORT SymmetricForcesDemonsRegistrationFunction :
00062 public PDEDeformableRegistrationFunction< TFixedImage,
00063 TMovingImage, TDeformationField>
00064 {
00065 public:
00066
00068 typedef SymmetricForcesDemonsRegistrationFunction Self;
00069 typedef PDEDeformableRegistrationFunction< TFixedImage,
00070 TMovingImage, TDeformationField >
00071 Superclass;
00072 typedef SmartPointer<Self> Pointer;
00073 typedef SmartPointer<const Self> ConstPointer;
00074
00076 itkNewMacro(Self);
00077
00079 itkTypeMacro( SymmetricForcesDemonsRegistrationFunction,
00080 PDEDeformableRegistrationFunction );
00081
00083 typedef typename Superclass::MovingImageType MovingImageType;
00084 typedef typename Superclass::MovingImagePointer MovingImagePointer;
00085
00087 typedef typename Superclass::FixedImageType FixedImageType;
00088 typedef typename Superclass::FixedImagePointer FixedImagePointer;
00089 typedef typename FixedImageType::IndexType IndexType;
00090 typedef typename FixedImageType::SizeType SizeType;
00091 typedef typename FixedImageType::SpacingType SpacingType;
00092
00094 typedef typename Superclass::DeformationFieldType DeformationFieldType;
00095 typedef typename Superclass::DeformationFieldTypePointer
00096 DeformationFieldTypePointer;
00097
00099 itkStaticConstMacro(ImageDimension, unsigned int,Superclass::ImageDimension);
00100
00102 typedef typename Superclass::PixelType PixelType;
00103 typedef typename Superclass::RadiusType RadiusType;
00104 typedef typename Superclass::NeighborhoodType NeighborhoodType;
00105 typedef typename Superclass::FloatOffsetType FloatOffsetType;
00106 typedef typename Superclass::TimeStepType TimeStepType;
00107
00109 typedef double CoordRepType;
00110 typedef InterpolateImageFunction<MovingImageType,CoordRepType>
00111 InterpolatorType;
00112 typedef typename InterpolatorType::Pointer InterpolatorPointer;
00113 typedef typename InterpolatorType::PointType PointType;
00114 typedef LinearInterpolateImageFunction<MovingImageType,CoordRepType>
00115 DefaultInterpolatorType;
00116
00118 typedef CovariantVector<double,itkGetStaticConstMacro(ImageDimension)> CovariantVectorType;
00119
00121 typedef CentralDifferenceImageFunction<FixedImageType> GradientCalculatorType;
00122 typedef typename GradientCalculatorType::Pointer GradientCalculatorPointer;
00123
00125 void SetMovingImageInterpolator( InterpolatorType * ptr )
00126 { m_MovingImageInterpolator = ptr; }
00127
00129 InterpolatorType * GetMovingImageInterpolator(void)
00130 { return m_MovingImageInterpolator; }
00131
00133 virtual TimeStepType ComputeGlobalTimeStep(void * itkNotUsed(GlobalData)) const
00134 { return m_TimeStep; }
00135
00138 virtual void *GetGlobalDataPointer() const
00139 {
00140 GlobalDataStruct *global = new GlobalDataStruct();
00141 global->m_SumOfSquaredDifference = 0.0;
00142 global->m_NumberOfPixelsProcessed = 0L;
00143 global->m_SumOfSquaredChange = 0;
00144 return global;
00145 }
00147
00149 virtual void ReleaseGlobalDataPointer( void *GlobalData ) const;
00150
00152 virtual void InitializeIteration();
00153
00156 virtual PixelType ComputeUpdate(const NeighborhoodType &neighborhood,
00157 void *globalData,
00158 const FloatOffsetType &offset = FloatOffsetType(0.0));
00159
00163 virtual double GetMetric() const
00164 { return m_Metric; }
00165
00167 virtual const double &GetRMSChange() const
00168 { return m_RMSChange; }
00169
00174 virtual void SetIntensityDifferenceThreshold(double);
00175 virtual double GetIntensityDifferenceThreshold() const;
00177
00178 protected:
00179 SymmetricForcesDemonsRegistrationFunction();
00180 ~SymmetricForcesDemonsRegistrationFunction() {}
00181 void PrintSelf(std::ostream& os, Indent indent) const;
00182
00184 typedef ConstNeighborhoodIterator<FixedImageType> FixedImageNeighborhoodIteratorType;
00185
00188 struct GlobalDataStruct
00189 {
00190 double m_SumOfSquaredDifference;
00191 unsigned long m_NumberOfPixelsProcessed;
00192 double m_SumOfSquaredChange;
00193 };
00194
00195 private:
00196 SymmetricForcesDemonsRegistrationFunction(const Self&);
00197 void operator=(const Self&);
00198
00200 SpacingType m_FixedImageSpacing;
00201 PointType m_FixedImageOrigin;
00202 double m_Normalizer;
00203
00205 GradientCalculatorPointer m_FixedImageGradientCalculator;
00206
00208 InterpolatorPointer m_MovingImageInterpolator;
00209
00211 TimeStepType m_TimeStep;
00212
00214 double m_DenominatorThreshold;
00215
00217 double m_IntensityDifferenceThreshold;
00218
00222 mutable double m_Metric;
00223 mutable double m_SumOfSquaredDifference;
00224 mutable unsigned long m_NumberOfPixelsProcessed;
00225 mutable double m_RMSChange;
00226 mutable double m_SumOfSquaredChange;
00227
00229 mutable SimpleFastMutexLock m_MetricCalculationLock;
00230
00231 };
00232
00233
00234 }
00235
00236 #ifndef ITK_MANUAL_INSTANTIATION
00237 #include "itkSymmetricForcesDemonsRegistrationFunction.txx"
00238 #endif
00239
00240 #endif
00241