]> Creatis software - clitk.git/blob - registration/clitkBLUTDIRGenericFilter.cxx
Use itkLabelGeometryImageFilter instead of itkLabelStatisticsImageFilter it
[clitk.git] / registration / clitkBLUTDIRGenericFilter.cxx
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 #ifndef clitkBLUTDIRGenericFilter_cxx
19 #define clitkBLUTDIRGenericFilter_cxx
20
21 /* =================================================
22  * @file   clitkBLUTDIRGenericFilter.cxx
23  * @author
24  * @date
25  *
26  * @brief
27  *
28  ===================================================*/
29
30 #include "clitkBLUTDIRGenericFilter.h"
31
32 namespace clitk
33 {
34
35   //==============================================================================
36   // Creating an observer class that allows output at each iteration
37   //==============================================================================
38   class CommandIterationUpdate : public itk::Command
39   {
40     public:
41       typedef  CommandIterationUpdate   Self;
42       typedef  itk::Command             Superclass;
43       typedef  itk::SmartPointer<Self>  Pointer;
44       itkNewMacro( Self );
45     protected:
46       CommandIterationUpdate() {};
47     public:
48       typedef   clitk::GenericOptimizer<args_info_clitkBLUTDIR>     OptimizerType;
49       typedef   const OptimizerType   *           OptimizerPointer;
50
51       // Set the generic optimizer
52       void SetOptimizer(OptimizerPointer o){m_Optimizer=o;}
53
54       // Execute
55       void Execute(itk::Object *caller, const itk::EventObject & event)
56       {
57         Execute( (const itk::Object *)caller, event);
58       }
59
60       void Execute(const itk::Object * object, const itk::EventObject & event)
61       {
62         if( !(itk::IterationEvent().CheckEvent( &event )) )
63         {
64           return;
65         }
66
67         m_Optimizer->OutputIterationInfo();
68       }
69
70       OptimizerPointer m_Optimizer;
71   };
72
73   //===========================================================================//
74   //Constructor
75   //==========================================================================//
76   BLUTDIRGenericFilter::BLUTDIRGenericFilter():
77     ImageToImageGenericFilter<Self>("Register DIR")
78   {
79     InitializeImageType<2>();
80     InitializeImageType<3>();
81     m_Verbose=true;
82   }
83
84   //=========================================================================//
85   //SetArgsInfo
86   //==========================================================================//
87   void BLUTDIRGenericFilter::SetArgsInfo(const args_info_clitkBLUTDIR & a){
88     m_ArgsInfo=a;
89     if (m_ArgsInfo.reference_given) AddInputFilename(m_ArgsInfo.reference_arg);
90
91     if (m_ArgsInfo.target_given) {
92       AddInputFilename(m_ArgsInfo.target_arg);
93     }
94
95     if (m_ArgsInfo.output_given) SetOutputFilename(m_ArgsInfo.output_arg);
96   }
97
98   //=========================================================================//
99   //===========================================================================//
100   template<unsigned int Dim>
101     void BLUTDIRGenericFilter::InitializeImageType()
102     {
103       ADD_DEFAULT_IMAGE_TYPES(3);
104     }
105   //--------------------------------------------------------------------
106
107   //==============================================================================
108   //Creating an observer class that allows us to change parameters at subsequent levels
109   //==============================================================================
110   template <typename TRegistration,class args_info_clitkBLUTDIR>
111     class RegistrationInterfaceCommand : public itk::Command
112   {
113     public:
114       typedef RegistrationInterfaceCommand   Self;
115       typedef itk::Command             Superclass;
116       typedef itk::SmartPointer<Self>  Pointer;
117       itkNewMacro( Self );
118     protected:
119       RegistrationInterfaceCommand() { };
120     public:
121
122       // Registration
123       typedef   TRegistration                              RegistrationType;
124       typedef   RegistrationType *                         RegistrationPointer;
125
126       // Transform
127       typedef typename RegistrationType::FixedImageType FixedImageType;
128       typedef typename FixedImageType::RegionType RegionType;
129       itkStaticConstMacro(ImageDimension, unsigned int,FixedImageType::ImageDimension);
130       typedef clitk::BSplineDeformableTransform<double, ImageDimension, ImageDimension> TransformType;
131       typedef clitk::BSplineDeformableTransformInitializer<TransformType, FixedImageType> InitializerType;
132       typedef typename InitializerType::CoefficientImageType CoefficientImageType;
133       typedef itk::CastImageFilter<CoefficientImageType, CoefficientImageType> CastImageFilterType;
134       typedef typename TransformType::ParametersType ParametersType;
135       typedef typename InitializerType::Pointer InitializerPointer;
136       typedef   CommandIterationUpdate::Pointer CommandIterationUpdatePointer;
137
138       // Optimizer
139       typedef clitk::GenericOptimizer<args_info_clitkBLUTDIR> GenericOptimizerType;
140       typedef typename GenericOptimizerType::Pointer GenericOptimizerPointer;
141
142       // Metric
143       typedef typename RegistrationType::FixedImageType    InternalImageType;
144       typedef clitk::GenericMetric<args_info_clitkBLUTDIR, InternalImageType, InternalImageType> GenericMetricType;
145       typedef typename GenericMetricType::Pointer GenericMetricPointer;
146
147       // Two arguments are passed to the Execute() method: the first
148       // is the pointer to the object which invoked the event and the
149       // second is the event that was invoked.
150       void Execute(itk::Object * object, const itk::EventObject & event)
151       {
152         if( !(itk::IterationEvent().CheckEvent( &event )) )
153         {
154           return;
155         }
156
157         // Get the levels
158         RegistrationPointer registration = dynamic_cast<RegistrationPointer>( object );
159         unsigned int numberOfLevels=registration->GetNumberOfLevels();
160         unsigned int currentLevel=registration->GetCurrentLevel()+1;
161
162         // Output the levels
163         std::cout<<std::endl;
164         std::cout<<"========================================"<<std::endl;
165         std::cout<<"Starting resolution level "<<currentLevel<<" of "<<numberOfLevels<<"..."<<std::endl;
166         std::cout<<"========================================"<<std::endl;
167         std::cout<<std::endl;
168
169         // Higher level?
170         if (currentLevel>1)
171         {
172           // fixed image region pyramid
173           typedef clitk::MultiResolutionPyramidRegionFilter<InternalImageType> FixedImageRegionPyramidType;
174           typename FixedImageRegionPyramidType::Pointer fixedImageRegionPyramid=FixedImageRegionPyramidType::New();
175           fixedImageRegionPyramid->SetRegion(m_MetricRegion);
176           fixedImageRegionPyramid->SetSchedule(registration->GetFixedImagePyramid()->GetSchedule());
177
178           // Reinitialize the metric (!= number of samples)
179           m_GenericMetric= GenericMetricType::New();
180           m_GenericMetric->SetArgsInfo(m_ArgsInfo);
181           m_GenericMetric->SetFixedImage(registration->GetFixedImagePyramid()->GetOutput(registration->GetCurrentLevel()));
182           if (m_ArgsInfo.referenceMask_given)  m_GenericMetric->SetFixedImageMask(registration->GetMetric()->GetFixedImageMask());
183           m_GenericMetric->SetFixedImageRegion(fixedImageRegionPyramid->GetOutput(registration->GetCurrentLevel()));
184           typedef itk::ImageToImageMetric< InternalImageType, InternalImageType >  MetricType;
185           typename  MetricType::Pointer metric=m_GenericMetric->GetMetricPointer();
186           registration->SetMetric(metric);
187
188           // Get the current coefficient image and make a COPY
189           typename itk::ImageDuplicator<CoefficientImageType>::Pointer caster=itk::ImageDuplicator<CoefficientImageType>::New();
190           caster->SetInputImage(m_Initializer->GetTransform()->GetCoefficientImage());
191           caster->Update();
192           typename CoefficientImageType::Pointer currentCoefficientImage=caster->GetOutput();
193
194           // Write the intermediate result?
195           if (m_ArgsInfo.intermediate_given>=numberOfLevels)
196             writeImage<CoefficientImageType>(currentCoefficientImage, m_ArgsInfo.intermediate_arg[currentLevel-2], m_ArgsInfo.verbose_flag);
197
198           // Set the new transform properties
199           m_Initializer->SetImage(registration->GetFixedImagePyramid()->GetOutput(currentLevel-1));
200           if( m_Initializer->m_ControlPointSpacingIsGiven)
201             m_Initializer->SetControlPointSpacing(m_Initializer->m_ControlPointSpacingArray[registration->GetCurrentLevel()]);
202           if( m_Initializer->m_NumberOfControlPointsIsGiven)
203             m_Initializer->SetNumberOfControlPointsInsideTheImage(m_Initializer->m_NumberOfControlPointsInsideTheImageArray[registration->GetCurrentLevel()]);
204
205           // Reinitialize the transform
206           if (m_ArgsInfo.verbose_flag) std::cout<<"Initializing transform for level "<<currentLevel<<" of "<<numberOfLevels<<"..."<<std::endl;
207           m_Initializer->InitializeTransform();
208           ParametersType* newParameters= new typename TransformType::ParametersType(m_Initializer->GetTransform()->GetNumberOfParameters());
209
210           // DS : if we want to skip the last pyramid level, force to only 1 iteration
211           DD(m_ArgsInfo.skipLastPyramidLevel_flag);
212           if ((currentLevel == numberOfLevels) && (m_ArgsInfo.skipLastPyramidLevel_flag)) {
213             DD(m_ArgsInfo.maxIt_arg);
214             std::cout << "I skip the last pyramid level : set max iteration to 0" << std::endl;
215             m_ArgsInfo.maxIt_arg = 0;
216             DD(m_ArgsInfo.maxIt_arg);
217           }
218
219           // Reinitialize an Optimizer (!= number of parameters)
220           m_GenericOptimizer = GenericOptimizerType::New();
221           m_GenericOptimizer->SetArgsInfo(m_ArgsInfo);
222           m_GenericOptimizer->SetMaximize(m_Maximize);
223           m_GenericOptimizer->SetNumberOfParameters(m_Initializer->GetTransform()->GetNumberOfParameters());
224
225
226           typedef itk::SingleValuedNonLinearOptimizer OptimizerType;
227           OptimizerType::Pointer optimizer = m_GenericOptimizer->GetOptimizerPointer();
228           optimizer->AddObserver( itk::IterationEvent(), m_CommandIterationUpdate);
229           registration->SetOptimizer(optimizer);
230           m_CommandIterationUpdate->SetOptimizer(m_GenericOptimizer);
231
232           // Set the previous transform parameters to the registration
233           // if(m_Initializer->m_Parameters!=NULL )delete m_Initializer->m_Parameters;
234           m_Initializer->SetInitialParameters(currentCoefficientImage,*newParameters);
235           registration->SetInitialTransformParametersOfNextLevel(*newParameters);
236         }
237       }
238
239       void Execute(const itk::Object * , const itk::EventObject & )
240       { return; }
241
242
243       // Members
244       void SetInitializer(InitializerPointer i){m_Initializer=i;}
245       InitializerPointer m_Initializer;
246
247       void SetArgsInfo(args_info_clitkBLUTDIR a){m_ArgsInfo=a;}
248       args_info_clitkBLUTDIR m_ArgsInfo;
249
250       void SetCommandIterationUpdate(CommandIterationUpdatePointer c){m_CommandIterationUpdate=c;};
251       CommandIterationUpdatePointer m_CommandIterationUpdate;
252
253       GenericOptimizerPointer m_GenericOptimizer;
254       void SetMaximize(bool b){m_Maximize=b;}
255       bool m_Maximize;
256
257       GenericMetricPointer m_GenericMetric;
258       void SetMetricRegion(RegionType i){m_MetricRegion=i;}
259       RegionType m_MetricRegion;
260
261
262   };
263
264   //==============================================================================
265   // Update with the number of dimensions and pixeltype
266   //==============================================================================
267   template<class InputImageType>
268     void BLUTDIRGenericFilter::UpdateWithInputImageType()
269     {
270       //=============================================================================
271       //Input
272       //=============================================================================
273       bool threadsGiven=m_ArgsInfo.threads_given;
274       int threads=m_ArgsInfo.threads_arg;
275       typedef typename InputImageType::PixelType PixelType;
276
277       typedef double TCoordRep;
278
279       typename InputImageType::Pointer fixedImage = this->template GetInput<InputImageType>(0);
280
281       typename InputImageType::Pointer inputFixedImage = this->template GetInput<InputImageType>(0);
282
283       // typedef input2
284       typename InputImageType::Pointer movingImage = this->template GetInput<InputImageType>(1);
285
286       typename InputImageType::Pointer inputMovingImage = this->template GetInput<InputImageType>(1);
287
288       typedef itk::Image< PixelType,InputImageType::ImageDimension >  FixedImageType;
289       typedef itk::Image< PixelType, InputImageType::ImageDimension>  MovingImageType;
290       const unsigned int SpaceDimension = InputImageType::ImageDimension;
291       //Whatever the pixel type, internally we work with an image represented in float
292       //Reading reference image
293       if (m_Verbose) std::cout<<"Reading images..."<<std::endl;
294       //=======================================================
295       //Input
296       //=======================================================
297       typename FixedImageType::Pointer croppedFixedImage=fixedImage;
298       //=======================================================
299       // Regions
300       //=======================================================
301       // The original input region
302       typename FixedImageType::RegionType fixedImageRegion = fixedImage->GetLargestPossibleRegion();
303
304       // The transform region with respect to the input region:
305       // where should the transform be DEFINED (depends on mask)
306       typename FixedImageType::RegionType transformRegion = fixedImage->GetLargestPossibleRegion();
307       typename FixedImageType::RegionType::SizeType transformRegionSize=transformRegion.GetSize();
308       typename FixedImageType::RegionType::IndexType transformRegionIndex=transformRegion.GetIndex();
309       typename FixedImageType::PointType transformRegionOrigin=fixedImage->GetOrigin();
310
311       // The metric region with respect to the extracted transform region:
312       // where should the metric be CALCULATED (depends on transform)
313       typename FixedImageType::RegionType metricRegion = fixedImage->GetLargestPossibleRegion();
314       typename FixedImageType::RegionType::SizeType metricRegionSize=metricRegion.GetSize();
315       typename FixedImageType::RegionType::IndexType metricRegionIndex=metricRegion.GetIndex();
316       typename FixedImageType::PointType metricRegionOrigin=fixedImage->GetOrigin();
317
318
319       //===========================================================================
320       // If given, we connect a mask to reference or target
321       //============================================================================
322       typedef itk::ImageMaskSpatialObject< InputImageType::ImageDimension >   MaskType;
323       typename MaskType::Pointer  fixedMask=NULL;
324       if (m_ArgsInfo.referenceMask_given)
325       {
326         fixedMask= MaskType::New();
327         typedef itk::Image< unsigned char,InputImageType::ImageDimension >   ImageMaskType;
328         typedef itk::ImageFileReader< ImageMaskType >    MaskReaderType;
329         typename MaskReaderType::Pointer  maskReader = MaskReaderType::New();
330         maskReader->SetFileName(m_ArgsInfo.referenceMask_arg);
331         try
332         {
333           maskReader->Update();
334         }
335         catch( itk::ExceptionObject & err )
336         {
337           std::cerr << "ExceptionObject caught while reading mask !" << std::endl;
338           std::cerr << err << std::endl;
339           return;
340         }
341         if (m_Verbose)std::cout <<"Reference image mask was read..." <<std::endl;
342
343         // Set the image to the spatialObject
344         fixedMask->SetImage( maskReader->GetOutput() );
345
346         // Find the bounding box of the "inside" label
347         typedef itk::LabelGeometryImageFilter<ImageMaskType> GeometryImageFilterType;
348         typename GeometryImageFilterType::Pointer geometryImageFilter=GeometryImageFilterType::New();
349         geometryImageFilter->SetInput(maskReader->GetOutput());
350         geometryImageFilter->Update();
351         typename GeometryImageFilterType::BoundingBoxType boundingBox = geometryImageFilter->GetBoundingBox(1);
352
353         // Limit the transform region to the mask
354         for (unsigned int i=0; i<InputImageType::ImageDimension; i++)
355         {
356           transformRegionIndex[i]=boundingBox[2*i];
357           transformRegionSize[i]=boundingBox[2*i+1]-boundingBox[2*i]+1;
358         }
359         transformRegion.SetSize(transformRegionSize);
360         transformRegion.SetIndex(transformRegionIndex);
361         fixedImage->TransformIndexToPhysicalPoint(transformRegion.GetIndex(), transformRegionOrigin);
362
363         // Crop the fixedImage to the bounding box to facilitate multi-resolution
364         typedef itk::ExtractImageFilter<FixedImageType,FixedImageType> ExtractImageFilterType;
365         typename ExtractImageFilterType::Pointer extractImageFilter=ExtractImageFilterType::New();
366         extractImageFilter->SetInput(fixedImage);
367         extractImageFilter->SetExtractionRegion(transformRegion);
368         extractImageFilter->Update();
369         croppedFixedImage=extractImageFilter->GetOutput();
370
371         // Update the metric region
372         metricRegion = croppedFixedImage->GetLargestPossibleRegion();
373         metricRegionIndex=metricRegion.GetIndex();
374         croppedFixedImage->TransformIndexToPhysicalPoint(metricRegionIndex, metricRegionOrigin);
375
376         // Set start index to zero (with respect to croppedFixedImage/transform region)
377         metricRegionIndex.Fill(0);
378         metricRegion.SetIndex(metricRegionIndex);
379         croppedFixedImage->SetRegions(metricRegion);
380         croppedFixedImage->SetOrigin(metricRegionOrigin);
381       }
382
383       typedef itk::ImageMaskSpatialObject< InputImageType::ImageDimension >   MaskType;
384       typename MaskType::Pointer  movingMask=NULL;
385       if (m_ArgsInfo.targetMask_given)
386       {
387         movingMask= MaskType::New();
388         typedef itk::Image< unsigned char, InputImageType::ImageDimension >   ImageMaskType;
389         typedef itk::ImageFileReader< ImageMaskType >    MaskReaderType;
390         typename MaskReaderType::Pointer  maskReader = MaskReaderType::New();
391         maskReader->SetFileName(m_ArgsInfo.targetMask_arg);
392         try
393         {
394           maskReader->Update();
395         }
396         catch( itk::ExceptionObject & err )
397         {
398           std::cerr << "ExceptionObject caught !" << std::endl;
399           std::cerr << err << std::endl;
400         }
401         if (m_Verbose)std::cout <<"Target image mask was read..." <<std::endl;
402
403         movingMask->SetImage( maskReader->GetOutput() );
404       }
405
406
407       //=======================================================
408       // Output Regions
409       //=======================================================
410
411       if (m_Verbose)
412       {
413         // Fixed image region
414         std::cout<<"The fixed image has its origin at "<<fixedImage->GetOrigin()<<std::endl
415           <<"The fixed image region starts at index "<<fixedImageRegion.GetIndex()<<std::endl
416           <<"The fixed image region has size "<< fixedImageRegion.GetSize()<<std::endl;
417
418         // Transform region
419         std::cout<<"The transform has its origin at "<<transformRegionOrigin<<std::endl
420           <<"The transform region will start at index "<<transformRegion.GetIndex()<<std::endl
421           <<"The transform region has size "<< transformRegion.GetSize()<<std::endl;
422
423         // Metric region
424         std::cout<<"The metric region has its origin at "<<metricRegionOrigin<<std::endl
425           <<"The metric region will start at index "<<metricRegion.GetIndex()<<std::endl
426           <<"The metric region has size "<< metricRegion.GetSize()<<std::endl;
427
428       }
429
430
431       //=======================================================
432       // Pyramids (update them for transform initializer)
433       //=======================================================
434       typedef itk::RecursiveMultiResolutionPyramidImageFilter< FixedImageType, FixedImageType>    FixedImagePyramidType;
435       typedef itk::RecursiveMultiResolutionPyramidImageFilter< MovingImageType, MovingImageType>    MovingImagePyramidType;
436       typename FixedImagePyramidType::Pointer fixedImagePyramid = FixedImagePyramidType::New();
437       typename MovingImagePyramidType::Pointer movingImagePyramid = MovingImagePyramidType::New();
438       fixedImagePyramid->SetUseShrinkImageFilter(false);
439       fixedImagePyramid->SetInput(croppedFixedImage);
440       fixedImagePyramid->SetNumberOfLevels(m_ArgsInfo.levels_arg);
441       movingImagePyramid->SetUseShrinkImageFilter(false);
442       movingImagePyramid->SetInput(movingImage);
443       movingImagePyramid->SetNumberOfLevels(m_ArgsInfo.levels_arg);
444       if (m_Verbose) std::cout<<"Creating the image pyramid..."<<std::endl;
445       fixedImagePyramid->Update();
446       movingImagePyramid->Update();
447       typedef clitk::MultiResolutionPyramidRegionFilter<FixedImageType> FixedImageRegionPyramidType;
448       typename FixedImageRegionPyramidType::Pointer fixedImageRegionPyramid=FixedImageRegionPyramidType::New();
449       fixedImageRegionPyramid->SetRegion(metricRegion);
450       fixedImageRegionPyramid->SetSchedule(fixedImagePyramid->GetSchedule());
451
452
453       //=======================================================
454       // Rigid or Affine Transform
455       //=======================================================
456       typedef itk::AffineTransform <double,3> RigidTransformType;
457       RigidTransformType::Pointer rigidTransform=NULL;
458       if (m_ArgsInfo.initMatrix_given)
459       {
460         if(m_Verbose) std::cout<<"Reading the prior transform matrix "<< m_ArgsInfo.initMatrix_arg<<"..."<<std::endl;
461         rigidTransform=RigidTransformType::New();
462         itk::Matrix<double,4,4> rigidTransformMatrix=clitk::ReadMatrix3D(m_ArgsInfo.initMatrix_arg);
463
464         //Set the rotation
465         itk::Matrix<double,3,3> finalRotation = clitk::GetRotationalPartMatrix3D(rigidTransformMatrix);
466         rigidTransform->SetMatrix(finalRotation);
467
468         //Set the translation
469         itk::Vector<double,3> finalTranslation = clitk::GetTranslationPartMatrix3D(rigidTransformMatrix);
470         rigidTransform->SetTranslation(finalTranslation);
471       }
472
473
474       //=======================================================
475       // B-LUT FFD Transform
476       //=======================================================
477       typedef  clitk::BSplineDeformableTransform<TCoordRep,InputImageType::ImageDimension, SpaceDimension > TransformType;
478       typename TransformType::Pointer transform= TransformType::New();
479       if (fixedMask) transform->SetMask( fixedMask );
480       if (rigidTransform) transform->SetBulkTransform( rigidTransform );
481
482       //-------------------------------------------------------------------------
483       // The transform initializer
484       //-------------------------------------------------------------------------
485       typedef clitk::BSplineDeformableTransformInitializer< TransformType,FixedImageType> InitializerType;
486       typename InitializerType::Pointer initializer = InitializerType::New();
487       initializer->SetVerbose(m_Verbose);
488       initializer->SetImage(fixedImagePyramid->GetOutput(0));
489       initializer->SetTransform(transform);
490
491       //-------------------------------------------------------------------------
492       // Order
493       //-------------------------------------------------------------------------
494       typename FixedImageType::RegionType::SizeType splineOrders ;
495       splineOrders.Fill(3);
496       if (m_ArgsInfo.order_given)
497         for(unsigned int i=0; i<InputImageType::ImageDimension;i++)
498           splineOrders[i]=m_ArgsInfo.order_arg[i];
499       if (m_Verbose) std::cout<<"Setting the spline orders  to "<<splineOrders<<"..."<<std::endl;
500       initializer->SetSplineOrders(splineOrders);
501
502       //-------------------------------------------------------------------------
503       // Levels
504       //-------------------------------------------------------------------------
505
506       // Spacing
507       if (m_ArgsInfo.spacing_given)
508       {
509         initializer->m_ControlPointSpacingArray.resize(m_ArgsInfo.levels_arg);
510         initializer->SetControlPointSpacing(m_ArgsInfo.spacing_arg);
511         initializer->m_ControlPointSpacingArray[m_ArgsInfo.levels_arg-1]=initializer->m_ControlPointSpacing;
512         if (m_Verbose) std::cout<<"Using a control point spacing of "<<initializer->m_ControlPointSpacingArray[m_ArgsInfo.levels_arg-1]
513           <<" at level "<<m_ArgsInfo.levels_arg<<" of "<<m_ArgsInfo.levels_arg<<"..."<<std::endl;
514
515         for (int i=1; i<m_ArgsInfo.levels_arg; i++ )
516         {
517           initializer->m_ControlPointSpacingArray[m_ArgsInfo.levels_arg-1-i]=initializer->m_ControlPointSpacingArray[m_ArgsInfo.levels_arg-i]*2;
518           if (m_Verbose) std::cout<<"Using a control point spacing of "<<initializer->m_ControlPointSpacingArray[m_ArgsInfo.levels_arg-1-i]
519             <<" at level "<<m_ArgsInfo.levels_arg-i<<" of "<<m_ArgsInfo.levels_arg<<"..."<<std::endl;
520         }
521
522       }
523
524       // Control
525       if (m_ArgsInfo.control_given)
526       {
527         initializer->m_NumberOfControlPointsInsideTheImageArray.resize(m_ArgsInfo.levels_arg);
528         initializer->SetNumberOfControlPointsInsideTheImage(m_ArgsInfo.control_arg);
529         initializer->m_NumberOfControlPointsInsideTheImageArray[m_ArgsInfo.levels_arg-1]=initializer->m_NumberOfControlPointsInsideTheImage;
530         if (m_Verbose) std::cout<<"Using "<< initializer->m_NumberOfControlPointsInsideTheImageArray[m_ArgsInfo.levels_arg-1]<<"control points inside the image"
531           <<" at level "<<m_ArgsInfo.levels_arg<<" of "<<m_ArgsInfo.levels_arg<<"..."<<std::endl;
532
533         for (int i=1; i<m_ArgsInfo.levels_arg; i++ )
534         {
535           for(unsigned int j=0;j<InputImageType::ImageDimension;j++)
536             initializer->m_NumberOfControlPointsInsideTheImageArray[m_ArgsInfo.levels_arg-1-i][j]=ceil ((double)initializer->m_NumberOfControlPointsInsideTheImageArray[m_ArgsInfo.levels_arg-i][j]/2.);
537           //        initializer->m_NumberOfControlPointsInsideTheImageArray[m_ArgsInfo.levels_arg-1-i]=ceil ((double)initializer->m_NumberOfControlPointsInsideTheImageArray[m_ArgsInfo.levels_arg-i]/2.);
538           if (m_Verbose) std::cout<<"Using "<< initializer->m_NumberOfControlPointsInsideTheImageArray[m_ArgsInfo.levels_arg-1-i]<<"control points inside the image"
539             <<" at level "<<m_ArgsInfo.levels_arg<<" of "<<m_ArgsInfo.levels_arg<<"..."<<std::endl;
540
541         }
542       }
543
544       // Inialize on the first level
545       if (m_ArgsInfo.verbose_flag) std::cout<<"Initializing transform for level 1 of "<<m_ArgsInfo.levels_arg<<"..."<<std::endl;
546       if (m_ArgsInfo.spacing_given) initializer->SetControlPointSpacing(        initializer->m_ControlPointSpacingArray[0]);
547       if (m_ArgsInfo.control_given) initializer->SetNumberOfControlPointsInsideTheImage(initializer->m_NumberOfControlPointsInsideTheImageArray[0]);
548       if (m_ArgsInfo.samplingFactor_given) initializer->SetSamplingFactors(m_ArgsInfo.samplingFactor_arg);
549
550       // Initialize
551       initializer->InitializeTransform();
552
553       //-------------------------------------------------------------------------
554       // Initial parameters (passed by reference)
555       //-------------------------------------------------------------------------
556       typedef typename TransformType::ParametersType     ParametersType;
557       const unsigned int numberOfParameters =    transform->GetNumberOfParameters();
558       ParametersType parameters(numberOfParameters);
559       parameters.Fill( 0.0 );
560       transform->SetParameters( parameters );
561       if (m_ArgsInfo.initCoeff_given) initializer->SetInitialParameters(m_ArgsInfo.initCoeff_arg, parameters);
562
563
564       //=======================================================
565       // Interpolator
566       //=======================================================
567       typedef clitk::GenericInterpolator<args_info_clitkBLUTDIR, FixedImageType, TCoordRep > GenericInterpolatorType;
568       typename   GenericInterpolatorType::Pointer genericInterpolator=GenericInterpolatorType::New();
569       genericInterpolator->SetArgsInfo(m_ArgsInfo);
570       typedef itk::InterpolateImageFunction< FixedImageType, TCoordRep >  InterpolatorType;
571       typename  InterpolatorType::Pointer interpolator=genericInterpolator->GetInterpolatorPointer();
572
573
574       //=======================================================
575       // Metric
576       //=======================================================
577       typedef clitk::GenericMetric<args_info_clitkBLUTDIR, FixedImageType,MovingImageType > GenericMetricType;
578       typename GenericMetricType::Pointer genericMetric=GenericMetricType::New();
579       genericMetric->SetArgsInfo(m_ArgsInfo);
580       genericMetric->SetFixedImage(fixedImagePyramid->GetOutput(0));
581       if (fixedMask) genericMetric->SetFixedImageMask(fixedMask);
582       genericMetric->SetFixedImageRegion(fixedImageRegionPyramid->GetOutput(0));
583       typedef itk::ImageToImageMetric< FixedImageType, MovingImageType >  MetricType;
584       typename  MetricType::Pointer metric=genericMetric->GetMetricPointer();
585       if (movingMask) metric->SetMovingImageMask(movingMask);
586
587 #ifdef ITK_USE_OPTIMIZED_REGISTRATION_METHODS
588       if (threadsGiven) {
589         metric->SetNumberOfThreads( threads );
590         if (m_Verbose) std::cout<< "Using " << threads << " threads." << std::endl;
591       }
592 #else
593       if (m_Verbose) std::cout<<"Not setting the number of threads (not compiled with USE_OPTIMIZED_REGISTRATION_METHODS)..."<<std::endl;
594 #endif
595
596
597       //=======================================================
598       // Optimizer
599       //=======================================================
600       typedef clitk::GenericOptimizer<args_info_clitkBLUTDIR> GenericOptimizerType;
601       GenericOptimizerType::Pointer genericOptimizer = GenericOptimizerType::New();
602       genericOptimizer->SetArgsInfo(m_ArgsInfo);
603       genericOptimizer->SetMaximize(genericMetric->GetMaximize());
604       genericOptimizer->SetNumberOfParameters(transform->GetNumberOfParameters());
605       typedef itk::SingleValuedNonLinearOptimizer OptimizerType;
606       OptimizerType::Pointer optimizer = genericOptimizer->GetOptimizerPointer();
607
608
609       //=======================================================
610       // Registration
611       //=======================================================
612       typedef itk::MultiResolutionImageRegistrationMethod<  FixedImageType, MovingImageType >    RegistrationType;
613       typename RegistrationType::Pointer   registration  = RegistrationType::New();
614       registration->SetMetric(        metric        );
615       registration->SetOptimizer(     optimizer     );
616       registration->SetInterpolator(  interpolator  );
617       registration->SetTransform (transform);
618       if(threadsGiven) {
619         registration->SetNumberOfThreads(threads);
620         if (m_Verbose) std::cout<< "Using " << threads << " threads." << std::endl;
621       }
622       registration->SetFixedImage(  croppedFixedImage   );
623       registration->SetMovingImage(  movingImage   );
624       registration->SetFixedImageRegion( metricRegion );
625       registration->SetFixedImagePyramid( fixedImagePyramid );
626       registration->SetMovingImagePyramid( movingImagePyramid );
627       registration->SetInitialTransformParameters( transform->GetParameters() );
628       registration->SetNumberOfLevels( m_ArgsInfo.levels_arg );
629       if (m_Verbose) std::cout<<"Setting the number of resolution levels to "<<m_ArgsInfo.levels_arg<<"..."<<std::endl;
630
631
632       //================================================================================================
633       // Observers
634       //================================================================================================
635       if (m_Verbose)
636       {
637         // Output iteration info
638         CommandIterationUpdate::Pointer observer = CommandIterationUpdate::New();
639         observer->SetOptimizer(genericOptimizer);
640         optimizer->AddObserver( itk::IterationEvent(), observer );
641
642         // Output level info
643         typedef RegistrationInterfaceCommand<RegistrationType,args_info_clitkBLUTDIR> CommandType;
644         typename CommandType::Pointer command = CommandType::New();
645         command->SetInitializer(initializer);
646         command->SetArgsInfo(m_ArgsInfo);
647         command->SetCommandIterationUpdate(observer);
648         command->SetMaximize(genericMetric->GetMaximize());
649         command->SetMetricRegion(metricRegion);
650         registration->AddObserver( itk::IterationEvent(), command );
651       }
652
653
654       //=======================================================
655       // Let's go
656       //=======================================================
657       if (m_Verbose) std::cout << std::endl << "Starting Registration" << std::endl;
658
659       try
660       {
661         registration->StartRegistration();
662       }
663       catch( itk::ExceptionObject & err )
664       {
665         std::cerr << "ExceptionObject caught while registering!" << std::endl;
666         std::cerr << err << std::endl;
667         return;
668       }
669
670
671       //=======================================================
672       // Get the result
673       //=======================================================
674       OptimizerType::ParametersType finalParameters =  registration->GetLastTransformParameters();
675       transform->SetParameters( finalParameters );
676       if (m_Verbose)
677       {
678         std::cout<<"Stop condition description: "
679           <<registration->GetOptimizer()->GetStopConditionDescription()<<std::endl;
680       }
681
682
683       //=======================================================
684       // Get the BSpline coefficient images and write them
685       //=======================================================
686       if (m_ArgsInfo.coeff_given)
687       {
688         typedef typename TransformType::CoefficientImageType CoefficientImageType;
689         typename CoefficientImageType::Pointer coefficientImage =transform->GetCoefficientImage();
690         typedef itk::ImageFileWriter<CoefficientImageType> CoeffWriterType;
691         typename CoeffWriterType::Pointer coeffWriter=CoeffWriterType::New();
692         coeffWriter->SetInput(coefficientImage);
693         coeffWriter->SetFileName(m_ArgsInfo.coeff_arg);
694         coeffWriter->Update();
695       }
696
697
698
699       //=======================================================
700       // Compute the DVF (only deformable transform)
701       //=======================================================
702       typedef itk::Vector< float, SpaceDimension >  DisplacementType;
703       typedef itk::Image< DisplacementType, InputImageType::ImageDimension >  DeformationFieldType;
704       typedef itk::TransformToDeformationFieldSource<DeformationFieldType, double> ConvertorType;
705       typename ConvertorType::Pointer filter= ConvertorType::New();
706       filter->SetNumberOfThreads(1);
707       transform->SetBulkTransform(NULL);
708       filter->SetTransform(transform);
709       filter->SetOutputParametersFromImage(fixedImage);
710       filter->Update();
711       typename DeformationFieldType::Pointer field = filter->GetOutput();
712
713
714       //=======================================================
715       // Write the DVF
716       //=======================================================
717       typedef itk::ImageFileWriter< DeformationFieldType >  FieldWriterType;
718       typename FieldWriterType::Pointer fieldWriter = FieldWriterType::New();
719       fieldWriter->SetFileName( m_ArgsInfo.vf_arg );
720       fieldWriter->SetInput( field );
721       try
722       {
723         fieldWriter->Update();
724       }
725       catch( itk::ExceptionObject & excp )
726       {
727         std::cerr << "Exception thrown writing the DVF" << std::endl;
728         std::cerr << excp << std::endl;
729         return;
730       }
731
732
733       //=======================================================
734       // Resample the moving image
735       //=======================================================
736       typedef itk::ResampleImageFilter< MovingImageType, FixedImageType >    ResampleFilterType;
737       typename ResampleFilterType::Pointer resampler = ResampleFilterType::New();
738       if (rigidTransform) transform->SetBulkTransform(rigidTransform);
739       resampler->SetTransform( transform );
740       resampler->SetInput( movingImage);
741       resampler->SetOutputParametersFromImage(fixedImage);
742       resampler->Update();
743       typename FixedImageType::Pointer result=resampler->GetOutput();
744
745       //     typedef itk::WarpImageFilter< MovingImageType, FixedImageType, DeformationFieldType >    WarpFilterType;
746       //     typename WarpFilterType::Pointer warp = WarpFilterType::New();
747
748       //     warp->SetDeformationField( field );
749       //     warp->SetInput( movingImageReader->GetOutput() );
750       //     warp->SetOutputOrigin(  fixedImage->GetOrigin() );
751       //     warp->SetOutputSpacing( fixedImage->GetSpacing() );
752       //     warp->SetOutputDirection( fixedImage->GetDirection() );
753       //     warp->SetEdgePaddingValue( 0.0 );
754       //     warp->Update();
755
756
757       //=======================================================
758       // Write the warped image
759       //=======================================================
760       typedef itk::ImageFileWriter< FixedImageType >  WriterType;
761       typename WriterType::Pointer      writer =  WriterType::New();
762       writer->SetFileName( m_ArgsInfo.output_arg );
763       writer->SetInput( result    );
764
765       try
766       {
767         writer->Update();
768       }
769       catch( itk::ExceptionObject & err )
770       {
771         std::cerr << "ExceptionObject caught !" << std::endl;
772         std::cerr << err << std::endl;
773         return;
774       }
775
776
777       //=======================================================
778       // Calculate the difference after the deformable transform
779       //=======================================================
780       typedef clitk::DifferenceImageFilter<  FixedImageType, FixedImageType> DifferenceFilterType;
781       if (m_ArgsInfo.after_given)
782       {
783         typename DifferenceFilterType::Pointer difference = DifferenceFilterType::New();
784         difference->SetValidInput( fixedImage );
785         difference->SetTestInput( result );
786
787         try
788         {
789           difference->Update();
790         }
791         catch( itk::ExceptionObject & err )
792         {
793           std::cerr << "ExceptionObject caught calculating the difference !" << std::endl;
794           std::cerr << err << std::endl;
795           return;
796         }
797
798         typename WriterType::Pointer differenceWriter=WriterType::New();
799         differenceWriter->SetInput(difference->GetOutput());
800         differenceWriter->SetFileName(m_ArgsInfo.after_arg);
801         differenceWriter->Update();
802
803       }
804
805
806       //=======================================================
807       // Calculate the difference before the deformable transform
808       //=======================================================
809       if( m_ArgsInfo.before_given )
810       {
811
812         typename FixedImageType::Pointer moving=FixedImageType::New();
813         if (m_ArgsInfo.initMatrix_given)
814         {
815           typedef itk::ResampleImageFilter<MovingImageType, FixedImageType> ResamplerType;
816           typename ResamplerType::Pointer resampler=ResamplerType::New();
817           resampler->SetInput(movingImage);
818           resampler->SetOutputOrigin(fixedImage->GetOrigin());
819           resampler->SetSize(fixedImage->GetLargestPossibleRegion().GetSize());
820           resampler->SetOutputSpacing(fixedImage->GetSpacing());
821           resampler->SetDefaultPixelValue( 0. );
822           if (rigidTransform ) resampler->SetTransform(rigidTransform);
823           resampler->Update();
824           moving=resampler->GetOutput();
825         }
826         else
827           moving=movingImage;
828
829         typename DifferenceFilterType::Pointer difference = DifferenceFilterType::New();
830         difference->SetValidInput( fixedImage );
831         difference->SetTestInput( moving );
832
833         try
834         {
835           difference->Update();
836         }
837         catch( itk::ExceptionObject & err )
838         {
839           std::cerr << "ExceptionObject caught calculating the difference !" << std::endl;
840           std::cerr << err << std::endl;
841           return;
842         }
843
844         typename WriterType::Pointer differenceWriter=WriterType::New();
845         writer->SetFileName( m_ArgsInfo.before_arg  );
846         writer->SetInput( difference->GetOutput()  );
847         writer->Update( );
848       }
849
850       return;
851
852     }
853 }//end clitk
854
855 #endif // #define clitkBLUTDIRGenericFilter_txx