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