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>
29 //-----------------------------------------------------------
31 //-----------------------------------------------------------
32 template<class args_info_type>
33 AffineTransformGenericFilter<args_info_type>::AffineTransformGenericFilter()
38 //-------------------------------------------------------------------
41 //-----------------------------------------------------------
43 //-----------------------------------------------------------
44 template<class args_info_type>
45 void AffineTransformGenericFilter<args_info_type>::Update()
47 // Read the Dimension and PixelType
48 int Dimension, Components;
49 std::string PixelType;
50 ReadImageDimensionAndPixelType(m_InputFileName, Dimension, PixelType, Components);
53 if(Dimension==2) UpdateWithDim<2>(PixelType, Components);
55 if(Dimension==3) UpdateWithDim<3>(PixelType, Components);
56 else if (Dimension==4)UpdateWithDim<4>(PixelType, Components);
58 std::cout<<"Error, Only for 2, 3 or 4 Dimensions!!!"<<std::endl ;
62 //-------------------------------------------------------------------
65 //-------------------------------------------------------------------
66 // Update with the number of dimensions
67 //-------------------------------------------------------------------
68 template<class args_info_type>
69 template<unsigned int Dimension>
71 AffineTransformGenericFilter<args_info_type>::UpdateWithDim(std::string PixelType, int Components)
73 if (m_Verbose) std::cout << "Image was detected to be "<<Dimension<<"D and "<<Components<<" component(s) of "<< PixelType<<"..."<<std::endl;
76 if(PixelType == "short") {
77 if (m_Verbose) std::cout << "Launching filter in "<< Dimension <<"D and signed short..." << std::endl;
78 UpdateWithDimAndPixelType<Dimension, signed short>();
80 // else if(PixelType == "unsigned_short"){
81 // if (m_Verbose) std::cout << "Launching filter in "<< Dimension <<"D and unsigned_short..." << std::endl;
82 // UpdateWithDimAndPixelType<Dimension, unsigned short>();
85 else if (PixelType == "unsigned_char") {
86 if (m_Verbose) std::cout << "Launching filter in "<< Dimension <<"D and unsigned_char..." << std::endl;
87 UpdateWithDimAndPixelType<Dimension, unsigned char>();
90 // else if (PixelType == "char"){
91 // if (m_Verbose) std::cout << "Launching filter in "<< Dimension <<"D and signed_char..." << std::endl;
92 // UpdateWithDimAndPixelType<Dimension, signed char>();
95 if (m_Verbose) std::cout << "Launching filter in "<< Dimension <<"D and float..." << std::endl;
96 UpdateWithDimAndPixelType<Dimension, float>();
100 else if (Components==3) {
101 if (m_Verbose) std::cout << "Launching transform in "<< Dimension <<"D and 3D float (DVF)" << std::endl;
102 UpdateWithDimAndVectorType<Dimension, itk::Vector<float, Dimension> >();
105 else std::cerr<<"Number of components is "<<Components<<", not supported!"<<std::endl;
108 //-------------------------------------------------------------------
111 //-------------------------------------------------------------------
112 // Update with the number of dimensions and the pixeltype
113 //-------------------------------------------------------------------
114 template<class args_info_type>
115 template <unsigned int Dimension, class PixelType>
117 AffineTransformGenericFilter<args_info_type>::UpdateWithDimAndPixelType()
121 typedef itk::Image<PixelType, Dimension> InputImageType;
122 typedef itk::Image<PixelType, Dimension> OutputImageType;
125 typedef itk::ImageFileReader<InputImageType> InputReaderType;
126 typename InputReaderType::Pointer reader = InputReaderType::New();
127 reader->SetFileName( m_InputFileName);
129 typename InputImageType::Pointer input= reader->GetOutput();
132 typedef itk::ResampleImageFilter< InputImageType,OutputImageType > ResampleFilterType;
133 typename ResampleFilterType::Pointer resampler = ResampleFilterType::New();
136 typename itk::Matrix<double, Dimension+1, Dimension+1> matrix;
137 if (m_ArgsInfo.rotate_given || m_ArgsInfo.translate_given)
139 if (m_ArgsInfo.matrix_given)
141 std::cerr << "You must use either rotate/translate or matrix options" << std::cout;
144 itk::Array<double> transformParameters(2 * Dimension);
145 transformParameters.Fill(0.0);
146 if (m_ArgsInfo.rotate_given)
149 transformParameters[0] = m_ArgsInfo.rotate_arg[0];
151 for (unsigned int i = 0; i < 3; i++)
152 transformParameters[i] = m_ArgsInfo.rotate_arg[i];
154 if (m_ArgsInfo.translate_given)
159 for (unsigned int i = 0; i < Dimension && i < 3; i++)
160 transformParameters[pos++] = m_ArgsInfo.translate_arg[i];
164 matrix.SetIdentity();
165 itk::Matrix<double, 4, 4> tmp = GetForwardAffineMatrix3D(transformParameters);
166 for (unsigned int i = 0; i < 3; ++i)
167 for (unsigned int j = 0; j < 3; ++j)
168 matrix[i][j] = tmp[i][j];
169 for (unsigned int i = 0; i < 3; ++i)
170 matrix[i][4] = tmp[i][3];
173 matrix = GetForwardAffineMatrix<Dimension>(transformParameters);
177 if (m_ArgsInfo.matrix_given)
179 matrix= clitk::ReadMatrix<Dimension>(m_ArgsInfo.matrix_arg);
180 if (m_Verbose) std::cout << "Reading the matrix..." << std::endl;
183 if (m_ArgsInfo.elastix_given) {
184 matrix = createMatrixFromElastixFile<Dimension,PixelType>(m_ArgsInfo.elastix_arg);
187 matrix.SetIdentity();
191 std::cout << "Using the following matrix:" << std::endl
192 << matrix << std::endl;
193 typename itk::Matrix<double, Dimension, Dimension> rotationMatrix = clitk::GetRotationalPartMatrix(matrix);
194 typename itk::Vector<double, Dimension> translationPart = clitk::GetTranslationPartMatrix(matrix);
197 typedef itk::AffineTransform<double, Dimension> AffineTransformType;
198 typename AffineTransformType::Pointer affineTransform=AffineTransformType::New();
199 affineTransform->SetMatrix(rotationMatrix);
200 affineTransform->SetTranslation(translationPart);
203 typedef clitk::GenericInterpolator<args_info_type, InputImageType, double> GenericInterpolatorType;
204 typename GenericInterpolatorType::Pointer genericInterpolator=GenericInterpolatorType::New();
205 genericInterpolator->SetArgsInfo(m_ArgsInfo);
208 if (m_ArgsInfo.like_given) {
209 typename InputReaderType::Pointer likeReader=InputReaderType::New();
210 likeReader->SetFileName(m_ArgsInfo.like_arg);
211 likeReader->Update();
212 resampler->SetOutputParametersFromImage(likeReader->GetOutput());
213 } else if(m_ArgsInfo.transform_grid_flag) {
214 typename itk::Matrix<double, Dimension+1, Dimension+1> invMatrix( matrix.GetInverse() );
215 typename itk::Matrix<double, Dimension, Dimension> invRotMatrix( clitk::GetRotationalPartMatrix(invMatrix) );
216 typename itk::Vector<double,Dimension> invTrans = clitk::GetTranslationPartMatrix(invMatrix);
219 if (m_ArgsInfo.spacing_given)
220 std::cout << "Warning --spacing ignored (because --transform_grid_flag)" << std::endl;
221 if (m_ArgsInfo.origin_given)
222 std::cout << "Warning --origin ignored (because --transform_grid_flag)" << std::endl;
224 // Spacing is influenced by affine transform matrix and input direction
225 typename InputImageType::SpacingType outputSpacing;
226 outputSpacing = invRotMatrix *
227 input->GetDirection() *
230 // Origin is influenced by translation but not by input direction
231 typename InputImageType::PointType outputOrigin;
232 outputOrigin = invRotMatrix *
236 // Size is influenced by affine transform matrix and input direction
237 // Size is converted to double, transformed and converted back to size type.
238 vnl_vector<double> vnlOutputSize(Dimension);
239 for(unsigned int i=0; i< Dimension; i++) {
240 vnlOutputSize[i] = input->GetLargestPossibleRegion().GetSize()[i];
242 vnlOutputSize = invRotMatrix *
243 input->GetDirection().GetVnlMatrix() *
245 typename OutputImageType::SizeType outputSize;
246 for(unsigned int i=0; i< Dimension; i++) {
247 // If the size is negative, we have a flip and we must modify
248 // the origin and the spacing accordingly.
249 if(vnlOutputSize[i]<0.) {
250 vnlOutputSize[i] *= -1.;
251 outputOrigin[i] = outputOrigin[i] + outputSpacing[i] * (vnlOutputSize[i]-1);
252 outputSpacing[i] *= -1.;
254 outputSize[i] = lrint(vnlOutputSize[i]);
256 resampler->SetSize( outputSize );
257 resampler->SetOutputSpacing( outputSpacing );
258 resampler->SetOutputOrigin( outputOrigin );
261 typename OutputImageType::SizeType outputSize;
262 if (m_ArgsInfo.size_given) {
263 for(unsigned int i=0; i< Dimension; i++)
264 outputSize[i]=m_ArgsInfo.size_arg[i];
265 } else outputSize=input->GetLargestPossibleRegion().GetSize();
268 typename OutputImageType::SpacingType outputSpacing;
269 if (m_ArgsInfo.spacing_given) {
270 for(unsigned int i=0; i< Dimension; i++)
271 outputSpacing[i]=m_ArgsInfo.spacing_arg[i];
272 } else outputSpacing=input->GetSpacing();
275 typename OutputImageType::PointType outputOrigin;
276 if (m_ArgsInfo.origin_given) {
277 for(unsigned int i=0; i< Dimension; i++)
278 outputOrigin[i]=m_ArgsInfo.origin_arg[i];
279 } else outputOrigin=input->GetOrigin();
282 resampler->SetSize( outputSize );
283 resampler->SetOutputSpacing( outputSpacing );
284 resampler->SetOutputOrigin( outputOrigin );
288 if (m_ArgsInfo.spacinglike_given) {
289 typename InputReaderType::Pointer likeReader=InputReaderType::New();
290 likeReader->SetFileName(m_ArgsInfo.spacinglike_arg);
291 likeReader->Update();
293 // set the support like the image
294 if (m_ArgsInfo.like_given) {
295 typename OutputImageType::SizeType outputSize;
296 outputSize[0] = ceil(resampler->GetSize()[0]*resampler->GetOutputSpacing()[0]
297 /likeReader->GetOutput()->GetSpacing()[0]);
298 outputSize[1] = ceil(resampler->GetSize()[1]*resampler->GetOutputSpacing()[1]
299 /likeReader->GetOutput()->GetSpacing()[1]);
300 outputSize[2] = ceil(resampler->GetSize()[2]*resampler->GetOutputSpacing()[2]
301 /likeReader->GetOutput()->GetSpacing()[2]);
302 if (m_ArgsInfo.verbose_flag) {
303 std::cout << "Compute the number of pixels such as the support is like " << m_ArgsInfo.like_arg << std::endl;
305 resampler->SetSize( outputSize );
308 resampler->SetOutputSpacing( likeReader->GetOutput()->GetSpacing() );
311 if (m_ArgsInfo.verbose_flag) {
312 std::cout << "Setting the output size to " << resampler->GetSize() << "..." << std::endl;
313 std::cout << "Setting the output spacing to " << resampler->GetOutputSpacing() << "..." << std::endl;
314 std::cout << "Setting the output origin to " << resampler->GetOutputOrigin() << "..." << std::endl;
317 resampler->SetInput( input );
318 resampler->SetTransform( affineTransform );
319 resampler->SetInterpolator( genericInterpolator->GetInterpolatorPointer());
320 resampler->SetDefaultPixelValue( static_cast<PixelType>(m_ArgsInfo.pad_arg) );
324 } catch(itk::ExceptionObject) {
325 std::cerr<<"Error resampling the image"<<std::endl;
328 typename OutputImageType::Pointer output = resampler->GetOutput();
331 typedef itk::ImageFileWriter<OutputImageType> WriterType;
332 typename WriterType::Pointer writer = WriterType::New();
333 writer->SetFileName(m_ArgsInfo.output_arg);
334 writer->SetInput(output);
338 //-------------------------------------------------------------------
341 //-------------------------------------------------------------------
342 // Update with the number of dimensions and the pixeltype (components)
343 //-------------------------------------------------------------------
344 template<class args_info_type>
345 template<unsigned int Dimension, class PixelType>
346 void AffineTransformGenericFilter<args_info_type>::UpdateWithDimAndVectorType()
349 typedef itk::Image<PixelType, Dimension> InputImageType;
350 typedef itk::Image<PixelType, Dimension> OutputImageType;
353 typedef itk::ImageFileReader<InputImageType> InputReaderType;
354 typename InputReaderType::Pointer reader = InputReaderType::New();
355 reader->SetFileName( m_InputFileName);
357 typename InputImageType::Pointer input= reader->GetOutput();
360 typedef itk::VectorResampleImageFilter< InputImageType,OutputImageType, double > ResampleFilterType;
361 typename ResampleFilterType::Pointer resampler = ResampleFilterType::New();
364 typename itk::Matrix<double, Dimension+1, Dimension+1> matrix;
365 if (m_ArgsInfo.rotate_given || m_ArgsInfo.translate_given)
367 if (m_ArgsInfo.matrix_given)
369 std::cerr << "You must use either rotate/translate or matrix options" << std::cout;
372 itk::Array<double> transformParameters(2 * Dimension);
373 transformParameters.Fill(0.0);
374 if (m_ArgsInfo.rotate_given)
377 transformParameters[0] = m_ArgsInfo.rotate_arg[0];
379 for (unsigned int i = 0; i < 3; i++)
380 transformParameters[i] = m_ArgsInfo.rotate_arg[i];
382 if (m_ArgsInfo.translate_given)
387 for (unsigned int i = 0; i < Dimension && i < 3; i++)
388 transformParameters[pos++] = m_ArgsInfo.translate_arg[i];
392 matrix.SetIdentity();
393 itk::Matrix<double, 4, 4> tmp = GetForwardAffineMatrix3D(transformParameters);
394 for (unsigned int i = 0; i < 3; ++i)
395 for (unsigned int j = 0; j < 3; ++j)
396 matrix[i][j] = tmp[i][j];
397 for (unsigned int i = 0; i < 3; ++i)
398 matrix[i][4] = tmp[i][3];
401 matrix = GetForwardAffineMatrix<Dimension>(transformParameters);
405 if (m_ArgsInfo.matrix_given)
407 matrix= clitk::ReadMatrix<Dimension>(m_ArgsInfo.matrix_arg);
408 if (m_Verbose) std::cout << "Reading the matrix..." << std::endl;
411 matrix.SetIdentity();
414 std::cout << "Using the following matrix:" << std::endl
415 << matrix << std::endl;
416 typename itk::Matrix<double, Dimension, Dimension> rotationMatrix = clitk::GetRotationalPartMatrix(matrix);
417 typename itk::Vector<double, Dimension> translationPart = clitk::GetTranslationPartMatrix(matrix);
420 typedef itk::AffineTransform<double, Dimension> AffineTransformType;
421 typename AffineTransformType::Pointer affineTransform=AffineTransformType::New();
422 affineTransform->SetMatrix(rotationMatrix);
423 affineTransform->SetTranslation(translationPart);
426 typedef clitk::GenericVectorInterpolator<args_info_type, InputImageType, double> GenericInterpolatorType;
427 typename GenericInterpolatorType::Pointer genericInterpolator=GenericInterpolatorType::New();
428 genericInterpolator->SetArgsInfo(m_ArgsInfo);
431 if (m_ArgsInfo.like_given) {
432 typename InputReaderType::Pointer likeReader=InputReaderType::New();
433 likeReader->SetFileName(m_ArgsInfo.like_arg);
434 likeReader->Update();
435 resampler->SetSize( likeReader->GetOutput()->GetLargestPossibleRegion().GetSize() );
436 resampler->SetOutputSpacing( likeReader->GetOutput()->GetSpacing() );
437 resampler->SetOutputOrigin( likeReader->GetOutput()->GetOrigin() );
440 typename OutputImageType::SizeType outputSize;
441 if (m_ArgsInfo.size_given) {
442 for(unsigned int i=0; i< Dimension; i++)
443 outputSize[i]=m_ArgsInfo.size_arg[i];
444 } else outputSize=input->GetLargestPossibleRegion().GetSize();
445 std::cout<<"Setting the size to "<<outputSize<<"..."<<std::endl;
448 typename OutputImageType::SpacingType outputSpacing;
449 if (m_ArgsInfo.spacing_given) {
450 for(unsigned int i=0; i< Dimension; i++)
451 outputSpacing[i]=m_ArgsInfo.spacing_arg[i];
452 } else outputSpacing=input->GetSpacing();
453 std::cout<<"Setting the spacing to "<<outputSpacing<<"..."<<std::endl;
456 typename OutputImageType::PointType outputOrigin;
457 if (m_ArgsInfo.origin_given) {
458 for(unsigned int i=0; i< Dimension; i++)
459 outputOrigin[i]=m_ArgsInfo.origin_arg[i];
460 } else outputOrigin=input->GetOrigin();
461 std::cout<<"Setting the origin to "<<outputOrigin<<"..."<<std::endl;
464 resampler->SetSize( outputSize );
465 resampler->SetOutputSpacing( outputSpacing );
466 resampler->SetOutputOrigin( outputOrigin );
470 resampler->SetInput( input );
471 resampler->SetTransform( affineTransform );
472 resampler->SetInterpolator( genericInterpolator->GetInterpolatorPointer());
473 resampler->SetDefaultPixelValue( static_cast<PixelType>(m_ArgsInfo.pad_arg) );
477 } catch(itk::ExceptionObject) {
478 std::cerr<<"Error resampling the image"<<std::endl;
481 typename OutputImageType::Pointer output = resampler->GetOutput();
484 typedef itk::ImageFileWriter<OutputImageType> WriterType;
485 typename WriterType::Pointer writer = WriterType::New();
486 writer->SetFileName(m_ArgsInfo.output_arg);
487 writer->SetInput(output);
491 //-------------------------------------------------------------------
494 //-------------------------------------------------------------------
495 template<class args_info_type>
496 template<unsigned int Dimension, class PixelType>
497 typename itk::Matrix<double, Dimension+1, Dimension+1>
498 AffineTransformGenericFilter<args_info_type>::createMatrixFromElastixFile(std::string filename)
500 if (Dimension != 3) {
501 FATAL("Only 3D yet" << std::endl);
503 typename itk::Matrix<double, Dimension+1, Dimension+1> matrix;
507 clitk::openFileForReading(is, filename);
511 bool b = GetElastixValueFromTag(is, "Transform ", s);
513 FATAL("Error must read 'Transform' in " << filename << std::endl);
515 if (s != "EulerTransform") {
516 FATAL("Sorry only 'EulerTransform'" << std::endl);
520 // (InitialTransformParametersFileName "NoInitialTransform")
522 // Get CenterOfRotationPoint
523 GetElastixValueFromTag(is, "CenterOfRotationPoint ", s); // space is needed
525 FATAL("Error must read 'CenterOfRotationPoint' in " << filename << std::endl);
527 std::vector<std::string> cor;
528 GetValuesFromValue(s, cor);
530 // Get Transformparameters
531 GetElastixValueFromTag(is, "TransformParameters ", s); // space is needed
533 FATAL("Error must read 'TransformParameters' in " << filename << std::endl);
535 std::vector<std::string> results;
536 GetValuesFromValue(s, results);
538 // construct a stream from the string
539 itk::CenteredEuler3DTransform<double>::Pointer mat = itk::CenteredEuler3DTransform<double>::New();
540 itk::CenteredEuler3DTransform<double>::ParametersType p;
542 for(uint i=0; i<3; i++)
543 p[i] = atof(results[i].c_str()); // Rotation
544 for(uint i=0; i<3; i++)
545 p[i+3] = atof(cor[i].c_str()); // Centre of rotation
546 for(uint i=0; i<3; i++)
547 p[i+6] = atof(results[i+3].c_str()); // Translation
548 mat->SetParameters(p);
551 std::cout << "Rotation (deg) : " << rad2deg(p[0]) << " " << rad2deg(p[1]) << " " << rad2deg(p[2]) << std::endl;
552 std::cout << "Translation (phy) : " << p[3] << " " << p[4] << " " << p[5] << std::endl;
553 std::cout << "Center of rot (phy) : " << p[6] << " " << p[7] << " " << p[8] << std::endl;
556 for(uint i=0; i<3; i++)
557 for(uint j=0; j<3; j++)
558 matrix[i][j] = mat->GetMatrix()[i][j];
559 // Offset is -Rc + t + c
560 matrix[0][3] = mat->GetOffset()[0];
561 matrix[1][3] = mat->GetOffset()[1];
562 matrix[2][3] = mat->GetOffset()[2];
568 //-------------------------------------------------------------------
569 template<class args_info_type>
571 AffineTransformGenericFilter<args_info_type>::GetElastixValueFromTag(std::ifstream & is,
576 is.seekg (0, is.beg);
577 while(std::getline(is, line)) {
578 unsigned pos = line.find(tag);
579 if (pos<line.size()) {
580 value=line.substr(pos+tag.size(),line.size()-2);// remove last ')'
581 value.erase (std::remove (value.begin(), value.end(), '"'), value.end());
582 value.erase (std::remove (value.begin(), value.end(), ')'), value.end());
588 //-------------------------------------------------------------------
591 //-------------------------------------------------------------------
592 template<class args_info_type>
594 AffineTransformGenericFilter<args_info_type>::GetValuesFromValue(const std::string & s,
595 std::vector<std::string> & values)
597 std::stringstream strstr(s);
598 std::istream_iterator<std::string> it(strstr);
599 std::istream_iterator<std::string> end;
600 std::vector<std::string> results(it, end);
602 values.resize(results.size());
603 for(uint i=0; i<results.size(); i++) values[i] = results[i];
605 //-------------------------------------------------------------------
610 #endif //#define clitkAffineTransformGenericFilter_txx