]> Creatis software - clitk.git/blob - registration/clitkAffineRegistrationGenericFilter.txx
*** empty log message ***
[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 // Update with the number of dimensions and pixeltype
154 //==============================================================================
155 template<class args_info_clitkAffineRegistration>
156 template<class InputImageType>
157 void AffineRegistrationGenericFilter<args_info_clitkAffineRegistration>::UpdateWithInputImageType()
158 {
159   //=============================================================================
160   //Input
161   //=============================================================================
162
163   typedef typename  InputImageType::PixelType PixelType;
164 //typedef typename InputImageType::ImageDimension Dimension;
165
166
167   bool threadsGiven=m_ArgsInfo.threads_given;
168   int threads=m_ArgsInfo.threads_arg;
169
170   //Coordinate Representation
171   typedef double TCoordRep;
172
173   //The pixeltype of the fixed image will be used for output
174   typedef itk::Image< PixelType, InputImageType::ImageDimension > FixedImageType;
175
176   //Whatever the pixel type, internally we work with an image represented in float
177   typedef typename  InputImageType::PixelType  InternalPixelType;
178   typedef itk::Image< InternalPixelType, InputImageType::ImageDimension > InternalImageType;
179
180   //Read in the reference/fixed image
181   typedef itk::ImageFileReader< InternalImageType > ReaderType;
182   typename ReaderType::Pointer  fixedImageReader  = ReaderType::New();
183   fixedImageReader->SetFileName( m_ArgsInfo.reference_arg);
184
185
186   //Read in the object/moving image
187   typename ReaderType::Pointer movingImageReader = ReaderType::New();
188   movingImageReader->SetFileName( m_ArgsInfo.target_arg );
189   if (m_Verbose) std::cout<<"Reading images..."<<std::endl;
190   fixedImageReader->Update();
191   movingImageReader->Update();
192
193   if (m_Verbose) std::cout  << "Reading images... " << std::endl;
194
195   //we connect pointers to these internal images
196   typename InternalImageType::Pointer fixedImage= fixedImageReader->GetOutput();
197   typename InternalImageType::Pointer movingImage= movingImageReader->GetOutput();
198
199   //We keep the images used for input for possible output
200   typename InternalImageType::Pointer inputFixedImage= fixedImageReader->GetOutput();
201   typename InternalImageType::Pointer inputMovingImage= movingImageReader->GetOutput();
202
203
204   //============================================================================
205   // Preprocessing
206   //============================================================================
207
208   //If given, the intensities of both images are first normalized to a zero mean and SD of 1
209   // (usefull for MI, not necessary for Mattes' MI but performed anyway for the ouput)
210   if ( m_ArgsInfo.normalize_flag ) {
211     typedef itk::NormalizeImageFilter< InternalImageType,InternalImageType >  NormalizeFilterType;
212
213     typename  NormalizeFilterType::Pointer  fixedNormalizeFilter = NormalizeFilterType::New();
214     typename  NormalizeFilterType::Pointer  movingNormalizeFilter = NormalizeFilterType::New();
215
216     fixedNormalizeFilter->SetInput( fixedImage );
217     movingNormalizeFilter->SetInput( movingImage );
218
219     fixedNormalizeFilter->Update();
220     movingNormalizeFilter->Update();
221
222     //We keep the images used for input for possible output
223     inputFixedImage= fixedNormalizeFilter->GetOutput();
224     inputMovingImage= movingNormalizeFilter->GetOutput();
225
226     //the pointers are reconnected for further output
227     fixedImage=fixedNormalizeFilter->GetOutput();
228     movingImage=movingNormalizeFilter->GetOutput();
229
230     if (m_Verbose)  std::cout <<  "Normalizing image intensities to zero mean and SD of 1..." << std::endl;
231   }
232
233
234   //If given, the images are blurred before processing
235   if ( m_ArgsInfo.blur_arg!= 0) {
236     typedef itk::DiscreteGaussianImageFilter<InternalImageType,InternalImageType> GaussianFilterType;
237     typename GaussianFilterType::Pointer fixedSmoother = GaussianFilterType::New();
238     typename GaussianFilterType::Pointer movingSmoother = GaussianFilterType::New();
239     fixedSmoother->SetVariance( m_ArgsInfo.blur_arg );
240     movingSmoother->SetVariance(m_ArgsInfo.blur_arg );
241
242     fixedSmoother->SetInput( fixedImage );
243     movingSmoother->SetInput( movingImage );
244
245     fixedSmoother->Update();
246     movingSmoother->Update();
247
248     fixedImage=fixedSmoother->GetOutput();
249     movingImage=movingSmoother->GetOutput();
250
251     if (m_Verbose)  std::cout <<  "Blurring images with a Gaussian with standard deviation of " << m_ArgsInfo.blur_arg <<"..." << std::endl;
252   }
253
254
255   //============================================================================
256   // Setting up the moving image in a reference system
257   //============================================================================
258   const itk::Vector<double, InputImageType::ImageDimension> movingResolution = movingImage->GetSpacing();
259   typename InternalImageType::RegionType movingRegion = movingImage->GetLargestPossibleRegion();
260   typename InternalImageType::RegionType::SizeType  movingSize = movingRegion.GetSize();
261
262   // Print the parameters of the moving image
263   if (m_Verbose) {
264     std::cout << "Object or Moving image:"<<std::endl;
265     std::cout << "Size: " << movingSize[0] << ", " << movingSize[1];
266     if (InputImageType::ImageDimension==3) std::cout<<", " << movingSize[2];
267     std::cout << std::endl;
268
269     std::cout<< "Resolution: "<< movingResolution[0] << ", " << movingResolution[1];
270     if (InputImageType::ImageDimension==3) std::cout<< ", " << movingResolution[2];
271     std::cout << std::endl;
272   }
273
274
275   //============================================================================
276   // Setting up the fixed image in a reference system
277   //============================================================================
278   const itk::Vector<double, InputImageType::ImageDimension> fixedResolution = fixedImage->GetSpacing();
279   typename InternalImageType::RegionType fixedRegion = fixedImage->GetLargestPossibleRegion();
280   typename InternalImageType::RegionType::SizeType fixedSize = fixedRegion.GetSize();
281
282   // Print the parameters of the moving image and the transform
283   if (m_Verbose) {
284     std::cout << "Target or Moving image:"<<std::endl;
285     std::cout << "Size: " << fixedSize[0] << ", " << fixedSize[1];
286     if (InputImageType::ImageDimension==3) std::cout<<", " << fixedSize[2];
287     std::cout << std::endl;
288
289     std::cout<< "Resolution: "<< fixedResolution[0] << ", " << fixedResolution[1];
290     if (InputImageType::ImageDimension==3) std::cout<< ", " << fixedResolution[2];
291     std::cout << std::endl;
292   }
293
294
295
296   //===========================================================================
297   // If given, we connect a mask to reference or target
298   //============================================================================
299   typedef itk::ImageMaskSpatialObject<  InputImageType::ImageDimension >   MaskType;
300   typename MaskType::Pointer  fixedMask=NULL;
301   if (m_ArgsInfo.referenceMask_given) {
302     fixedMask= MaskType::New();
303     typedef itk::Image< unsigned char, InputImageType::ImageDimension >   ImageMaskType;
304     typedef itk::ImageFileReader< ImageMaskType >    MaskReaderType;
305     typename MaskReaderType::Pointer  maskReader = MaskReaderType::New();
306     maskReader->SetFileName(m_ArgsInfo.referenceMask_arg);
307     try {
308       maskReader->Update();
309     } catch ( itk::ExceptionObject & err ) {
310       std::cerr << "ExceptionObject caught while reading mask !" << std::endl;
311       std::cerr << err << std::endl;
312       return;
313     }
314     if (m_Verbose)std::cout <<"Reference image mask was read..." <<std::endl;
315
316     // Set the image to the spatialObject
317     fixedMask->SetImage( maskReader->GetOutput() );
318   }
319
320   typedef itk::ImageMaskSpatialObject<  InputImageType::ImageDimension >   MaskType;
321   typename MaskType::Pointer  movingMask=NULL;
322   if (m_ArgsInfo.targetMask_given) {
323     movingMask= MaskType::New();
324     typedef itk::Image< unsigned char, InputImageType::ImageDimension >   ImageMaskType;
325     typedef itk::ImageFileReader< ImageMaskType >    MaskReaderType;
326     typename MaskReaderType::Pointer  maskReader = MaskReaderType::New();
327     maskReader->SetFileName(m_ArgsInfo.targetMask_arg);
328     try {
329       maskReader->Update();
330     } catch ( itk::ExceptionObject & err ) {
331       std::cerr << "ExceptionObject caught !" << std::endl;
332       std::cerr << err << std::endl;
333     }
334     if (m_Verbose)std::cout <<"Target image mask was read..." <<std::endl;
335
336     movingMask->SetImage( maskReader->GetOutput() );
337   }
338
339
340   //============================================================================
341   // The image pyramids
342   //============================================================================
343   typedef itk::RecursiveMultiResolutionPyramidImageFilter<InternalImageType,InternalImageType >  FixedImagePyramidType;
344   typedef itk::RecursiveMultiResolutionPyramidImageFilter<InternalImageType,InternalImageType >  MovingImagePyramidType;
345   typename FixedImagePyramidType::Pointer fixedImagePyramid = FixedImagePyramidType::New();
346   typename MovingImagePyramidType::Pointer movingImagePyramid = MovingImagePyramidType::New();
347   fixedImagePyramid->SetUseShrinkImageFilter(false);
348   fixedImagePyramid->SetInput(fixedImage);
349   fixedImagePyramid->SetNumberOfLevels(m_ArgsInfo.levels_arg);
350   movingImagePyramid->SetUseShrinkImageFilter(false);
351   movingImagePyramid->SetInput(movingImage);
352   movingImagePyramid->SetNumberOfLevels(m_ArgsInfo.levels_arg);
353   if (m_Verbose) std::cout<<"Creating the image pyramid..."<<std::endl;
354   fixedImagePyramid->Update();
355   movingImagePyramid->Update();
356
357
358
359   //============================================================================
360   // We retrieve the type of metric from the command line
361   //============================================================================
362   typedef clitk::GenericMetric<args_info_clitkAffineRegistration, InternalImageType, InternalImageType> GenericMetricType;
363   typename GenericMetricType::Pointer genericMetric=GenericMetricType::New();
364   genericMetric->SetArgsInfo(m_ArgsInfo);
365   genericMetric->SetFixedImage(fixedImagePyramid->GetOutput(0));
366   if (fixedMask) genericMetric->SetFixedImageMask(fixedMask);
367   typedef itk::ImageToImageMetric< InternalImageType, InternalImageType >  MetricType;
368   typename  MetricType::Pointer metric=genericMetric->GetMetricPointer();
369   if (movingMask) metric->SetMovingImageMask(movingMask);
370
371 #ifdef ITK_USE_OPTIMIZED_REGISTRATION_METHODS
372   if (threadsGiven) metric->SetNumberOfThreads( threads );
373 #else
374   if (m_Verbose) std::cout<<"Not setting the number of threads (not compiled with USE_OPTIMIZED_REGISTRATION_METHODS)..."<<std::endl;
375 #endif
376
377   //============================================================================
378   // Initialize using image moments.
379   //============================================================================
380   if (m_ArgsInfo.moment_flag) {
381     typedef itk::ImageMomentsCalculator< InternalImageType > CalculatorType;
382     typename CalculatorType::Pointer fixedCalculator= CalculatorType::New();
383
384     typename InternalImageType::Pointer fixedThresh;
385     if (m_ArgsInfo.intThreshold_given) {
386       typedef itk::ThresholdImageFilter<InternalImageType> ThresholdImageFilterType;
387       typename ThresholdImageFilterType::Pointer thresholder = ThresholdImageFilterType::New();
388       thresholder->SetInput(fixedImage);
389       thresholder->SetLower(m_ArgsInfo.intThreshold_arg);
390       thresholder->Update();
391       fixedThresh=thresholder->GetOutput();
392     } else fixedThresh=fixedImage;
393
394     fixedCalculator->SetImage(fixedThresh);
395     fixedCalculator->Compute();
396     Vector<double, InputImageType::ImageDimension> fixedCenter=fixedCalculator->GetCenterOfGravity();
397     if (m_Verbose)std::cout<<"The fixed center of gravity is "<<fixedCenter<<"..."<<std::endl;
398
399     typedef itk::ImageMomentsCalculator< InternalImageType > CalculatorType;
400     typename CalculatorType::Pointer movingCalculator= CalculatorType::New();
401
402     typename InternalImageType::Pointer movingThresh;
403     if (m_ArgsInfo.intThreshold_given) {
404       typedef itk::ThresholdImageFilter<InternalImageType> ThresholdImageFilterType;
405       typename ThresholdImageFilterType::Pointer thresholder = ThresholdImageFilterType::New();
406       thresholder->SetInput(movingImage);
407       thresholder->SetLower(m_ArgsInfo.intThreshold_arg);
408       thresholder->Update();
409       movingThresh=thresholder->GetOutput();
410     } else movingThresh=movingImage;
411
412     movingCalculator->SetImage(movingThresh);
413     movingCalculator->Compute();
414     Vector<double, InputImageType::ImageDimension> movingCenter=movingCalculator->GetCenterOfGravity();
415     if (m_Verbose)std::cout<<"The moving center of gravity is "<<movingCenter<<"..."<<std::endl;
416
417     Vector<double, InputImageType::ImageDimension> shift= movingCenter-fixedCenter;
418     if (m_Verbose)std::cout<<"The initial shift applied is "<<shift<<"..."<<std::endl;
419
420     m_ArgsInfo.transX_arg= shift [0];
421     m_ArgsInfo.transY_arg= shift [1];
422     if (InputImageType::ImageDimension==3) m_ArgsInfo.transZ_arg=shift [2];
423   }
424
425   //============================================================================
426   // Transform
427   //============================================================================
428   typedef clitk::GenericAffineTransform<args_info_clitkAffineRegistration, TCoordRep, InputImageType::ImageDimension > GenericAffineTransformType;
429   typename GenericAffineTransformType::Pointer genericAffineTransform = GenericAffineTransformType::New();
430   genericAffineTransform->SetArgsInfo(m_ArgsInfo);
431   typedef itk::Transform< double, InputImageType::ImageDimension, InputImageType::ImageDimension > TransformType;
432   typename TransformType::Pointer transform = genericAffineTransform->GetTransform();
433
434
435   //=======================================================
436   // Interpolator
437   //=======================================================
438   typedef clitk::GenericInterpolator<args_info_clitkAffineRegistration, InternalImageType,TCoordRep > GenericInterpolatorType;
439   typename GenericInterpolatorType::Pointer genericInterpolator=GenericInterpolatorType::New();
440   genericInterpolator->SetArgsInfo(m_ArgsInfo);
441   typedef itk::InterpolateImageFunction< InternalImageType, TCoordRep >  InterpolatorType;
442   typename  InterpolatorType::Pointer interpolator=genericInterpolator->GetInterpolatorPointer();
443
444
445   //============================================================================
446   // Optimizer
447   //============================================================================
448   typedef clitk::GenericOptimizer<args_info_clitkAffineRegistration> GenericOptimizerType;
449   unsigned int nParam = transform->GetNumberOfParameters();
450   typename GenericOptimizerType::Pointer genericOptimizer=GenericOptimizerType::New();
451   genericOptimizer->SetArgsInfo(m_ArgsInfo);
452   genericOptimizer->SetOutputIteration(m_Verbose);
453   genericOptimizer->SetOutputPosition(m_Verbose);
454   genericOptimizer->SetOutputValue(m_Verbose);
455   genericOptimizer->SetOutputGradient(m_ArgsInfo.gradient_flag);
456   genericOptimizer->SetMaximize(genericMetric->GetMaximize());
457   genericOptimizer->SetNumberOfParameters(nParam);
458   typedef itk::SingleValuedNonLinearOptimizer OptimizerType;
459   OptimizerType::Pointer optimizer=genericOptimizer->GetOptimizerPointer();
460
461   // Scales
462   itk::Optimizer::ScalesType scales( nParam );
463   for (unsigned int i=nParam-InputImageType::ImageDimension; i<nParam; i++) //Translations
464     scales[i] = m_ArgsInfo.tWeight_arg;
465   for (unsigned int i=0; i<nParam-InputImageType::ImageDimension; i++)      //Rest
466     scales[i] = m_ArgsInfo.rWeight_arg*180./M_PI;
467   optimizer->SetScales(scales);
468   //============================================================================
469   // Multiresolution registration
470   //============================================================================
471   typedef itk::MultiResolutionImageRegistrationMethod< InternalImageType,InternalImageType >  RegistrationType;
472   typename  RegistrationType::Pointer registration = RegistrationType::New();
473   registration->SetFixedImage( fixedImage  );
474   registration->SetFixedImageRegion(fixedImage->GetLargestPossibleRegion());
475   registration->SetMovingImage(  movingImage );
476   registration->SetFixedImagePyramid( fixedImagePyramid );
477   registration->SetMovingImagePyramid( movingImagePyramid );
478   registration->SetTransform( transform );
479   registration->SetInitialTransformParameters( transform->GetParameters() );
480   registration->SetInterpolator( interpolator );
481   registration->SetMetric(metric);
482   registration->SetOptimizer(optimizer);
483   registration->SetNumberOfLevels( m_ArgsInfo.levels_arg );
484   if (m_Verbose) std::cout << "Setting "<< m_ArgsInfo.levels_arg <<" resolution levels..." << std::endl;
485   if (m_Verbose) std::cout << "Initial Transform: "<< registration->GetInitialTransformParameters()<<std::endl;
486
487   //============================================================================
488   // Connecting the commander to the registration to monitor it
489   //============================================================================
490   if (m_Verbose) {
491
492     // Output iteration info
493     CommandIterationUpdate::Pointer observer = CommandIterationUpdate::New();
494     observer->SetOptimizer(genericOptimizer);
495     optimizer->AddObserver( itk::IterationEvent(), observer );
496
497
498     // Output level info
499     typedef RegistrationInterfaceCommand<RegistrationType, args_info_clitkAffineRegistration> CommandType;
500     typename CommandType::Pointer command = CommandType::New();
501     command->SetArgsInfo(m_ArgsInfo);
502     registration->AddObserver( itk::IterationEvent(), command );
503
504   }
505
506
507   //============================================================================
508   // Finally we can start the registration with the given amount of multiresolution levels
509   //============================================================================
510   if (m_Verbose) std::cout << "Starting the registration now..." << std::endl;
511
512   try {
513     registration->StartRegistration();
514   } catch ( itk::ExceptionObject & err ) {
515     std::cerr << "ExceptionObject caught !" << std::endl;
516     std::cerr << err << std::endl;
517   }
518
519
520   //============================================================================
521   // Processing the result of the registration
522   //============================================================================
523   OptimizerType::ParametersType finalParameters = registration->GetLastTransformParameters();
524   std::cout<< "Result : " <<std::setprecision(12)<<std::endl;
525
526   for (unsigned int i=nParam-InputImageType::ImageDimension; i<nParam; i++) //Translations
527     std::cout << " Translation " << i << " = " << finalParameters[i];
528   for (unsigned int i=0; i<nParam-InputImageType::ImageDimension; i++)      //Rest
529     std::cout << " Other parameter " << i << " = " << finalParameters[i];
530
531
532   itk::Matrix<double,InputImageType::ImageDimension+1,InputImageType::ImageDimension+1> matrix;
533   if (m_ArgsInfo.transform_arg == 3) {
534     for (unsigned int i=0; i<InputImageType::ImageDimension; i++) {
535       matrix[i][3] =  finalParameters[nParam-InputImageType::ImageDimension+i];
536       for (unsigned int j=0; j<InputImageType::ImageDimension; j++) {
537         matrix[i][j] = finalParameters[i*3+j];
538       }
539       matrix[3][3] = 1.0;
540     }
541   } else {
542     matrix = clitk::GetBackwardAffineMatrix<InputImageType::ImageDimension>(finalParameters);
543   }
544
545   std::cout << " Affine transform matrix =" << std::endl;
546   std::cout << matrix <<std::setprecision(6)<< std::endl;
547
548   // Write matrix to a file
549   if (m_ArgsInfo.matrix_given) {
550     std::ofstream mFile;
551     mFile.open(m_ArgsInfo.matrix_arg);
552     mFile<<std::setprecision(12)<<matrix<< std::setprecision(6)<<std::endl;
553     mFile.close();
554   }
555
556   //============================================================================
557   // Prepare the resampling filter in order to transform the moving image.
558   //============================================================================
559   if (m_ArgsInfo.output_given || m_ArgsInfo.checker_after_given || m_ArgsInfo.after_given ) {
560     transform->SetParameters( finalParameters );
561     typedef itk::ResampleImageFilter< InternalImageType,InternalImageType >    ResampleFilterType;
562     typename    ResampleFilterType::Pointer resampler = ResampleFilterType::New();
563
564     resampler->SetTransform( transform );
565     resampler->SetInput( movingImageReader->GetOutput() );
566     resampler->SetSize( fixedImageReader ->GetOutput()->GetLargestPossibleRegion().GetSize() );
567     resampler->SetOutputOrigin(  fixedImageReader ->GetOutput()->GetOrigin() );
568     resampler->SetOutputSpacing( fixedImageReader ->GetOutput()->GetSpacing() );
569     resampler->SetDefaultPixelValue( 0 );
570
571     //Output?
572     if (m_ArgsInfo.output_given) {
573       //We write an output in the same pixeltype then the input
574       typedef itk::ImageFileWriter< FixedImageType >  WriterType;
575       typename WriterType::Pointer outputWriter =  WriterType::New();
576       outputWriter->SetFileName(m_ArgsInfo.output_arg );
577       outputWriter->SetInput( resampler->GetOutput()   );
578       outputWriter->Update();
579     }
580
581
582     //============================================================================
583     // Checker after?
584     //============================================================================
585     if (m_ArgsInfo.checker_after_given) {
586       //To display correctly the checkerboard image, the intensities must lie in the same range (normalized)
587       //We write the image in the internal image type
588       typedef itk::ResampleImageFilter< InternalImageType,InternalImageType >    ResampleFilterType;
589       typename    ResampleFilterType::Pointer internalResampler = ResampleFilterType::New();
590       internalResampler->SetTransform( transform );
591       internalResampler->SetInput( inputMovingImage );
592       internalResampler->SetSize( fixedImage->GetLargestPossibleRegion().GetSize() );
593       internalResampler->SetOutputOrigin(  fixedImage->GetOrigin() );
594       internalResampler->SetOutputSpacing( fixedImage->GetSpacing() );
595       internalResampler->SetDefaultPixelValue( 0 );
596
597       //We pass the normalized images to the checker filter
598       typedef itk::CheckerBoardImageFilter< InternalImageType > CheckerBoardFilterType;
599       typename CheckerBoardFilterType::Pointer checkerFilter= CheckerBoardFilterType::New();
600
601       checkerFilter->SetInput1(inputFixedImage);
602       checkerFilter->SetInput2(internalResampler->GetOutput());
603       typedef itk::ImageFileWriter< InternalImageType >  InternalWriterType;
604       typename  InternalWriterType::Pointer checkerWriter =  InternalWriterType::New();
605       checkerWriter->SetFileName(m_ArgsInfo.checker_after_arg);
606       checkerWriter->SetInput( checkerFilter->GetOutput() );
607       checkerWriter->Update();
608     }
609
610
611     //============================================================================
612     // Checker before?
613     //============================================================================
614     if (m_ArgsInfo.checker_before_given) {
615       //To display correctly the checkerboard image, the intensities must lie in the same range (normalized)
616       //We write the image in the internal image type
617       //We pass the normalized images to the checker filter
618       typedef itk::CheckerBoardImageFilter< InternalImageType > CheckerBoardFilterType;
619       typename CheckerBoardFilterType::Pointer checkerFilter= CheckerBoardFilterType::New();
620
621       checkerFilter->SetInput1(inputFixedImage);
622       checkerFilter->SetInput2(inputMovingImage);
623       typedef itk::ImageFileWriter< InternalImageType >  InternalWriterType;
624       typename  InternalWriterType::Pointer checkerWriter =  InternalWriterType::New();
625       checkerWriter->SetFileName(m_ArgsInfo.checker_before_arg);
626       checkerWriter->SetInput( checkerFilter->GetOutput() );
627       checkerWriter->Update();
628     }
629
630
631     //============================================================================
632     // Difference After?
633     //============================================================================
634     if (m_ArgsInfo.after_given) {
635       typedef itk::SubtractImageFilter< InternalImageType, FixedImageType,FixedImageType > DifferenceImageFilterType;
636       typename DifferenceImageFilterType::Pointer differenceAfterFilter= DifferenceImageFilterType::New();
637
638       differenceAfterFilter->SetInput1(fixedImageReader ->GetOutput());
639       differenceAfterFilter->SetInput2(resampler->GetOutput());
640
641       // Prepare a writer to write the difference image
642       typedef itk::ImageFileWriter< FixedImageType >  WriterType;
643       typename   WriterType::Pointer     differenceAfterWriter =  WriterType::New();
644       differenceAfterWriter->SetFileName(m_ArgsInfo.after_arg );
645       differenceAfterWriter->SetInput( differenceAfterFilter->GetOutput()   );
646       differenceAfterWriter->Update();
647     }
648   }
649
650   //============================================================================
651   // Difference Before?
652   //============================================================================
653   if (m_ArgsInfo.before_given) {
654     typedef itk::CastImageFilter< InternalImageType,FixedImageType > CastFilterType;
655     typename    CastFilterType::Pointer  caster =  CastFilterType::New();
656     caster->SetInput( movingImageReader->GetOutput() );
657
658     typedef itk::SubtractImageFilter< InternalImageType, FixedImageType, FixedImageType > DifferenceImageFilterType;
659     typename DifferenceImageFilterType::Pointer differenceBeforeFilter= DifferenceImageFilterType::New();
660
661
662     differenceBeforeFilter->SetInput1(fixedImageReader ->GetOutput());
663     differenceBeforeFilter->SetInput2(caster->GetOutput());
664
665     // Prepare a writer to write the difference image
666     typedef itk::ImageFileWriter< FixedImageType >  WriterType;
667     typename WriterType::Pointer     differenceBeforeWriter =  WriterType::New();
668     differenceBeforeWriter->SetFileName(m_ArgsInfo.before_arg);
669     differenceBeforeWriter->SetInput( differenceBeforeFilter->GetOutput()   );
670     differenceBeforeWriter->Update();
671   }
672
673 }
674 }
675 #endif  //#define  CLITKAFFINEREGISTRATIONGENERICFILTER_TXX
676