]> Creatis software - clitk.git/blob - tools/clitkAffineTransformGenericFilter.txx
Merge branch 'master' of git.creatis.insa-lyon.fr:clitk
[clitk.git] / tools / clitkAffineTransformGenericFilter.txx
1 /*=========================================================================
2   Program:   vv                     http://www.creatis.insa-lyon.fr/rio/vv
3
4   Authors belong to:
5   - University of LYON              http://www.universite-lyon.fr/
6   - Léon Bérard cancer center       http://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 clitkAffineTransformGenericFilter_txx
19 #define clitkAffineTransformGenericFilter_txx
20
21 /* =================================================
22  * @file   clitkAffineTransformGenericFilter.txx
23  * @author
24  * @date
25  *
26  * @brief
27  *
28  ===================================================*/
29
30
31 namespace clitk
32 {
33
34 //-----------------------------------------------------------
35 // Constructor
36 //-----------------------------------------------------------
37 template<class args_info_type>
38 AffineTransformGenericFilter<args_info_type>::AffineTransformGenericFilter()
39 {
40   m_Verbose=false;
41   m_InputFileName="";
42 }
43
44
45 //-----------------------------------------------------------
46 // Update
47 //-----------------------------------------------------------
48 template<class args_info_type>
49 void AffineTransformGenericFilter<args_info_type>::Update()
50 {
51   // Read the Dimension and PixelType
52   int Dimension, Components;
53   std::string PixelType;
54   ReadImageDimensionAndPixelType(m_InputFileName, Dimension, PixelType, Components);
55
56
57   // Call UpdateWithDim
58   if(Dimension==2) UpdateWithDim<2>(PixelType, Components);
59   else if(Dimension==3) UpdateWithDim<3>(PixelType, Components);
60   else if (Dimension==4)UpdateWithDim<4>(PixelType, Components);
61   else {
62     std::cout<<"Error, Only for 2, 3 or 4  Dimensions!!!"<<std::endl ;
63     return;
64   }
65 }
66
67 //-------------------------------------------------------------------
68 // Update with the number of dimensions
69 //-------------------------------------------------------------------
70 template<class args_info_type>
71 template<unsigned int Dimension>
72 void
73 AffineTransformGenericFilter<args_info_type>::UpdateWithDim(std::string PixelType, int Components)
74 {
75   if (m_Verbose) std::cout << "Image was detected to be "<<Dimension<<"D and "<<Components<<" component(s) of "<<  PixelType<<"..."<<std::endl;
76
77   if (Components==1) {
78     if(PixelType == "short") {
79       if (m_Verbose) std::cout << "Launching filter in "<< Dimension <<"D and signed short..." << std::endl;
80       UpdateWithDimAndPixelType<Dimension, signed short>();
81     }
82     //    else if(PixelType == "unsigned_short"){
83     //       if (m_Verbose) std::cout  << "Launching filter in "<< Dimension <<"D and unsigned_short..." << std::endl;
84     //       UpdateWithDimAndPixelType<Dimension, unsigned short>();
85     //     }
86
87     else if (PixelType == "unsigned_char") {
88       if (m_Verbose) std::cout  << "Launching filter in "<< Dimension <<"D and unsigned_char..." << std::endl;
89       UpdateWithDimAndPixelType<Dimension, unsigned char>();
90     }
91
92     //     else if (PixelType == "char"){
93     //       if (m_Verbose) std::cout  << "Launching filter in "<< Dimension <<"D and signed_char..." << std::endl;
94     //       UpdateWithDimAndPixelType<Dimension, signed char>();
95     //     }
96     else {
97       if (m_Verbose) std::cout  << "Launching filter in "<< Dimension <<"D and float..." << std::endl;
98       UpdateWithDimAndPixelType<Dimension, float>();
99     }
100   }
101
102   else if (Components==3) {
103     if (m_Verbose) std::cout  << "Launching transform in "<< Dimension <<"D and 3D float (DVF)" << std::endl;
104     UpdateWithDimAndVectorType<Dimension, itk::Vector<float, Dimension> >();
105   }
106
107   else std::cerr<<"Number of components is "<<Components<<", not supported!"<<std::endl;
108
109 }
110
111
112 //-------------------------------------------------------------------
113 // Update with the number of dimensions and the pixeltype
114 //-------------------------------------------------------------------
115 template<class args_info_type>
116 template <unsigned int Dimension, class  PixelType>
117 void
118 AffineTransformGenericFilter<args_info_type>::UpdateWithDimAndPixelType()
119 {
120
121   // ImageTypes
122   typedef itk::Image<PixelType, Dimension> InputImageType;
123   typedef itk::Image<PixelType, Dimension> OutputImageType;
124
125   // Read the input
126   typedef itk::ImageFileReader<InputImageType> InputReaderType;
127   typename InputReaderType::Pointer reader = InputReaderType::New();
128   reader->SetFileName( m_InputFileName);
129   reader->Update();
130   typename InputImageType::Pointer input= reader->GetOutput();
131
132   //Filter
133   typedef  itk::ResampleImageFilter< InputImageType,OutputImageType >  ResampleFilterType;
134   typename ResampleFilterType::Pointer resampler = ResampleFilterType::New();
135
136   // Matrix
137   typename itk::Matrix<double, Dimension+1, Dimension+1> matrix;
138   if (m_ArgsInfo.rotate_given || m_ArgsInfo.translate_given)
139   {
140     if (m_ArgsInfo.matrix_given)
141     {
142       std::cerr << "You must use either rotate/translate or matrix options" << std::cout;
143       return;
144     }
145     itk::Array<double> transformParameters(2 * Dimension);
146     transformParameters.Fill(0.0);
147     if (m_ArgsInfo.rotate_given)
148     {
149       if (Dimension == 2)
150         transformParameters[0] = m_ArgsInfo.rotate_arg[0];
151       else
152         for (unsigned int i = 0; i < 3; i++)
153           transformParameters[i] = m_ArgsInfo.rotate_arg[i];
154     }
155     if (m_ArgsInfo.translate_given)
156     {
157       int pos = 3;
158       if (Dimension == 2)
159         pos = 1;
160       for (unsigned int i = 0; i < Dimension && i < 3; i++)
161         transformParameters[pos++] = m_ArgsInfo.translate_arg[i];
162     }
163     if (Dimension == 4)
164     {
165       matrix.SetIdentity();
166       itk::Matrix<double, 4, 4> tmp = GetForwardAffineMatrix3D(transformParameters);
167       for (unsigned int i = 0; i < 3; ++i)
168         for (unsigned int j = 0; j < 3; ++j)
169           matrix[i][j] = tmp[i][j];
170       for (unsigned int i = 0; i < 3; ++i)
171         matrix[i][4] = tmp[i][3];
172     }
173     else
174       matrix = GetForwardAffineMatrix<Dimension>(transformParameters);
175   }
176   else
177   {
178     if (m_ArgsInfo.matrix_given)
179     {
180       matrix= clitk::ReadMatrix<Dimension>(m_ArgsInfo.matrix_arg);
181       if (m_Verbose) std::cout << "Reading the matrix..." << std::endl;
182     }
183     else
184       matrix.SetIdentity();
185   }
186   if (m_Verbose)
187     std::cout << "Using the following matrix:" << std::endl
188               << matrix << std::endl;
189   typename itk::Matrix<double, Dimension, Dimension> rotationMatrix = clitk::GetRotationalPartMatrix(matrix);
190   typename itk::Vector<double, Dimension> translationPart = clitk::GetTranslationPartMatrix(matrix);
191
192   // Transform
193   typedef itk::AffineTransform<double, Dimension> AffineTransformType;
194   typename AffineTransformType::Pointer affineTransform=AffineTransformType::New();
195   affineTransform->SetMatrix(rotationMatrix);
196   affineTransform->SetTranslation(translationPart);
197
198   // Interp
199   typedef clitk::GenericInterpolator<args_info_type, InputImageType, double> GenericInterpolatorType;
200   typename GenericInterpolatorType::Pointer genericInterpolator=GenericInterpolatorType::New();
201   genericInterpolator->SetArgsInfo(m_ArgsInfo);
202
203   // Properties
204   if (m_ArgsInfo.like_given) {
205     typename InputReaderType::Pointer likeReader=InputReaderType::New();
206     likeReader->SetFileName(m_ArgsInfo.like_arg);
207     likeReader->Update();
208     resampler->SetOutputParametersFromImage(likeReader->GetOutput());
209   } else if(m_ArgsInfo.transform_grid_flag) {
210     typename itk::Matrix<double, Dimension+1, Dimension+1> invMatrix( matrix.GetInverse() );
211     typename itk::Matrix<double, Dimension, Dimension> invRotMatrix( clitk::GetRotationalPartMatrix(invMatrix) );
212     typename itk::Vector<double,Dimension> invTrans =  clitk::GetTranslationPartMatrix(invMatrix);
213
214     // Display warning
215     if (m_ArgsInfo.spacing_given)
216       std::cout << "Warning --spacing ignored (because --transform_grid_flag)" << std::endl;
217     if (m_ArgsInfo.origin_given)
218       std::cout << "Warning --origin ignored (because --transform_grid_flag)" << std::endl;
219
220     // Spacing is influenced by affine transform matrix and input direction
221     typename InputImageType::SpacingType outputSpacing;
222     outputSpacing = invRotMatrix *
223                     input->GetDirection() *
224                     input->GetSpacing();
225
226     // Origin is influenced by translation but not by input direction
227     typename InputImageType::PointType outputOrigin;
228     outputOrigin = invRotMatrix *
229                    input->GetOrigin() +
230                    invTrans;
231
232     // Size is influenced by affine transform matrix and input direction
233     // Size is converted to double, transformed and converted back to size type.
234     vnl_vector<double> vnlOutputSize(Dimension);
235     for(unsigned int i=0; i< Dimension; i++) {
236       vnlOutputSize[i] = input->GetLargestPossibleRegion().GetSize()[i];
237     }
238     vnlOutputSize = invRotMatrix *
239                     input->GetDirection().GetVnlMatrix() *
240                     vnlOutputSize;
241     typename OutputImageType::SizeType outputSize;
242     for(unsigned int i=0; i< Dimension; i++) {
243       // If the size is negative, we have a flip and we must modify
244       // the origin and the spacing accordingly.
245       if(vnlOutputSize[i]<0.) {
246         vnlOutputSize[i] *= -1.;
247         outputOrigin[i] = outputOrigin[i] + outputSpacing[i] * (vnlOutputSize[i]-1);
248         outputSpacing[i] *= -1.;
249       }
250       outputSize[i] = lrint(vnlOutputSize[i]);
251     }
252     resampler->SetSize( outputSize );
253     resampler->SetOutputSpacing( outputSpacing );
254     resampler->SetOutputOrigin( outputOrigin );
255   } else {
256     //Size
257     typename OutputImageType::SizeType outputSize;
258     if (m_ArgsInfo.size_given) {
259       for(unsigned int i=0; i< Dimension; i++)
260         outputSize[i]=m_ArgsInfo.size_arg[i];
261     } else outputSize=input->GetLargestPossibleRegion().GetSize();
262
263     //Spacing
264     typename OutputImageType::SpacingType outputSpacing;
265     if (m_ArgsInfo.spacing_given) {
266       for(unsigned int i=0; i< Dimension; i++)
267         outputSpacing[i]=m_ArgsInfo.spacing_arg[i];
268     } else outputSpacing=input->GetSpacing();
269
270     //Origin
271     typename OutputImageType::PointType outputOrigin;
272     if (m_ArgsInfo.origin_given) {
273       for(unsigned int i=0; i< Dimension; i++)
274         outputOrigin[i]=m_ArgsInfo.origin_arg[i];
275     } else outputOrigin=input->GetOrigin();
276
277     // Set
278     resampler->SetSize( outputSize );
279     resampler->SetOutputSpacing( outputSpacing );
280     resampler->SetOutputOrigin(  outputOrigin );
281
282   }
283
284   if (m_ArgsInfo.verbose_flag) {
285     std::cout << "Setting the output size to " << resampler->GetSize() << "..." << std::endl;
286     std::cout << "Setting the output spacing to " << resampler->GetOutputSpacing() << "..." << std::endl;
287     std::cout << "Setting the output origin to " << resampler->GetOutputOrigin() << "..." << std::endl;
288   }
289
290   resampler->SetInput( input );
291   resampler->SetTransform( affineTransform );
292   resampler->SetInterpolator( genericInterpolator->GetInterpolatorPointer());
293   resampler->SetDefaultPixelValue( static_cast<PixelType>(m_ArgsInfo.pad_arg) );
294
295   try {
296     resampler->Update();
297   } catch(itk::ExceptionObject) {
298     std::cerr<<"Error resampling the image"<<std::endl;
299   }
300
301   typename OutputImageType::Pointer output = resampler->GetOutput();
302
303   // Output
304   typedef itk::ImageFileWriter<OutputImageType> WriterType;
305   typename WriterType::Pointer writer = WriterType::New();
306   writer->SetFileName(m_ArgsInfo.output_arg);
307   writer->SetInput(output);
308   writer->Update();
309
310 }
311
312 //-------------------------------------------------------------------
313 // Update with the number of dimensions and the pixeltype (components)
314 //-------------------------------------------------------------------
315 template<class args_info_type>
316 template<unsigned int Dimension, class PixelType>
317 void AffineTransformGenericFilter<args_info_type>::UpdateWithDimAndVectorType()
318 {
319   // ImageTypes
320   typedef itk::Image<PixelType, Dimension> InputImageType;
321   typedef itk::Image<PixelType, Dimension> OutputImageType;
322
323   // Read the input
324   typedef itk::ImageFileReader<InputImageType> InputReaderType;
325   typename InputReaderType::Pointer reader = InputReaderType::New();
326   reader->SetFileName( m_InputFileName);
327   reader->Update();
328   typename InputImageType::Pointer input= reader->GetOutput();
329
330   //Filter
331   typedef  itk::VectorResampleImageFilter< InputImageType,OutputImageType, double >  ResampleFilterType;
332   typename ResampleFilterType::Pointer resampler = ResampleFilterType::New();
333
334   // Matrix
335   typename itk::Matrix<double, Dimension+1, Dimension+1> matrix;
336   if (m_ArgsInfo.rotate_given || m_ArgsInfo.translate_given)
337   {
338     if (m_ArgsInfo.matrix_given)
339     {
340       std::cerr << "You must use either rotate/translate or matrix options" << std::cout;
341       return;
342     }
343     itk::Array<double> transformParameters(2 * Dimension);
344     transformParameters.Fill(0.0);
345     if (m_ArgsInfo.rotate_given)
346     {
347       if (Dimension == 2)
348         transformParameters[0] = m_ArgsInfo.rotate_arg[0];
349       else
350         for (unsigned int i = 0; i < 3; i++)
351           transformParameters[i] = m_ArgsInfo.rotate_arg[i];
352     }
353     if (m_ArgsInfo.translate_given)
354     {
355       int pos = 3;
356       if (Dimension == 2)
357         pos = 1;
358       for (unsigned int i = 0; i < Dimension && i < 3; i++)
359         transformParameters[pos++] = m_ArgsInfo.translate_arg[i];
360     }
361     if (Dimension == 4)
362     {
363       matrix.SetIdentity();
364       itk::Matrix<double, 4, 4> tmp = GetForwardAffineMatrix3D(transformParameters);
365       for (unsigned int i = 0; i < 3; ++i)
366         for (unsigned int j = 0; j < 3; ++j)
367           matrix[i][j] = tmp[i][j];
368       for (unsigned int i = 0; i < 3; ++i)
369         matrix[i][4] = tmp[i][3];
370     }
371     else
372       matrix = GetForwardAffineMatrix<Dimension>(transformParameters);
373   }
374   else
375   {
376     if (m_ArgsInfo.matrix_given)
377     {
378       matrix= clitk::ReadMatrix<Dimension>(m_ArgsInfo.matrix_arg);
379       if (m_Verbose) std::cout << "Reading the matrix..." << std::endl;
380     }
381     else
382       matrix.SetIdentity();
383   }
384   if (m_Verbose)
385     std::cout << "Using the following matrix:" << std::endl
386               << matrix << std::endl;
387   typename itk::Matrix<double, Dimension, Dimension> rotationMatrix = clitk::GetRotationalPartMatrix(matrix);
388   typename itk::Vector<double, Dimension> translationPart = clitk::GetTranslationPartMatrix(matrix);
389
390   // Transform
391   typedef itk::AffineTransform<double, Dimension> AffineTransformType;
392   typename AffineTransformType::Pointer affineTransform=AffineTransformType::New();
393   affineTransform->SetMatrix(rotationMatrix);
394   affineTransform->SetTranslation(translationPart);
395
396   // Interp
397   typedef clitk::GenericVectorInterpolator<args_info_type, InputImageType, double> GenericInterpolatorType;
398   typename GenericInterpolatorType::Pointer genericInterpolator=GenericInterpolatorType::New();
399   genericInterpolator->SetArgsInfo(m_ArgsInfo);
400
401   // Properties
402   if (m_ArgsInfo.like_given) {
403     typename InputReaderType::Pointer likeReader=InputReaderType::New();
404     likeReader->SetFileName(m_ArgsInfo.like_arg);
405     likeReader->Update();
406     resampler->SetSize( likeReader->GetOutput()->GetLargestPossibleRegion().GetSize() );
407     resampler->SetOutputSpacing( likeReader->GetOutput()->GetSpacing() );
408     resampler->SetOutputOrigin(  likeReader->GetOutput()->GetOrigin() );
409   } else {
410     //Size
411     typename OutputImageType::SizeType outputSize;
412     if (m_ArgsInfo.size_given) {
413       for(unsigned int i=0; i< Dimension; i++)
414         outputSize[i]=m_ArgsInfo.size_arg[i];
415     } else outputSize=input->GetLargestPossibleRegion().GetSize();
416     std::cout<<"Setting the size to "<<outputSize<<"..."<<std::endl;
417
418     //Spacing
419     typename OutputImageType::SpacingType outputSpacing;
420     if (m_ArgsInfo.spacing_given) {
421       for(unsigned int i=0; i< Dimension; i++)
422         outputSpacing[i]=m_ArgsInfo.spacing_arg[i];
423     } else outputSpacing=input->GetSpacing();
424     std::cout<<"Setting the spacing to "<<outputSpacing<<"..."<<std::endl;
425
426     //Origin
427     typename OutputImageType::PointType outputOrigin;
428     if (m_ArgsInfo.origin_given) {
429       for(unsigned int i=0; i< Dimension; i++)
430         outputOrigin[i]=m_ArgsInfo.origin_arg[i];
431     } else outputOrigin=input->GetOrigin();
432     std::cout<<"Setting the origin to "<<outputOrigin<<"..."<<std::endl;
433
434     // Set
435     resampler->SetSize( outputSize );
436     resampler->SetOutputSpacing( outputSpacing );
437     resampler->SetOutputOrigin(  outputOrigin );
438
439   }
440
441   resampler->SetInput( input );
442   resampler->SetTransform( affineTransform );
443   resampler->SetInterpolator( genericInterpolator->GetInterpolatorPointer());
444   resampler->SetDefaultPixelValue( static_cast<PixelType>(m_ArgsInfo.pad_arg) );
445
446   try {
447     resampler->Update();
448   } catch(itk::ExceptionObject) {
449     std::cerr<<"Error resampling the image"<<std::endl;
450   }
451
452   typename OutputImageType::Pointer output = resampler->GetOutput();
453
454   // Output
455   typedef itk::ImageFileWriter<OutputImageType> WriterType;
456   typename WriterType::Pointer writer = WriterType::New();
457   writer->SetFileName(m_ArgsInfo.output_arg);
458   writer->SetInput(output);
459   writer->Update();
460
461 }
462
463
464 } //end clitk
465
466 #endif //#define clitkAffineTransformGenericFilter_txx