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

itkFEMRegistrationFilter.h

Go to the documentation of this file.
00001 /*=========================================================================
00002 
00003   Program:   Insight Segmentation & Registration Toolkit
00004   Module:    $RCSfile: itkFEMRegistrationFilter.h,v $
00005   Language:  C++
00006   Date:      $Date: 2003/09/10 14:28:30 $
00007   Version:   $Revision: 1.20 $
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 
00018 #ifndef _itkFEMRegistrationFilter_h_
00019 #define _itkFEMRegistrationFilter_h_
00020 
00021 #include "itkFEMLinearSystemWrapperItpack.h"
00022 #include "itkFEMLinearSystemWrapperDenseVNL.h"
00023 #include "itkFEMGenerateMesh.h"
00024 #include "itkFEMSolverCrankNicolson.h"
00025 #include "itkFEMMaterialLinearElasticity.h"
00026 #include "itkFEMImageMetricLoad.h"
00027 #include "itkFEMFiniteDifferenceFunctionLoad.h"
00028 
00029 #include "itkImage.h"
00030 #include "itkVector.h"
00031 #include "itkImageRegionIteratorWithIndex.h"
00032 #include "itkVectorCastImageFilter.h"
00033 #include "itkVectorIndexSelectionCastImageFilter.h"
00034 #include "itkWarpImageFilter.h"
00035 #include "itkImageToImageMetric.h"
00036 #include "itkTranslationTransform.h"
00037 #include "itkVectorExpandImageFilter.h"
00038 
00039 #include "itkRecursiveMultiResolutionPyramidImageFilter.h"
00040 #include "itkFEMLoadLandmark.h"
00041 
00042 #include "vnl/vnl_vector.h"
00043 #include "vnl/vnl_math.h"
00044 #include "vnl/vnl_vector_fixed.h"
00045 
00046 
00047 #include <iostream>
00048 #include <string>
00049 
00050 
00051 
00052 namespace itk {
00053 namespace fem {
00054 
00101 template<class TMovingImage,class TFixedImage> 
00102 class  ITK_EXPORT  FEMRegistrationFilter : public ImageToImageFilter<TMovingImage, TFixedImage>
00103 {
00104 public:
00105   typedef FEMRegistrationFilter                              Self;
00106   typedef ImageToImageFilter<TMovingImage, TFixedImage> Superclass;
00107   typedef SmartPointer<Self> Pointer;
00108   typedef SmartPointer<const Self> ConstPointer;
00109 
00111   itkNewMacro(Self);
00112   
00114   itkTypeMacro(FEMRegistrationFilter, ImageToImageFilter );
00115   
00116   typedef TMovingImage                              MovingImageType;
00117   typedef TFixedImage                               FixedImageType;
00118   typedef typename FixedImageType::PixelType             PixelType;
00119   typedef typename FixedImageType::SizeType              ImageSizeType;
00120 
00122   itkStaticConstMacro(ImageDimension, unsigned int,
00123                       FixedImageType::ImageDimension);
00124 
00125   typedef Image< float, itkGetStaticConstMacro(ImageDimension) >            FloatImageType;
00126   typedef LinearSystemWrapperItpack                 LinearSystemSolverType;
00127   typedef SolverCrankNicolson                       SolverType;
00128   enum Sign { positive = 1, negative = -1 };
00129   typedef double                                    Float;
00130   typedef Load::ArrayType LoadArray;
00131     
00132 
00133   typedef std::vector<typename LoadLandmark::Pointer> LandmarkArrayType;
00134   typedef itk::Vector<float,itkGetStaticConstMacro(ImageDimension)>         VectorType;
00135   typedef itk::Image<VectorType,itkGetStaticConstMacro(ImageDimension)>     FieldType;
00136   typedef itk::WarpImageFilter<MovingImageType,FixedImageType, FieldType> WarperType;
00137   typedef MaterialLinearElasticity                  MaterialType;
00138   typedef itk::ImageRegionIteratorWithIndex<FixedImageType>         ImageIterator; 
00139   typedef itk::ImageRegionIteratorWithIndex<FloatImageType>         FloatImageIterator; 
00140   typedef itk::ImageRegionIteratorWithIndex<FieldType>         FieldIterator; 
00141   typedef itk::VectorIndexSelectionCastImageFilter<FieldType,FloatImageType> IndexSelectCasterType;
00142 
00143 
00145   typedef double CoordRepType;
00146   typedef VectorInterpolateImageFunction<FieldType,CoordRepType> 
00147     InterpolatorType;
00148   typedef typename InterpolatorType::Pointer InterpolatorPointer;
00149   typedef VectorLinearInterpolateImageFunction<FieldType,CoordRepType> 
00150     DefaultInterpolatorType;
00151 
00153   itkSetObjectMacro( Interpolator, InterpolatorType );
00154 
00156   itkGetObjectMacro( Interpolator, InterpolatorType );
00157 
00158 
00159   typedef itk::VectorExpandImageFilter<FieldType,FieldType> ExpanderType;
00160   typedef typename ExpanderType::ExpandFactorsType ExpandFactorsType;
00161 
00162   typedef itk::RecursiveMultiResolutionPyramidImageFilter<FixedImageType,FixedImageType>
00163     FixedPyramidType;
00164 
00166 //#define USEIMAGEMETRIC
00167 #ifdef  USEIMAGEMETRIC
00168   typedef ImageToImageMetric<ImageType,FixedImageType>   MetricBaseType;
00169   typedef  ImageMetricLoad<ImageType,ImageType>  ImageMetricLoadType;
00170 #else
00171   typedef  FiniteDifferenceFunctionLoad<MovingImageType,FixedImageType>  ImageMetricLoadType;
00172   typedef PDEDeformableRegistrationFunction<FixedImageType,MovingImageType,FieldType>   MetricBaseType;
00173 #endif
00174   typedef typename MetricBaseType::Pointer          MetricBaseTypePointer;
00175   /* Main functions */
00176  
00178   bool      ReadConfigFile(const char*);
00179   
00180 
00182   void      RunRegistration(); 
00183   
00187   void      WriteWarpedImage(const char* fn);
00188 
00189 
00191   void      IterativeSolve(SolverType& S);  
00192 
00194   void      MultiResSolve();
00195 
00197   void      WarpImage(const MovingImageType * R);      
00198 
00200   int       WriteDisplacementField(unsigned int index);
00201 
00203   void      SetMovingFile(const char* r) {m_MovingFileName=r;}
00204 
00205   std::string GetMovingFile() {return m_MovingFileName;}
00206   
00207   void      SetFixedFile(const char* t) {m_FixedFileName=t;}
00208   
00209   std::string GetFixedFile() {return m_FixedFileName;}
00210 
00211   
00215   void SetMovingImage(MovingImageType* R);
00216   
00218   void SetFixedImage(FixedImageType* T);
00219   
00220   MovingImageType* GetMovingImage(){return m_MovingImage;}
00221   MovingImageType* GetOriginalMovingImage(){return m_OriginalMovingImage;}
00222   
00223   FixedImageType* GetFixedImage(){return m_FixedImage;}
00224   
00225   
00228   FixedImageType* GetWarpedImage(){return m_WarpedImage;}
00229 
00231   void ComputeJacobian(float sign=1.0, FieldType* field=NULL , float smooth=0.0);
00232 
00234   FloatImageType* GetJacobianImage(){return m_FloatImage;}
00235   
00236 
00238   FieldType* GetDeformationField(){return m_Field;}
00240   void SetDeformationField(FieldType* F)
00241   {  
00242     m_FieldSize=F->GetLargestPossibleRegion().GetSize();
00243     m_Field=F;
00244   }
00245 
00248   void      SetLandmarkFile(const char* l) {m_LandmarkFileName=l; }
00249 
00251   void      UseLandmarks(bool b) {m_UseLandmarks=b;}
00252 
00253 
00261   void      EnforceDiffeomorphism(float thresh , SolverType& S ,  bool onlywriteimages);
00262 
00263 
00268   void      SetResultsFile(const char* r) {m_ResultsFileName=r;}
00269 
00270   void      SetResultsFileName (const char* f){m_ResultsFileName=f;}
00271 
00272   std::string GetResultsFileName () {return m_ResultsFileName;} 
00273 
00275   void      SetDisplacementsFile(const char* r) {m_DisplacementsFileName=r;}
00276   
00283   void      SetMeshPixelsPerElementAtEachResolution(unsigned int i,unsigned int which=0){ m_MeshPixelsPerElementAtEachResolution[which]=i;}
00284   
00288   void      SetNumberOfIntegrationPoints(unsigned int i,unsigned int which=0){ m_NumberOfIntegrationPoints[which]=i;}
00289   
00295   void      SetWidthOfMetricRegion(unsigned int i,unsigned int which=0) { m_MetricWidth[which]=i;}
00296   unsigned int      GetWidthOfMetricRegion(unsigned int which=0) { return m_MetricWidth[which];}
00297   
00302   void      SetMaximumIterations(unsigned int i,unsigned int which) { m_Maxiters[which]=i;}
00303   
00306   void      SetTimeStep(Float i) { m_dT=i;}
00307   
00309   void      SetAlpha(Float a) { m_Alpha=a;}
00310   
00313   void      SetEnergyReductionFactor(Float i) { m_EnergyReductionFactor=i;}
00314 
00316   void      SetElasticity(Float i,unsigned int which=0) { m_E[which]=i;} 
00317 
00319   Float     GetElasticity(unsigned int which=0) { return m_E[which];}
00320   
00322   void      SetRho(Float r,unsigned int which=0) { m_Rho[which]=r;} 
00323 
00325   void      SetGamma(Float r,unsigned int which=0) { m_Gamma[which]=r;} 
00326 
00328   void      SetDescentDirectionMinimize() { m_DescentDirection=positive;}
00329   
00331   void      SetDescentDirectionMaximize() { m_DescentDirection=negative;} 
00332 
00334   void      DoLineSearch(unsigned int b) { m_DoLineSearchOnImageEnergy=b; } 
00335 
00336 
00338   void      DoMultiRes(bool b) { m_DoMultiRes=b; } 
00339 
00341   void      EmployRegridding(unsigned int b) { m_EmployRegridding=b; } 
00342 
00344   void      SetLineSearchMaximumIterations(unsigned int f) { m_LineSearchMaximumIterations=f; } 
00345   
00347   void      SetWriteDisplacements(bool b) {m_WriteDisplacementField=b;}
00348   
00350   bool      GetWriteDisplacements() {return m_WriteDisplacementField;}
00351     
00354   void      SetConfigFileName (const char* f){m_ConfigFileName=f;}
00355   
00356   std::string GetConfigFileName () {return m_ConfigFileName; }
00357   
00358   ImageSizeType GetImageSize(){ return m_FullImageSize; }
00359 
00361   MetricBaseTypePointer    GetMetric() { return m_Metric; }
00362   void      SetMetric(MetricBaseTypePointer MP) { m_Metric=MP; }
00363   
00366   void      ChooseMetric( float whichmetric); 
00367   
00369   void      SetElement(Element::Pointer e) {m_Element=e;}
00370 
00372   void      SetMaterial(MaterialType::Pointer m) {m_Material=m;}
00373 
00374   void      PrintVectorField(unsigned int modnum=1000);
00375 
00376   void      SetNumLevels(unsigned int i) { m_NumLevels=i; }
00377   void      SetMaxLevel(unsigned int i) { m_MaxLevel=i; }
00378 
00379   void      SetTemp(Float i) { m_Temp=i; }
00380 
00382   FEMRegistrationFilter( ); 
00383   ~FEMRegistrationFilter(); 
00384 
00385 // HELPER FUNCTIONS
00386 protected :
00387 
00388   
00393   class FEMOF : public FEMObjectFactory<FEMLightObject>{
00394   protected:
00395     FEMOF();
00396     ~FEMOF();
00397     };
00398 
00399  
00401   void      CreateMesh(double ElementsPerSide, Solver& S, ImageSizeType sz);
00402 
00404   void      ApplyLoads(SolverType& S,ImageSizeType Isz,double* spacing=NULL); 
00405 
00407   void      ApplyImageLoads(SolverType& S, MovingImageType* i1, FixedImageType* i2); 
00408 
00409   
00412   void      CreateLinearSystemSolver();
00413 
00414   
00416   Float     EvaluateEnergy();
00417  
00422   void      InterpolateVectorField(SolverType& S); 
00423 
00426   FloatImageType*      GetMetricImage(FieldType* F); 
00427 
00428  
00430   typedef  typename FieldType::Pointer FieldPointer;
00431   FieldPointer ExpandVectorField(ExpandFactorsType* expandFactors, FieldType* f);
00432   
00433 
00435   void      SampleVectorFieldAtNodes(SolverType& S);
00436   
00437   Float EvaluateResidual(SolverType& mySolver,Float t);
00438 
00439   /* Finds a triplet that brackets the energy minimum.  From Numerical Recipes.*/
00440   void FindBracketingTriplet(SolverType& mySolver,Float* a, Float* b, Float* c);
00441   
00445   Float GoldenSection(SolverType& mySolver,Float tol=0.01,unsigned int MaxIters=25);
00446 
00448 //  itkSetMacro( Load, ImageMetricLoadType* );
00449   itkGetMacro( Load, ImageMetricLoadType* );
00450 
00451 
00452   void PrintSelf(std::ostream& os, Indent indent) const;
00453 
00454 private :
00455 
00456   void InitializeField();
00457 
00458   FEMRegistrationFilter(const Self&); //purposely not implemented
00459   void operator=(const Self&); //purposely not implemented
00460     
00461   std::string      m_ConfigFileName;
00462   std::string      m_ResultsFileName;
00463   std::string      m_MovingFileName;  
00464   std::string      m_FixedFileName;
00465   std::string      m_LandmarkFileName;
00466   std::string      m_DisplacementsFileName;
00467   std::string      m_MeshFileName;
00468 
00469   unsigned int     m_DoLineSearchOnImageEnergy; 
00470   unsigned int     m_LineSearchMaximumIterations;
00471 
00472   vnl_vector<unsigned int>     m_NumberOfIntegrationPoints;// resolution of integration
00473   vnl_vector<unsigned int>     m_MetricWidth;
00474   vnl_vector<unsigned int>     m_Maxiters; // max iterations
00475   unsigned int     m_TotalIterations;
00476   unsigned int     m_NumLevels; // Number of Resolution Levels
00477   unsigned int     m_MaxLevel;  // Maximum Level (NumLevels is original resolution).
00478   unsigned int     m_MeshLevels;// Number of Mesh Resolutions ( should be >= 1)
00479   unsigned int     m_MeshStep;  // Ratio Between Mesh Resolutions ( currently set to 2, should be >= 1)
00480   unsigned int     m_FileCount; // keeps track of number of files written
00481   unsigned int     m_CurrentLevel;
00482   typename FixedImageType::SizeType     m_CurrentLevelImageSize; 
00483   unsigned int     m_WhichMetric;
00484  
00487   vnl_vector<unsigned int> m_MeshPixelsPerElementAtEachResolution;
00488 
00489   Float     m_dT; // time step
00490   vnl_vector<Float>     m_E;  // elasticity 
00491   vnl_vector<Float>     m_Rho;   // mass matrix weight
00492   vnl_vector<Float>     m_Gamma;   // image similarity weight
00493   Float     m_Energy; // current value of energy
00494   Float     m_MinE;  // minimum recorded energy
00495   Float     m_MinJacobian;  // minimum recorded energy
00496   Float     m_Alpha; // difference parameter 
00498   Float     m_EnergyReductionFactor; 
00499   Float     m_Temp;
00500 
00501   bool  m_WriteDisplacementField;
00502   bool  m_DoMultiRes;
00503   bool  m_UseLandmarks;
00504   bool  m_ReadMeshFile;
00505   bool  m_UseMassMatrix;
00506   unsigned int m_EmployRegridding;
00507   Sign  m_DescentDirection;
00508 
00509   ImageSizeType     m_FullImageSize; // image size
00510   ImageSizeType     m_ImageOrigin; // image size
00512   ImageSizeType     m_ImageScaling; 
00513   ImageSizeType     m_CurrentImageScaling; 
00514   typename FieldType::RegionType   m_FieldRegion;
00515   typename FieldType::SizeType     m_FieldSize;
00516   typename FieldType::Pointer      m_Field;
00517   // only use TotalField if re-gridding is employed.
00518   typename FieldType::Pointer      m_TotalField;
00519   ImageMetricLoadType* m_Load; // Defines the load to use
00520    
00521   // define the warper
00522   typename WarperType::Pointer m_Warper; 
00523 
00524  // declare a new image to hold the warped  reference
00525   typename FixedImageType::Pointer      m_WarpedImage;
00526   typename FloatImageType::Pointer      m_FloatImage;
00527   typename FixedImageType::RegionType   m_Wregion; 
00528   typename FixedImageType::IndexType    m_Windex;
00529  
00530  // declare images for target and reference
00531   typename MovingImageType::Pointer      m_MovingImage;
00532   typename MovingImageType::Pointer      m_OriginalMovingImage;
00533   typename FixedImageType::Pointer      m_FixedImage;
00534 
00535   // element and metric pointers
00536   typename Element::Pointer        m_Element;
00537   typename MaterialType::Pointer   m_Material;
00538   MetricBaseTypePointer            m_Metric;
00539 
00540 
00541   // multi-resolution stuff
00542 //  typename FixedPyramidType::Pointer   m_FixedPyramid;
00543 //  typename FixedPyramidType::Pointer   m_MovingPyramid;
00544   
00545   LandmarkArrayType    m_LandmarkArray;
00546 
00547   InterpolatorPointer    m_Interpolator;
00548 
00549  
00550 
00551 };
00552 
00553 }} // end namespace itk::fem
00554 
00555 #ifndef ITK_MANUAL_INSTANTIATION
00556 #include "itkFEMRegistrationFilter.txx"
00557 #endif
00558 
00559 #endif
00560 

Generated at Tue Sep 16 11:31:58 2003 for ITK by doxygen 1.2.15 written by Dimitri van Heesch, © 1997-2000