]> Creatis software - clitk.git/blob - registration/clitkAffineRegistrationGenericFilter.txx
vv now builds even with all options turned off
[clitk.git] / registration / clitkAffineRegistrationGenericFilter.txx
1 /*=========================================================================
2   Program:   vv                     http://www.creatis.insa-lyon.fr/rio/vv
3
4   Authors belong to:
5   - University of LYON              http://www.universite-lyon.fr/
6   - Léon Bérard cancer center       http://oncora1.lyon.fnclcc.fr
7   - CREATIS CNRS laboratory         http://www.creatis.insa-lyon.fr
8
9   This software is distributed WITHOUT ANY WARRANTY; without even
10   the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
11   PURPOSE.  See the copyright notices for more information.
12
13   It is distributed under dual licence
14
15   - BSD        See included LICENSE.txt file
16   - CeCILL-B   http://www.cecill.info/licences/Licence_CeCILL-B_V1-en.html
17 ======================================================================-====*/
18
19 #ifndef CLITKAFFINEREGISTRATIONGENERICFILTER_TXX
20 #define CLITKAFFINEREGISTRATIONGENERICFILTER_TXX
21
22 namespace clitk
23 {
24
25 //==============================================================================
26 //Creating an observer class that allows us to monitor the registration
27 //================================================================================
28 class CommandIterationUpdate : public itk::Command
29 {
30 public:
31   typedef  CommandIterationUpdate   Self;
32   typedef  itk::Command             Superclass;
33   typedef  itk::SmartPointer<Self>  Pointer;
34   itkNewMacro( Self );
35 protected:
36   CommandIterationUpdate() {};
37 public:
38   typedef   clitk::GenericOptimizer<args_info_clitkAffineRegistration>     OptimizerType;
39   typedef   const OptimizerType   *           OptimizerPointer;
40
41   // Set the generic optimizer
42   void SetOptimizer(OptimizerPointer o) {
43     m_Optimizer=o;
44   }
45
46   // Execute
47   void Execute(itk::Object *caller, const itk::EventObject & event) {
48     Execute( (const itk::Object *)caller, event);
49   }
50
51   void Execute(const itk::Object * object, const itk::EventObject & event) {
52     if ( !(itk::IterationEvent().CheckEvent( &event )) ) {
53       return;
54     }
55
56     m_Optimizer->OutputIterationInfo();
57   }
58
59   OptimizerPointer m_Optimizer;
60 };
61 //==================================================================================================================================//
62 //Constructor
63 //===================================================================================================================================//
64
65 template<class args_info_clitkAffineRegistration>
66 AffineRegistrationGenericFilter<args_info_clitkAffineRegistration>::AffineRegistrationGenericFilter():
67     ImageToImageGenericFilter<Self>("Register")
68
69 {
70   InitializeImageType<2>();
71   InitializeImageType<3>();
72   m_Verbose=false;
73 }
74 //==========================================================================================================//
75 //============================================================================================================//
76 //--------------------------------------------------------------------
77 template<class args_info_type>
78 template<unsigned int Dim>
79 void AffineRegistrationGenericFilter<args_info_type>::InitializeImageType()
80 {
81   ADD_DEFAULT_IMAGE_TYPES(Dim);
82 }
83 //--------------------------------------------------------------------
84
85
86 //==============================================================================
87 //Creating an observer class that allows us to change parameters at subsequent levels
88 //==============================================================================
89 template <typename TRegistration, class args_info_clitkAffineRegistration>
90 class RegistrationInterfaceCommand : public itk::Command
91 {
92 public:
93   typedef  RegistrationInterfaceCommand   Self;
94   typedef  itk::Command             Superclass;
95   typedef itk::SmartPointer<Self>  Pointer;
96   itkNewMacro( Self );
97 protected:
98   RegistrationInterfaceCommand() {};
99 public:
100
101   // Registration
102   typedef   TRegistration                              RegistrationType;
103   typedef   RegistrationType *                         RegistrationPointer;
104
105   // Metric
106   typedef typename RegistrationType::FixedImageType    InternalImageType;
107   typedef clitk::GenericMetric<args_info_clitkAffineRegistration, InternalImageType, InternalImageType> GenericMetricType;
108
109   // Two arguments are passed to the Execute() method: the first
110   // is the pointer to the object which invoked the event and the
111   // second is the event that was invoked.
112   void Execute(itk::Object * object, const itk::EventObject & event) {
113     if ( !(itk::IterationEvent().CheckEvent( &event )) ) {
114       return;
115     }
116
117     // Get the levels
118     RegistrationPointer registration = dynamic_cast<RegistrationPointer>( object );
119     unsigned int numberOfLevels=registration->GetNumberOfLevels();
120     unsigned int currentLevel=registration->GetCurrentLevel()+1;
121
122     // Output the levels
123     std::cout<<std::endl;
124     std::cout<<"========================================"<<std::endl;
125     std::cout<<"Starting resolution level "<<currentLevel<<" of "<<numberOfLevels<<"..."<<std::endl;
126     std::cout<<"========================================"<<std::endl;
127     std::cout<<std::endl;
128
129     // Higher level?
130     if (currentLevel>1) {
131       // Reinitialize the metric (!= number of samples)
132       typename GenericMetricType::Pointer genericMetric= GenericMetricType::New();
133       genericMetric->SetArgsInfo(m_ArgsInfo);
134       genericMetric->SetFixedImage(registration->GetFixedImagePyramid()->GetOutput(registration->GetCurrentLevel()));
135       if (m_ArgsInfo.referenceMask_given)  genericMetric->SetFixedImageMask(registration->GetMetric()->GetFixedImageMask());
136       typedef itk::ImageToImageMetric< InternalImageType, InternalImageType >  MetricType;
137       typename  MetricType::Pointer metric=genericMetric->GetMetricPointer();
138       registration->SetMetric(metric);
139     }
140   }
141
142   void Execute(const itk::Object * , const itk::EventObject & ) {
143     return;
144   }
145
146  void SetArgsInfo(args_info_clitkAffineRegistration a) {
147     m_ArgsInfo=a;
148   }
149   args_info_clitkAffineRegistration m_ArgsInfo;
150 };
151
152 //==============================================================================================//
153 // ArgsInfo
154 //==============================================================================================//
155 template<class args_info_clitkAffineRegistration>
156 void AffineRegistrationGenericFilter<args_info_clitkAffineRegistration>::SetArgsInfo(const args_info_clitkAffineRegistration & a)
157 {
158   m_ArgsInfo=a;
159   if (m_ArgsInfo.reference_given) AddInputFilename(m_ArgsInfo.reference_arg);
160   if (m_ArgsInfo.target_given) {
161     AddInputFilename(m_ArgsInfo.target_arg);
162   }
163
164   if (m_ArgsInfo.output_given) SetOutputFilename(m_ArgsInfo.output_arg);
165 }
166 //==============================================================================
167 // Update with the number of dimensions and pixeltype
168 //==============================================================================
169 template<class args_info_clitkAffineRegistration>
170 template<class InputImageType>
171 void AffineRegistrationGenericFilter<args_info_clitkAffineRegistration>::UpdateWithInputImageType()
172 {
173   //=============================================================================
174   //Input
175   //=============================================================================
176
177   typedef typename  InputImageType::PixelType PixelType;
178 //typedef typename InputImageType::ImageDimension Dimension;
179
180
181   bool threadsGiven=m_ArgsInfo.threads_given;
182   int threads=m_ArgsInfo.threads_arg;
183
184   //Coordinate Representation
185   typedef double TCoordRep;
186
187
188   typename InputImageType::Pointer fixedImage = this->template GetInput<InputImageType>(0);
189
190   typename InputImageType::Pointer inputFixedImage = this->template GetInput<InputImageType>(0);
191
192   // typedef input2
193   typename InputImageType::Pointer movingImage = this->template GetInput<InputImageType>(1);
194
195   typename InputImageType::Pointer inputMovingImage = this->template GetInput<InputImageType>(1);
196
197
198
199   //The pixeltype of the fixed image will be used for output
200   typedef itk::Image< PixelType, InputImageType::ImageDimension > FixedImageType;
201
202   //Whatever the pixel type, internally we work with an image represented in float
203   typedef typename  InputImageType::PixelType  InternalPixelType;
204   typedef itk::Image< InternalPixelType, InputImageType::ImageDimension > InternalImageType;
205
206   //Read in the reference/fixed image
207 //  typedef itk::ImageFileReader< InternalImageType > ReaderType;
208 //  typename ReaderType::Pointer  fixedImageReader  = ReaderType::New();
209 //  fixedImageReader->SetFileName( m_ArgsInfo.reference_arg);
210
211
212   //Read in the object/moving image
213 //  typename ReaderType::Pointer movingImageReader = ReaderType::New();
214 //  movingImageReader->SetFileName( m_ArgsInfo.target_arg );
215   if (m_Verbose) std::cout<<"Reading images..."<<std::endl;
216 //  fixedImageReader->Update();
217 //  movingImageReader->Update();
218
219   if (m_Verbose) std::cout  << "Reading images... " << std::endl;
220
221   //we connect pointers to these internal images
222  // typedef typename  fixedImageReader fixedImage;
223  // typedef typename  movingImageReader movingImage;
224
225   //We keep the images used for input for possible output
226 // typedef typename   fixedImageReader inputFixedImage;
227 // typedef typename  movingImageReader inputMovingImage;
228
229
230   //============================================================================
231   // Preprocessing
232   //============================================================================
233
234   //If given, the intensities of both images are first normalized to a zero mean and SD of 1
235   // (usefull for MI, not necessary for Mattes' MI but performed anyway for the ouput)
236   if ( m_ArgsInfo.normalize_flag ) {
237     typedef itk::NormalizeImageFilter< InternalImageType,InternalImageType >  NormalizeFilterType;
238
239     typename  NormalizeFilterType::Pointer  fixedNormalizeFilter = NormalizeFilterType::New();
240     typename  NormalizeFilterType::Pointer  movingNormalizeFilter = NormalizeFilterType::New();
241
242     fixedNormalizeFilter->SetInput( fixedImage );
243     movingNormalizeFilter->SetInput( movingImage );
244
245     fixedNormalizeFilter->Update();
246     movingNormalizeFilter->Update();
247
248     //We keep the images used for input for possible output
249     inputFixedImage= fixedNormalizeFilter->GetOutput();
250     inputMovingImage= movingNormalizeFilter->GetOutput();
251
252     //the pointers are reconnected for further output
253     fixedImage=fixedNormalizeFilter->GetOutput();
254     movingImage=movingNormalizeFilter->GetOutput();
255
256     if (m_Verbose)  std::cout <<  "Normalizing image intensities to zero mean and SD of 1..." << std::endl;
257   }
258
259
260   //If given, the images are blurred before processing
261   if ( m_ArgsInfo.blur_arg!= 0) {
262     typedef itk::DiscreteGaussianImageFilter<InternalImageType,InternalImageType> GaussianFilterType;
263     typename GaussianFilterType::Pointer fixedSmoother = GaussianFilterType::New();
264     typename GaussianFilterType::Pointer movingSmoother = GaussianFilterType::New();
265     fixedSmoother->SetVariance( m_ArgsInfo.blur_arg );
266     movingSmoother->SetVariance(m_ArgsInfo.blur_arg );
267
268     fixedSmoother->SetInput( fixedImage );
269     movingSmoother->SetInput( movingImage );
270
271     fixedSmoother->Update();
272     movingSmoother->Update();
273
274     fixedImage=fixedSmoother->GetOutput();
275     movingImage=movingSmoother->GetOutput();
276
277     if (m_Verbose)  std::cout <<  "Blurring images with a Gaussian with standard deviation of " << m_ArgsInfo.blur_arg <<"..." << std::endl;
278   }
279
280
281   //============================================================================
282   // Setting up the moving image in a reference system
283   //============================================================================
284   const itk::Vector<double, InputImageType::ImageDimension> movingResolution = movingImage->GetSpacing();
285   typename InternalImageType::RegionType movingRegion = movingImage->GetLargestPossibleRegion();
286   typename InternalImageType::RegionType::SizeType  movingSize = movingRegion.GetSize();
287
288   // Print the parameters of the moving image
289   if (m_Verbose) {
290     std::cout << "Object or Moving image:"<<std::endl;
291     std::cout << "Size: " << movingSize[0] << ", " << movingSize[1];
292     if (InputImageType::ImageDimension==3) std::cout<<", " << movingSize[2];
293     std::cout << std::endl;
294
295     std::cout<< "Resolution: "<< movingResolution[0] << ", " << movingResolution[1];
296     if (InputImageType::ImageDimension==3) std::cout<< ", " << movingResolution[2];
297     std::cout << std::endl;
298   }
299
300
301   //============================================================================
302   // Setting up the fixed image in a reference system
303   //============================================================================
304   const itk::Vector<double, InputImageType::ImageDimension> fixedResolution = fixedImage->GetSpacing();
305   typename InternalImageType::RegionType fixedRegion = fixedImage->GetLargestPossibleRegion();
306   typename InternalImageType::RegionType::SizeType fixedSize = fixedRegion.GetSize();
307
308   // Print the parameters of the moving image and the transform
309   if (m_Verbose) {
310     std::cout << "Target or Moving image:"<<std::endl;
311     std::cout << "Size: " << fixedSize[0] << ", " << fixedSize[1];
312     if (InputImageType::ImageDimension==3) std::cout<<", " << fixedSize[2];
313     std::cout << std::endl;
314
315     std::cout<< "Resolution: "<< fixedResolution[0] << ", " << fixedResolution[1];
316     if (InputImageType::ImageDimension==3) std::cout<< ", " << fixedResolution[2];
317     std::cout << std::endl;
318   }
319
320
321
322   //===========================================================================
323   // If given, we connect a mask to reference or target
324   //============================================================================
325   typedef itk::ImageMaskSpatialObject<  InputImageType::ImageDimension >   MaskType;
326   typename MaskType::Pointer  fixedMask=NULL;
327   if (m_ArgsInfo.referenceMask_given) {
328     fixedMask= MaskType::New();
329     typedef itk::Image< unsigned char, InputImageType::ImageDimension >   ImageMaskType;
330     typedef itk::ImageFileReader< ImageMaskType >    MaskReaderType;
331     typename MaskReaderType::Pointer  maskReader = MaskReaderType::New();
332     maskReader->SetFileName(m_ArgsInfo.referenceMask_arg);
333     try {
334       maskReader->Update();
335     } catch ( itk::ExceptionObject & err ) {
336       std::cerr << "ExceptionObject caught while reading mask !" << std::endl;
337       std::cerr << err << std::endl;
338       return;
339     }
340     if (m_Verbose)std::cout <<"Reference image mask was read..." <<std::endl;
341
342     // Set the image to the spatialObject
343     fixedMask->SetImage( maskReader->GetOutput() );
344   }
345
346   typedef itk::ImageMaskSpatialObject<  InputImageType::ImageDimension >   MaskType;
347   typename MaskType::Pointer  movingMask=NULL;
348   if (m_ArgsInfo.targetMask_given) {
349     movingMask= MaskType::New();
350     typedef itk::Image< unsigned char, InputImageType::ImageDimension >   ImageMaskType;
351     typedef itk::ImageFileReader< ImageMaskType >    MaskReaderType;
352     typename MaskReaderType::Pointer  maskReader = MaskReaderType::New();
353     maskReader->SetFileName(m_ArgsInfo.targetMask_arg);
354     try {
355       maskReader->Update();
356     } catch ( itk::ExceptionObject & err ) {
357       std::cerr << "ExceptionObject caught !" << std::endl;
358       std::cerr << err << std::endl;
359     }
360     if (m_Verbose)std::cout <<"Target image mask was read..." <<std::endl;
361
362     movingMask->SetImage( maskReader->GetOutput() );
363   }
364
365
366   //============================================================================
367   // The image pyramids
368   //============================================================================
369   typedef itk::RecursiveMultiResolutionPyramidImageFilter<InternalImageType,InternalImageType >  FixedImagePyramidType;
370   typedef itk::RecursiveMultiResolutionPyramidImageFilter<InternalImageType,InternalImageType >  MovingImagePyramidType;
371   typename FixedImagePyramidType::Pointer fixedImagePyramid = FixedImagePyramidType::New();
372   typename MovingImagePyramidType::Pointer movingImagePyramid = MovingImagePyramidType::New();
373   fixedImagePyramid->SetUseShrinkImageFilter(false);
374   fixedImagePyramid->SetInput(fixedImage);
375   fixedImagePyramid->SetNumberOfLevels(m_ArgsInfo.levels_arg);
376   movingImagePyramid->SetUseShrinkImageFilter(false);
377   movingImagePyramid->SetInput(movingImage);
378   movingImagePyramid->SetNumberOfLevels(m_ArgsInfo.levels_arg);
379   if (m_Verbose) std::cout<<"Creating the image pyramid..."<<std::endl;
380   fixedImagePyramid->Update();
381   movingImagePyramid->Update();
382
383
384
385   //============================================================================
386   // We retrieve the type of metric from the command line
387   //============================================================================
388   typedef clitk::GenericMetric<args_info_clitkAffineRegistration, InternalImageType, InternalImageType> GenericMetricType;
389   typename GenericMetricType::Pointer genericMetric=GenericMetricType::New();
390   genericMetric->SetArgsInfo(m_ArgsInfo);
391   genericMetric->SetFixedImage(fixedImagePyramid->GetOutput(0));
392   if (fixedMask) genericMetric->SetFixedImageMask(fixedMask);
393   typedef itk::ImageToImageMetric< InternalImageType, InternalImageType >  MetricType;
394   typename  MetricType::Pointer metric=genericMetric->GetMetricPointer();
395   if (movingMask) metric->SetMovingImageMask(movingMask);
396
397 #ifdef ITK_USE_OPTIMIZED_REGISTRATION_METHODS
398   if (threadsGiven) metric->SetNumberOfThreads( threads );
399 #else
400   if (m_Verbose) std::cout<<"Not setting the number of threads (not compiled with USE_OPTIMIZED_REGISTRATION_METHODS)..."<<std::endl;
401 #endif
402
403   //============================================================================
404   // Initialize using image moments.
405   //============================================================================
406   if (m_ArgsInfo.moment_flag) {
407     typedef itk::ImageMomentsCalculator< InternalImageType > CalculatorType;
408     typename CalculatorType::Pointer fixedCalculator= CalculatorType::New();
409
410     typename InternalImageType::Pointer fixedThresh;
411     if (m_ArgsInfo.intThreshold_given) {
412       typedef itk::ThresholdImageFilter<InternalImageType> ThresholdImageFilterType;
413       typename ThresholdImageFilterType::Pointer thresholder = ThresholdImageFilterType::New();
414       thresholder->SetInput(fixedImage);
415       thresholder->SetLower(m_ArgsInfo.intThreshold_arg);
416       thresholder->Update();
417       fixedThresh=thresholder->GetOutput();
418     } else fixedThresh=fixedImage;
419
420     fixedCalculator->SetImage(fixedThresh);
421     fixedCalculator->Compute();
422     Vector<double, InputImageType::ImageDimension> fixedCenter=fixedCalculator->GetCenterOfGravity();
423     if (m_Verbose)std::cout<<"The fixed center of gravity is "<<fixedCenter<<"..."<<std::endl;
424
425     typedef itk::ImageMomentsCalculator< InternalImageType > CalculatorType;
426     typename CalculatorType::Pointer movingCalculator= CalculatorType::New();
427
428     typename InternalImageType::Pointer movingThresh;
429     if (m_ArgsInfo.intThreshold_given) {
430       typedef itk::ThresholdImageFilter<InternalImageType> ThresholdImageFilterType;
431       typename ThresholdImageFilterType::Pointer thresholder = ThresholdImageFilterType::New();
432       thresholder->SetInput(movingImage);
433       thresholder->SetLower(m_ArgsInfo.intThreshold_arg);
434       thresholder->Update();
435       movingThresh=thresholder->GetOutput();
436     } else movingThresh=movingImage;
437
438     movingCalculator->SetImage(movingThresh);
439     movingCalculator->Compute();
440     Vector<double, InputImageType::ImageDimension> movingCenter=movingCalculator->GetCenterOfGravity();
441     if (m_Verbose)std::cout<<"The moving center of gravity is "<<movingCenter<<"..."<<std::endl;
442
443     Vector<double, InputImageType::ImageDimension> shift= movingCenter-fixedCenter;
444     if (m_Verbose)std::cout<<"The initial shift applied is "<<shift<<"..."<<std::endl;
445
446     m_ArgsInfo.transX_arg= shift [0];
447     m_ArgsInfo.transY_arg= shift [1];
448     if (InputImageType::ImageDimension==3) m_ArgsInfo.transZ_arg=shift [2];
449   }
450
451   //============================================================================
452   // Transform
453   //============================================================================
454   typedef clitk::GenericAffineTransform<args_info_clitkAffineRegistration, TCoordRep, InputImageType::ImageDimension > GenericAffineTransformType;
455   typename GenericAffineTransformType::Pointer genericAffineTransform = GenericAffineTransformType::New();
456   genericAffineTransform->SetArgsInfo(m_ArgsInfo);
457   typedef itk::Transform< double, InputImageType::ImageDimension, InputImageType::ImageDimension > TransformType;
458   typename TransformType::Pointer transform = genericAffineTransform->GetTransform();
459
460
461   //=======================================================
462   // Interpolator
463   //=======================================================
464   typedef clitk::GenericInterpolator<args_info_clitkAffineRegistration, InternalImageType,TCoordRep > GenericInterpolatorType;
465   typename GenericInterpolatorType::Pointer genericInterpolator=GenericInterpolatorType::New();
466   genericInterpolator->SetArgsInfo(m_ArgsInfo);
467   typedef itk::InterpolateImageFunction< InternalImageType, TCoordRep >  InterpolatorType;
468   typename  InterpolatorType::Pointer interpolator=genericInterpolator->GetInterpolatorPointer();
469
470
471   //============================================================================
472   // Optimizer
473   //============================================================================
474   typedef clitk::GenericOptimizer<args_info_clitkAffineRegistration> GenericOptimizerType;
475   unsigned int nParam = transform->GetNumberOfParameters();
476   typename GenericOptimizerType::Pointer genericOptimizer=GenericOptimizerType::New();
477   genericOptimizer->SetArgsInfo(m_ArgsInfo);
478   genericOptimizer->SetOutputIteration(m_Verbose);
479   genericOptimizer->SetOutputPosition(m_Verbose);
480   genericOptimizer->SetOutputValue(m_Verbose);
481   genericOptimizer->SetOutputGradient(m_ArgsInfo.gradient_flag);
482   genericOptimizer->SetMaximize(genericMetric->GetMaximize());
483   genericOptimizer->SetNumberOfParameters(nParam);
484   typedef itk::SingleValuedNonLinearOptimizer OptimizerType;
485   OptimizerType::Pointer optimizer=genericOptimizer->GetOptimizerPointer();
486
487   // Scales
488   itk::Optimizer::ScalesType scales( nParam );
489   for (unsigned int i=nParam-InputImageType::ImageDimension; i<nParam; i++) //Translations
490     scales[i] = m_ArgsInfo.tWeight_arg;
491   for (unsigned int i=0; i<nParam-InputImageType::ImageDimension; i++)      //Rest
492     scales[i] = m_ArgsInfo.rWeight_arg*180./M_PI;
493   optimizer->SetScales(scales);
494   //============================================================================
495   // Multiresolution registration
496   //============================================================================
497   typedef itk::MultiResolutionImageRegistrationMethod< InternalImageType,InternalImageType >  RegistrationType;
498   typename  RegistrationType::Pointer registration = RegistrationType::New();
499   registration->SetFixedImage( fixedImage  );
500   registration->SetFixedImageRegion(fixedImage->GetLargestPossibleRegion());
501   registration->SetMovingImage(  movingImage );
502   registration->SetFixedImagePyramid( fixedImagePyramid );
503   registration->SetMovingImagePyramid( movingImagePyramid );
504   registration->SetTransform( transform );
505   registration->SetInitialTransformParameters( transform->GetParameters() );
506   registration->SetInterpolator( interpolator );
507   registration->SetMetric(metric);
508   registration->SetOptimizer(optimizer);
509   registration->SetNumberOfLevels( m_ArgsInfo.levels_arg );
510   if (m_Verbose) std::cout << "Setting "<< m_ArgsInfo.levels_arg <<" resolution levels..." << std::endl;
511   if (m_Verbose) std::cout << "Initial Transform: "<< registration->GetInitialTransformParameters()<<std::endl;
512
513   //============================================================================
514   // Connecting the commander to the registration to monitor it
515   //============================================================================
516   if (m_Verbose) {
517
518     // Output iteration info
519     CommandIterationUpdate::Pointer observer = CommandIterationUpdate::New();
520     observer->SetOptimizer(genericOptimizer);
521     optimizer->AddObserver( itk::IterationEvent(), observer );
522
523
524     // Output level info
525     typedef RegistrationInterfaceCommand<RegistrationType, args_info_clitkAffineRegistration> CommandType;
526     typename CommandType::Pointer command = CommandType::New();
527     command->SetArgsInfo(m_ArgsInfo);
528     registration->AddObserver( itk::IterationEvent(), command );
529
530   }
531
532
533   //============================================================================
534   // Finally we can start the registration with the given amount of multiresolution levels
535   //============================================================================
536   if (m_Verbose) std::cout << "Starting the registration now..." << std::endl;
537
538   try {
539     registration->StartRegistration();
540   } catch ( itk::ExceptionObject & err ) {
541     std::cerr << "ExceptionObject caught !" << std::endl;
542     std::cerr << err << std::endl;
543   }
544
545
546   //============================================================================
547   // Processing the result of the registration
548   //============================================================================
549   OptimizerType::ParametersType finalParameters = registration->GetLastTransformParameters();
550   std::cout<< "Result : " <<std::setprecision(12)<<std::endl;
551
552   for (unsigned int i=nParam-InputImageType::ImageDimension; i<nParam; i++) //Translations
553     std::cout << " Translation " << i << " = " << finalParameters[i];
554   for (unsigned int i=0; i<nParam-InputImageType::ImageDimension; i++)      //Rest
555     std::cout << " Other parameter " << i << " = " << finalParameters[i];
556
557
558   itk::Matrix<double,InputImageType::ImageDimension+1,InputImageType::ImageDimension+1> matrix;
559   if (m_ArgsInfo.transform_arg == 3) {
560     for (unsigned int i=0; i<InputImageType::ImageDimension; i++) {
561       matrix[i][3] =  finalParameters[nParam-InputImageType::ImageDimension+i];
562       for (unsigned int j=0; j<InputImageType::ImageDimension; j++) {
563         matrix[i][j] = finalParameters[i*3+j];
564       }
565       matrix[3][3] = 1.0;
566     }
567   } else {
568     matrix = clitk::GetBackwardAffineMatrix<InputImageType::ImageDimension>(finalParameters);
569   }
570
571   std::cout << " Affine transform matrix =" << std::endl;
572   std::cout << matrix <<std::setprecision(6)<< std::endl;
573
574   // Write matrix to a file
575   if (m_ArgsInfo.matrix_given) {
576     std::ofstream mFile;
577     mFile.open(m_ArgsInfo.matrix_arg);
578     mFile<<std::setprecision(12)<<matrix<< std::setprecision(6)<<std::endl;
579     mFile.close();
580   }
581
582   //============================================================================
583   // Prepare the resampling filter in order to transform the moving image.
584   //============================================================================
585   if (m_ArgsInfo.output_given || m_ArgsInfo.checker_after_given || m_ArgsInfo.after_given ) {
586     transform->SetParameters( finalParameters );
587     typedef itk::ResampleImageFilter< InternalImageType,InternalImageType >    ResampleFilterType;
588     typename    ResampleFilterType::Pointer resampler = ResampleFilterType::New();
589
590     resampler->SetTransform( transform );
591     resampler->SetInput( movingImage );
592     resampler->SetSize( fixedImage->GetLargestPossibleRegion().GetSize() );
593     resampler->SetOutputOrigin(  fixedImage->GetOrigin() );
594     resampler->SetOutputSpacing( fixedImage->GetSpacing() );
595     resampler->SetDefaultPixelValue( 0 );
596
597     //Output?
598     if (m_ArgsInfo.output_given) {
599       //We write an output in the same pixeltype then the input
600       typedef itk::ImageFileWriter< FixedImageType >  WriterType;
601       typename WriterType::Pointer outputWriter =  WriterType::New();
602       outputWriter->SetFileName(m_ArgsInfo.output_arg );
603       outputWriter->SetInput( resampler->GetOutput()   );
604       outputWriter->Update();
605     }
606
607
608     //============================================================================
609     // Checker after?
610     //============================================================================
611     if (m_ArgsInfo.checker_after_given) {
612       //To display correctly the checkerboard image, the intensities must lie in the same range (normalized)
613       //We write the image in the internal image type
614       typedef itk::ResampleImageFilter< InternalImageType,InternalImageType >    ResampleFilterType;
615       typename    ResampleFilterType::Pointer internalResampler = ResampleFilterType::New();
616       internalResampler->SetTransform( transform );
617       internalResampler->SetInput( inputMovingImage );
618       internalResampler->SetSize( fixedImage->GetLargestPossibleRegion().GetSize() );
619       internalResampler->SetOutputOrigin(  fixedImage->GetOrigin() );
620       internalResampler->SetOutputSpacing( fixedImage->GetSpacing() );
621       internalResampler->SetDefaultPixelValue( 0 );
622
623       //We pass the normalized images to the checker filter
624       typedef itk::CheckerBoardImageFilter< InternalImageType > CheckerBoardFilterType;
625       typename CheckerBoardFilterType::Pointer checkerFilter= CheckerBoardFilterType::New();
626
627       checkerFilter->SetInput1(inputFixedImage);
628       checkerFilter->SetInput2(internalResampler->GetOutput());
629       typedef itk::ImageFileWriter< InternalImageType >  InternalWriterType;
630       typename  InternalWriterType::Pointer checkerWriter =  InternalWriterType::New();
631       checkerWriter->SetFileName(m_ArgsInfo.checker_after_arg);
632       checkerWriter->SetInput( checkerFilter->GetOutput() );
633       checkerWriter->Update();
634     }
635
636
637     //============================================================================
638     // Checker before?
639     //============================================================================
640     if (m_ArgsInfo.checker_before_given) {
641       //To display correctly the checkerboard image, the intensities must lie in the same range (normalized)
642       //We write the image in the internal image type
643       //We pass the normalized images to the checker filter
644       typedef itk::CheckerBoardImageFilter< InternalImageType > CheckerBoardFilterType;
645       typename CheckerBoardFilterType::Pointer checkerFilter= CheckerBoardFilterType::New();
646
647       checkerFilter->SetInput1(inputFixedImage);
648       checkerFilter->SetInput2(inputMovingImage);
649       typedef itk::ImageFileWriter< InternalImageType >  InternalWriterType;
650       typename  InternalWriterType::Pointer checkerWriter =  InternalWriterType::New();
651       checkerWriter->SetFileName(m_ArgsInfo.checker_before_arg);
652       checkerWriter->SetInput( checkerFilter->GetOutput() );
653       checkerWriter->Update();
654     }
655
656
657     //============================================================================
658     // Difference After?
659     //============================================================================
660     if (m_ArgsInfo.after_given) {
661       typedef itk::SubtractImageFilter< InternalImageType, FixedImageType,FixedImageType > DifferenceImageFilterType;
662       typename DifferenceImageFilterType::Pointer differenceAfterFilter= DifferenceImageFilterType::New();
663
664       differenceAfterFilter->SetInput1(fixedImage);
665       differenceAfterFilter->SetInput2(resampler->GetOutput());
666
667       // Prepare a writer to write the difference image
668       typedef itk::ImageFileWriter< FixedImageType >  WriterType;
669       typename   WriterType::Pointer     differenceAfterWriter =  WriterType::New();
670       differenceAfterWriter->SetFileName(m_ArgsInfo.after_arg );
671       differenceAfterWriter->SetInput( differenceAfterFilter->GetOutput()   );
672       differenceAfterWriter->Update();
673     }
674   }
675
676   //============================================================================
677   // Difference Before?
678   //============================================================================
679   if (m_ArgsInfo.before_given) {
680     typedef itk::CastImageFilter< InternalImageType,FixedImageType > CastFilterType;
681     typename    CastFilterType::Pointer  caster =  CastFilterType::New();
682     caster->SetInput( movingImage );
683
684     typedef itk::SubtractImageFilter< InternalImageType, FixedImageType, FixedImageType > DifferenceImageFilterType;
685     typename DifferenceImageFilterType::Pointer differenceBeforeFilter= DifferenceImageFilterType::New();
686
687
688     differenceBeforeFilter->SetInput1(fixedImage);
689     differenceBeforeFilter->SetInput2(caster->GetOutput());
690
691     // Prepare a writer to write the difference image
692     typedef itk::ImageFileWriter< FixedImageType >  WriterType;
693     typename WriterType::Pointer     differenceBeforeWriter =  WriterType::New();
694     differenceBeforeWriter->SetFileName(m_ArgsInfo.before_arg);
695     differenceBeforeWriter->SetInput( differenceBeforeFilter->GetOutput()   );
696     differenceBeforeWriter->Update();
697   }
698
699 }
700 }
701 #endif  //#define  CLITKAFFINEREGISTRATIONGENERICFILTER_TXX
702