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/11/28 03:09:57 $ 00007 Version: $Revision: 1.21 $ 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 int WriteDisplacementFieldMultiComponent(); 00204 00206 void SetMovingFile(const char* r) {m_MovingFileName=r;} 00207 00208 std::string GetMovingFile() {return m_MovingFileName;} 00209 00210 void SetFixedFile(const char* t) {m_FixedFileName=t;} 00211 00212 std::string GetFixedFile() {return m_FixedFileName;} 00213 00214 00218 void SetMovingImage(MovingImageType* R); 00219 00221 void SetFixedImage(FixedImageType* T); 00222 00223 MovingImageType* GetMovingImage(){return m_MovingImage;} 00224 MovingImageType* GetOriginalMovingImage(){return m_OriginalMovingImage;} 00225 00226 FixedImageType* GetFixedImage(){return m_FixedImage;} 00227 00228 00231 FixedImageType* GetWarpedImage(){return m_WarpedImage;} 00232 00234 void ComputeJacobian(float sign=1.0, FieldType* field=NULL , float smooth=0.0); 00235 00237 FloatImageType* GetJacobianImage(){return m_FloatImage;} 00238 00239 00241 FieldType* GetDeformationField(){return m_Field;} 00243 void SetDeformationField(FieldType* F) 00244 { 00245 m_FieldSize=F->GetLargestPossibleRegion().GetSize(); 00246 m_Field=F; 00247 } 00248 00251 void SetLandmarkFile(const char* l) {m_LandmarkFileName=l; } 00252 00254 void UseLandmarks(bool b) {m_UseLandmarks=b;} 00255 00256 00264 void EnforceDiffeomorphism(float thresh , SolverType& S , bool onlywriteimages); 00265 00266 00271 void SetResultsFile(const char* r) {m_ResultsFileName=r;} 00272 00273 void SetResultsFileName (const char* f){m_ResultsFileName=f;} 00274 00275 std::string GetResultsFileName () {return m_ResultsFileName;} 00276 00278 void SetDisplacementsFile(const char* r) {m_DisplacementsFileName=r;} 00279 00286 void SetMeshPixelsPerElementAtEachResolution(unsigned int i,unsigned int which=0){ m_MeshPixelsPerElementAtEachResolution[which]=i;} 00287 00291 void SetNumberOfIntegrationPoints(unsigned int i,unsigned int which=0){ m_NumberOfIntegrationPoints[which]=i;} 00292 00298 void SetWidthOfMetricRegion(unsigned int i,unsigned int which=0) { m_MetricWidth[which]=i;} 00299 unsigned int GetWidthOfMetricRegion(unsigned int which=0) { return m_MetricWidth[which];} 00300 00305 void SetMaximumIterations(unsigned int i,unsigned int which) { m_Maxiters[which]=i;} 00306 00309 void SetTimeStep(Float i) { m_dT=i;} 00310 00312 void SetAlpha(Float a) { m_Alpha=a;} 00313 00316 void SetEnergyReductionFactor(Float i) { m_EnergyReductionFactor=i;} 00317 00319 void SetElasticity(Float i,unsigned int which=0) { m_E[which]=i;} 00320 00322 Float GetElasticity(unsigned int which=0) { return m_E[which];} 00323 00325 void SetRho(Float r,unsigned int which=0) { m_Rho[which]=r;} 00326 00328 void SetGamma(Float r,unsigned int which=0) { m_Gamma[which]=r;} 00329 00331 void SetDescentDirectionMinimize() { m_DescentDirection=positive;} 00332 00334 void SetDescentDirectionMaximize() { m_DescentDirection=negative;} 00335 00337 void DoLineSearch(unsigned int b) { m_DoLineSearchOnImageEnergy=b; } 00338 00339 00341 void DoMultiRes(bool b) { m_DoMultiRes=b; } 00342 00344 void EmployRegridding(unsigned int b) { m_EmployRegridding=b; } 00345 00347 void SetLineSearchMaximumIterations(unsigned int f) { m_LineSearchMaximumIterations=f; } 00348 00350 void SetWriteDisplacements(bool b) {m_WriteDisplacementField=b;} 00351 00353 bool GetWriteDisplacements() {return m_WriteDisplacementField;} 00354 00357 void SetConfigFileName (const char* f){m_ConfigFileName=f;} 00358 00359 std::string GetConfigFileName () {return m_ConfigFileName; } 00360 00361 ImageSizeType GetImageSize(){ return m_FullImageSize; } 00362 00364 MetricBaseTypePointer GetMetric() { return m_Metric; } 00365 void SetMetric(MetricBaseTypePointer MP) { m_Metric=MP; } 00366 00369 void ChooseMetric( float whichmetric); 00370 00372 void SetElement(Element::Pointer e) {m_Element=e;} 00373 00375 void SetMaterial(MaterialType::Pointer m) {m_Material=m;} 00376 00377 void PrintVectorField(unsigned int modnum=1000); 00378 00379 void SetNumLevels(unsigned int i) { m_NumLevels=i; } 00380 void SetMaxLevel(unsigned int i) { m_MaxLevel=i; } 00381 00382 void SetTemp(Float i) { m_Temp=i; } 00383 00385 FEMRegistrationFilter( ); 00386 ~FEMRegistrationFilter(); 00387 00388 // HELPER FUNCTIONS 00389 protected : 00390 00391 00396 class FEMOF : public FEMObjectFactory<FEMLightObject>{ 00397 protected: 00398 FEMOF(); 00399 ~FEMOF(); 00400 }; 00401 00402 00404 void CreateMesh(double ElementsPerSide, Solver& S, ImageSizeType sz); 00405 00407 void ApplyLoads(SolverType& S,ImageSizeType Isz,double* spacing=NULL); 00408 00410 void ApplyImageLoads(SolverType& S, MovingImageType* i1, FixedImageType* i2); 00411 00412 00415 void CreateLinearSystemSolver(); 00416 00417 00419 Float EvaluateEnergy(); 00420 00425 void InterpolateVectorField(SolverType& S); 00426 00429 FloatImageType* GetMetricImage(FieldType* F); 00430 00431 00433 typedef typename FieldType::Pointer FieldPointer; 00434 FieldPointer ExpandVectorField(ExpandFactorsType* expandFactors, FieldType* f); 00435 00436 00438 void SampleVectorFieldAtNodes(SolverType& S); 00439 00440 Float EvaluateResidual(SolverType& mySolver,Float t); 00441 00442 /* Finds a triplet that brackets the energy minimum. From Numerical Recipes.*/ 00443 void FindBracketingTriplet(SolverType& mySolver,Float* a, Float* b, Float* c); 00444 00448 Float GoldenSection(SolverType& mySolver,Float tol=0.01,unsigned int MaxIters=25); 00449 00451 // itkSetMacro( Load, ImageMetricLoadType* ); 00452 itkGetMacro( Load, ImageMetricLoadType* ); 00453 00454 00455 void PrintSelf(std::ostream& os, Indent indent) const; 00456 00457 private : 00458 00459 void InitializeField(); 00460 00461 FEMRegistrationFilter(const Self&); //purposely not implemented 00462 void operator=(const Self&); //purposely not implemented 00463 00464 std::string m_ConfigFileName; 00465 std::string m_ResultsFileName; 00466 std::string m_MovingFileName; 00467 std::string m_FixedFileName; 00468 std::string m_LandmarkFileName; 00469 std::string m_DisplacementsFileName; 00470 std::string m_MeshFileName; 00471 00472 unsigned int m_DoLineSearchOnImageEnergy; 00473 unsigned int m_LineSearchMaximumIterations; 00474 00475 vnl_vector<unsigned int> m_NumberOfIntegrationPoints;// resolution of integration 00476 vnl_vector<unsigned int> m_MetricWidth; 00477 vnl_vector<unsigned int> m_Maxiters; // max iterations 00478 unsigned int m_TotalIterations; 00479 unsigned int m_NumLevels; // Number of Resolution Levels 00480 unsigned int m_MaxLevel; // Maximum Level (NumLevels is original resolution). 00481 unsigned int m_MeshLevels;// Number of Mesh Resolutions ( should be >= 1) 00482 unsigned int m_MeshStep; // Ratio Between Mesh Resolutions ( currently set to 2, should be >= 1) 00483 unsigned int m_FileCount; // keeps track of number of files written 00484 unsigned int m_CurrentLevel; 00485 typename FixedImageType::SizeType m_CurrentLevelImageSize; 00486 unsigned int m_WhichMetric; 00487 00490 vnl_vector<unsigned int> m_MeshPixelsPerElementAtEachResolution; 00491 00492 Float m_dT; // time step 00493 vnl_vector<Float> m_E; // elasticity 00494 vnl_vector<Float> m_Rho; // mass matrix weight 00495 vnl_vector<Float> m_Gamma; // image similarity weight 00496 Float m_Energy; // current value of energy 00497 Float m_MinE; // minimum recorded energy 00498 Float m_MinJacobian; // minimum recorded energy 00499 Float m_Alpha; // difference parameter 00501 Float m_EnergyReductionFactor; 00502 Float m_Temp; 00503 00504 bool m_WriteDisplacementField; 00505 bool m_DoMultiRes; 00506 bool m_UseLandmarks; 00507 bool m_ReadMeshFile; 00508 bool m_UseMassMatrix; 00509 unsigned int m_EmployRegridding; 00510 Sign m_DescentDirection; 00511 00512 ImageSizeType m_FullImageSize; // image size 00513 ImageSizeType m_ImageOrigin; // image size 00515 ImageSizeType m_ImageScaling; 00516 ImageSizeType m_CurrentImageScaling; 00517 typename FieldType::RegionType m_FieldRegion; 00518 typename FieldType::SizeType m_FieldSize; 00519 typename FieldType::Pointer m_Field; 00520 // only use TotalField if re-gridding is employed. 00521 typename FieldType::Pointer m_TotalField; 00522 ImageMetricLoadType* m_Load; // Defines the load to use 00523 00524 // define the warper 00525 typename WarperType::Pointer m_Warper; 00526 00527 // declare a new image to hold the warped reference 00528 typename FixedImageType::Pointer m_WarpedImage; 00529 typename FloatImageType::Pointer m_FloatImage; 00530 typename FixedImageType::RegionType m_Wregion; 00531 typename FixedImageType::IndexType m_Windex; 00532 00533 // declare images for target and reference 00534 typename MovingImageType::Pointer m_MovingImage; 00535 typename MovingImageType::Pointer m_OriginalMovingImage; 00536 typename FixedImageType::Pointer m_FixedImage; 00537 00538 // element and metric pointers 00539 typename Element::Pointer m_Element; 00540 typename MaterialType::Pointer m_Material; 00541 MetricBaseTypePointer m_Metric; 00542 00543 00544 // multi-resolution stuff 00545 // typename FixedPyramidType::Pointer m_FixedPyramid; 00546 // typename FixedPyramidType::Pointer m_MovingPyramid; 00547 00548 LandmarkArrayType m_LandmarkArray; 00549 00550 InterpolatorPointer m_Interpolator; 00551 00552 00553 00554 }; 00555 00556 }} // end namespace itk::fem 00557 00558 #ifndef ITK_MANUAL_INSTANTIATION 00559 #include "itkFEMRegistrationFilter.txx" 00560 #endif 00561 00562 #endif 00563

Generated at Sun Apr 1 02:28:25 2007 for ITK by doxygen 1.3.8 written by Dimitri van Heesch, © 1997-2000