1 /*=========================================================================
2 Program: vv http://www.creatis.insa-lyon.fr/rio/vv
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
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.
13 It is distributed under dual licence
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
24 #include <itkCenteredEuler3DTransform.h>
25 #include "clitkElastix.h"
30 //-----------------------------------------------------------
32 //-----------------------------------------------------------
33 template<class args_info_type>
34 AffineTransformGenericFilter<args_info_type>::AffineTransformGenericFilter()
39 //-------------------------------------------------------------------
42 //-----------------------------------------------------------
44 //-----------------------------------------------------------
45 template<class args_info_type>
46 void AffineTransformGenericFilter<args_info_type>::Update()
48 // Read the Dimension and PixelType
49 int Dimension, Components;
50 std::string PixelType;
51 ReadImageDimensionAndPixelType(m_InputFileName, Dimension, PixelType, Components);
54 if(Dimension==2) UpdateWithDim<2>(PixelType, Components);
56 if(Dimension==3) UpdateWithDim<3>(PixelType, Components);
57 else if (Dimension==4)UpdateWithDim<4>(PixelType, Components);
59 std::cout<<"Error, Only for 2, 3 or 4 Dimensions!!!"<<std::endl ;
63 //-------------------------------------------------------------------
66 //-------------------------------------------------------------------
67 // Update with the number of dimensions
68 //-------------------------------------------------------------------
69 template<class args_info_type>
70 template<unsigned int Dimension>
72 AffineTransformGenericFilter<args_info_type>::UpdateWithDim(std::string PixelType, int Components)
74 if (m_Verbose) std::cout << "Image was detected to be "<<Dimension<<"D and "<<Components<<" component(s) of "<< PixelType<<"..."<<std::endl;
77 if(PixelType == "short") {
78 if (m_Verbose) std::cout << "Launching filter in "<< Dimension <<"D and signed short..." << std::endl;
79 UpdateWithDimAndPixelType<Dimension, signed short>();
81 else if(PixelType == "unsigned_short"){
82 if (m_Verbose) std::cout << "Launching filter in "<< Dimension <<"D and unsigned_short..." << std::endl;
83 UpdateWithDimAndPixelType<Dimension, unsigned short>();
86 else if (PixelType == "unsigned_char") {
87 if (m_Verbose) std::cout << "Launching filter in "<< Dimension <<"D and unsigned_char..." << std::endl;
88 UpdateWithDimAndPixelType<Dimension, unsigned char>();
91 // else if (PixelType == "char"){
92 // if (m_Verbose) std::cout << "Launching filter in "<< Dimension <<"D and signed_char..." << std::endl;
93 // UpdateWithDimAndPixelType<Dimension, signed char>();
96 if (m_Verbose) std::cout << "Launching filter in "<< Dimension <<"D and float..." << std::endl;
97 UpdateWithDimAndPixelType<Dimension, float>();
101 else if (Components==3) {
102 if (m_Verbose) std::cout << "Launching transform in "<< Dimension <<"D and 3D float (DVF)" << std::endl;
103 UpdateWithDimAndVectorType<Dimension, itk::Vector<float, Dimension> >();
106 else std::cerr<<"Number of components is "<<Components<<", not supported!"<<std::endl;
109 //-------------------------------------------------------------------
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>
118 AffineTransformGenericFilter<args_info_type>::UpdateWithDimAndPixelType()
122 typedef itk::Image<PixelType, Dimension> InputImageType;
123 typedef itk::Image<PixelType, Dimension> OutputImageType;
126 typedef itk::ImageFileReader<InputImageType> InputReaderType;
127 typename InputReaderType::Pointer reader = InputReaderType::New();
128 reader->SetFileName( m_InputFileName);
130 typename InputImageType::Pointer input= reader->GetOutput();
133 typedef itk::ResampleImageFilter< InputImageType,OutputImageType > ResampleFilterType;
134 typename ResampleFilterType::Pointer resampler = ResampleFilterType::New();
137 typename itk::Matrix<double, Dimension+1, Dimension+1> matrix;
138 if (m_ArgsInfo.rotate_given || m_ArgsInfo.translate_given)
140 if (m_ArgsInfo.matrix_given)
142 std::cerr << "You must use either rotate/translate or matrix options" << std::cout;
145 itk::Array<double> transformParameters(2 * Dimension);
146 transformParameters.Fill(0.0);
147 if (m_ArgsInfo.rotate_given)
150 transformParameters[0] = m_ArgsInfo.rotate_arg[0];
152 for (unsigned int i = 0; i < 3; i++)
153 transformParameters[i] = m_ArgsInfo.rotate_arg[i];
155 if (m_ArgsInfo.translate_given)
160 for (unsigned int i = 0; i < Dimension && i < 3; i++)
161 transformParameters[pos++] = m_ArgsInfo.translate_arg[i];
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];
174 matrix = GetForwardAffineMatrix<Dimension>(transformParameters);
178 if (m_ArgsInfo.matrix_given)
180 matrix= clitk::ReadMatrix<Dimension>(m_ArgsInfo.matrix_arg);
181 if (m_Verbose) std::cout << "Reading the matrix..." << std::endl;
184 if (m_ArgsInfo.elastix_given) {
185 std::string filename(m_ArgsInfo.elastix_arg);
186 matrix = createMatrixFromElastixFile<Dimension>(filename, m_Verbose);
189 matrix.SetIdentity();
193 std::cout << "Using the following matrix:" << std::endl
194 << matrix << std::endl;
195 typename itk::Matrix<double, Dimension, Dimension> rotationMatrix = clitk::GetRotationalPartMatrix(matrix);
196 typename itk::Vector<double, Dimension> translationPart = clitk::GetTranslationPartMatrix(matrix);
199 typedef itk::AffineTransform<double, Dimension> AffineTransformType;
200 typename AffineTransformType::Pointer affineTransform=AffineTransformType::New();
201 affineTransform->SetMatrix(rotationMatrix);
202 affineTransform->SetTranslation(translationPart);
205 typedef clitk::GenericInterpolator<args_info_type, InputImageType, double> GenericInterpolatorType;
206 typename GenericInterpolatorType::Pointer genericInterpolator=GenericInterpolatorType::New();
207 genericInterpolator->SetArgsInfo(m_ArgsInfo);
210 if (m_ArgsInfo.like_given) {
211 typename InputReaderType::Pointer likeReader=InputReaderType::New();
212 likeReader->SetFileName(m_ArgsInfo.like_arg);
213 likeReader->Update();
214 resampler->SetOutputParametersFromImage(likeReader->GetOutput());
215 resampler->SetOutputDirection(likeReader->GetOutput()->GetDirection());
216 } else if(m_ArgsInfo.transform_grid_flag) {
217 typename itk::Matrix<double, Dimension+1, Dimension+1> invMatrix( matrix.GetInverse() );
218 typename itk::Matrix<double, Dimension, Dimension> invRotMatrix( clitk::GetRotationalPartMatrix(invMatrix) );
219 typename itk::Vector<double,Dimension> invTrans = clitk::GetTranslationPartMatrix(invMatrix);
222 if (m_ArgsInfo.spacing_given)
223 std::cout << "Warning --spacing ignored (because --transform_grid_flag)" << std::endl;
224 if (m_ArgsInfo.origin_given)
225 std::cout << "Warning --origin ignored (because --transform_grid_flag)" << std::endl;
227 // Spacing is influenced by affine transform matrix and input direction
228 typename InputImageType::SpacingType outputSpacing;
229 outputSpacing = invRotMatrix *
230 input->GetDirection() *
233 // Origin is influenced by translation but not by input direction
234 typename InputImageType::PointType outputOrigin;
235 outputOrigin = invRotMatrix *
239 // Size is influenced by affine transform matrix and input direction
240 // Size is converted to double, transformed and converted back to size type.
241 vnl_vector<double> vnlOutputSize(Dimension);
242 for(unsigned int i=0; i< Dimension; i++) {
243 vnlOutputSize[i] = input->GetLargestPossibleRegion().GetSize()[i];
245 vnlOutputSize = invRotMatrix *
246 input->GetDirection().GetVnlMatrix() *
248 typename OutputImageType::SizeType outputSize;
249 for(unsigned int i=0; i< Dimension; i++) {
250 // If the size is negative, we have a flip and we must modify
251 // the origin and the spacing accordingly.
252 if(vnlOutputSize[i]<0.) {
253 vnlOutputSize[i] *= -1.;
254 outputOrigin[i] = outputOrigin[i] + outputSpacing[i] * (vnlOutputSize[i]-1);
255 outputSpacing[i] *= -1.;
257 outputSize[i] = lrint(vnlOutputSize[i]);
259 resampler->SetSize( outputSize );
260 resampler->SetOutputSpacing( outputSpacing );
261 resampler->SetOutputOrigin( outputOrigin );
264 typename OutputImageType::SizeType outputSize;
265 if (m_ArgsInfo.size_given) {
266 for(unsigned int i=0; i< Dimension; i++)
267 outputSize[i]=m_ArgsInfo.size_arg[i];
268 } else outputSize=input->GetLargestPossibleRegion().GetSize();
271 typename OutputImageType::SpacingType outputSpacing;
272 if (m_ArgsInfo.spacing_given) {
273 for(unsigned int i=0; i< Dimension; i++)
274 outputSpacing[i]=m_ArgsInfo.spacing_arg[i];
275 } else outputSpacing=input->GetSpacing();
278 typename OutputImageType::PointType outputOrigin;
279 if (m_ArgsInfo.origin_given) {
280 for(unsigned int i=0; i< Dimension; i++)
281 outputOrigin[i]=m_ArgsInfo.origin_arg[i];
282 } else outputOrigin=input->GetOrigin();
285 resampler->SetSize( outputSize );
286 resampler->SetOutputSpacing( outputSpacing );
287 resampler->SetOutputOrigin( outputOrigin );
291 if (m_ArgsInfo.spacinglike_given) {
292 typename InputReaderType::Pointer likeReader=InputReaderType::New();
293 likeReader->SetFileName(m_ArgsInfo.spacinglike_arg);
294 likeReader->Update();
296 // set the support like the image
297 if (m_ArgsInfo.like_given) {
298 typename OutputImageType::SizeType outputSize;
299 outputSize[0] = ceil(resampler->GetSize()[0]*resampler->GetOutputSpacing()[0]
300 /likeReader->GetOutput()->GetSpacing()[0]);
301 outputSize[1] = ceil(resampler->GetSize()[1]*resampler->GetOutputSpacing()[1]
302 /likeReader->GetOutput()->GetSpacing()[1]);
303 outputSize[2] = ceil(resampler->GetSize()[2]*resampler->GetOutputSpacing()[2]
304 /likeReader->GetOutput()->GetSpacing()[2]);
305 if (m_ArgsInfo.verbose_flag) {
306 std::cout << "Compute the number of pixels such as the support is like " << m_ArgsInfo.like_arg << std::endl;
308 resampler->SetSize( outputSize );
311 resampler->SetOutputSpacing( likeReader->GetOutput()->GetSpacing() );
314 if (m_ArgsInfo.verbose_flag) {
315 std::cout << "Setting the output size to " << resampler->GetSize() << "..." << std::endl;
316 std::cout << "Setting the output spacing to " << resampler->GetOutputSpacing() << "..." << std::endl;
317 std::cout << "Setting the output origin to " << resampler->GetOutputOrigin() << "..." << std::endl;
320 resampler->SetInput( input );
321 resampler->SetTransform( affineTransform );
322 resampler->SetInterpolator( genericInterpolator->GetInterpolatorPointer());
323 resampler->SetDefaultPixelValue( static_cast<PixelType>(m_ArgsInfo.pad_arg) );
327 } catch(itk::ExceptionObject) {
328 std::cerr<<"Error resampling the image"<<std::endl;
331 typename OutputImageType::Pointer output = resampler->GetOutput();
334 typedef itk::ImageFileWriter<OutputImageType> WriterType;
335 typename WriterType::Pointer writer = WriterType::New();
336 writer->SetFileName(m_ArgsInfo.output_arg);
337 writer->SetInput(output);
341 //-------------------------------------------------------------------
344 //-------------------------------------------------------------------
345 // Update with the number of dimensions and the pixeltype (components)
346 //-------------------------------------------------------------------
347 template<class args_info_type>
348 template<unsigned int Dimension, class PixelType>
349 void AffineTransformGenericFilter<args_info_type>::UpdateWithDimAndVectorType()
352 typedef itk::Image<PixelType, Dimension> InputImageType;
353 typedef itk::Image<PixelType, Dimension> OutputImageType;
356 typedef itk::ImageFileReader<InputImageType> InputReaderType;
357 typename InputReaderType::Pointer reader = InputReaderType::New();
358 reader->SetFileName( m_InputFileName);
360 typename InputImageType::Pointer input= reader->GetOutput();
363 typedef itk::VectorResampleImageFilter< InputImageType,OutputImageType, double > ResampleFilterType;
364 typename ResampleFilterType::Pointer resampler = ResampleFilterType::New();
367 typename itk::Matrix<double, Dimension+1, Dimension+1> matrix;
368 if (m_ArgsInfo.rotate_given || m_ArgsInfo.translate_given)
370 if (m_ArgsInfo.matrix_given)
372 std::cerr << "You must use either rotate/translate or matrix options" << std::cout;
375 itk::Array<double> transformParameters(2 * Dimension);
376 transformParameters.Fill(0.0);
377 if (m_ArgsInfo.rotate_given)
380 transformParameters[0] = m_ArgsInfo.rotate_arg[0];
382 for (unsigned int i = 0; i < 3; i++)
383 transformParameters[i] = m_ArgsInfo.rotate_arg[i];
385 if (m_ArgsInfo.translate_given)
390 for (unsigned int i = 0; i < Dimension && i < 3; i++)
391 transformParameters[pos++] = m_ArgsInfo.translate_arg[i];
395 matrix.SetIdentity();
396 itk::Matrix<double, 4, 4> tmp = GetForwardAffineMatrix3D(transformParameters);
397 for (unsigned int i = 0; i < 3; ++i)
398 for (unsigned int j = 0; j < 3; ++j)
399 matrix[i][j] = tmp[i][j];
400 for (unsigned int i = 0; i < 3; ++i)
401 matrix[i][4] = tmp[i][3];
404 matrix = GetForwardAffineMatrix<Dimension>(transformParameters);
408 if (m_ArgsInfo.matrix_given)
410 matrix= clitk::ReadMatrix<Dimension>(m_ArgsInfo.matrix_arg);
411 if (m_Verbose) std::cout << "Reading the matrix..." << std::endl;
414 matrix.SetIdentity();
417 std::cout << "Using the following matrix:" << std::endl
418 << matrix << std::endl;
419 typename itk::Matrix<double, Dimension, Dimension> rotationMatrix = clitk::GetRotationalPartMatrix(matrix);
420 typename itk::Vector<double, Dimension> translationPart = clitk::GetTranslationPartMatrix(matrix);
423 typedef itk::AffineTransform<double, Dimension> AffineTransformType;
424 typename AffineTransformType::Pointer affineTransform=AffineTransformType::New();
425 affineTransform->SetMatrix(rotationMatrix);
426 affineTransform->SetTranslation(translationPart);
429 typedef clitk::GenericVectorInterpolator<args_info_type, InputImageType, double> GenericInterpolatorType;
430 typename GenericInterpolatorType::Pointer genericInterpolator=GenericInterpolatorType::New();
431 genericInterpolator->SetArgsInfo(m_ArgsInfo);
434 if (m_ArgsInfo.like_given) {
435 typename InputReaderType::Pointer likeReader=InputReaderType::New();
436 likeReader->SetFileName(m_ArgsInfo.like_arg);
437 likeReader->Update();
438 resampler->SetSize( likeReader->GetOutput()->GetLargestPossibleRegion().GetSize() );
439 resampler->SetOutputSpacing( likeReader->GetOutput()->GetSpacing() );
440 resampler->SetOutputOrigin( likeReader->GetOutput()->GetOrigin() );
443 typename OutputImageType::SizeType outputSize;
444 if (m_ArgsInfo.size_given) {
445 for(unsigned int i=0; i< Dimension; i++)
446 outputSize[i]=m_ArgsInfo.size_arg[i];
447 } else outputSize=input->GetLargestPossibleRegion().GetSize();
448 std::cout<<"Setting the size to "<<outputSize<<"..."<<std::endl;
451 typename OutputImageType::SpacingType outputSpacing;
452 if (m_ArgsInfo.spacing_given) {
453 for(unsigned int i=0; i< Dimension; i++)
454 outputSpacing[i]=m_ArgsInfo.spacing_arg[i];
455 } else outputSpacing=input->GetSpacing();
456 std::cout<<"Setting the spacing to "<<outputSpacing<<"..."<<std::endl;
459 typename OutputImageType::PointType outputOrigin;
460 if (m_ArgsInfo.origin_given) {
461 for(unsigned int i=0; i< Dimension; i++)
462 outputOrigin[i]=m_ArgsInfo.origin_arg[i];
463 } else outputOrigin=input->GetOrigin();
464 std::cout<<"Setting the origin to "<<outputOrigin<<"..."<<std::endl;
467 resampler->SetSize( outputSize );
468 resampler->SetOutputSpacing( outputSpacing );
469 resampler->SetOutputOrigin( outputOrigin );
473 resampler->SetInput( input );
474 resampler->SetTransform( affineTransform );
475 resampler->SetInterpolator( genericInterpolator->GetInterpolatorPointer());
476 resampler->SetDefaultPixelValue( static_cast<PixelType>(m_ArgsInfo.pad_arg) );
480 } catch(itk::ExceptionObject) {
481 std::cerr<<"Error resampling the image"<<std::endl;
484 typename OutputImageType::Pointer output = resampler->GetOutput();
487 typedef itk::ImageFileWriter<OutputImageType> WriterType;
488 typename WriterType::Pointer writer = WriterType::New();
489 writer->SetFileName(m_ArgsInfo.output_arg);
490 writer->SetInput(output);
494 //-------------------------------------------------------------------
498 #endif //#define clitkAffineTransformGenericFilter_txx