]> Creatis software - clitk.git/commitdiff
refactoring...
authorRomulo Pinho <romulo.pinho@lyon.unicancer.fr>
Thu, 19 Jan 2012 16:51:00 +0000 (17:51 +0100)
committerRomulo Pinho <romulo.pinho@lyon.unicancer.fr>
Thu, 19 Jan 2012 16:51:00 +0000 (17:51 +0100)
- Put coeff to dvf function in separate file, as a global function
- Changed affected files/classes

common/clitkCoeffsToDVF.h [new file with mode: 0644]
tools/clitkComposeVFGenericFilter.h
tools/clitkComposeVFGenericFilter.txx
tools/clitkInvertVFGenericFilter.h
tools/clitkInvertVFGenericFilter.txx
tools/clitkWarpImageGenericFilter.h
tools/clitkWarpImageGenericFilter.txx

diff --git a/common/clitkCoeffsToDVF.h b/common/clitkCoeffsToDVF.h
new file mode 100644 (file)
index 0000000..0fe30bb
--- /dev/null
@@ -0,0 +1,86 @@
+/*=========================================================================
+  Program:   vv                     http://www.creatis.insa-lyon.fr/rio/vv
+
+  Authors belong to:
+  - University of LYON              http://www.universite-lyon.fr/
+  - Léon Bérard cancer center       http://www.centreleonberard.fr
+  - CREATIS CNRS laboratory         http://www.creatis.insa-lyon.fr
+
+  This software is distributed WITHOUT ANY WARRANTY; without even
+  the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
+  PURPOSE.  See the copyright notices for more information.
+
+  It is distributed under dual licence
+
+  - BSD        See included LICENSE.txt file
+  - CeCILL-B   http://www.cecill.info/licences/Licence_CeCILL-B_V1-en.html
+===========================================================================**/
+#ifndef clitkCoeffsToDVF_h
+#define clitkCoeffsToDVF_h
+
+#include "clitkBSplineDeformableTransform.h"
+#if ITK_VERSION_MAJOR >= 4
+#include "itkTransformToDisplacementFieldSource.h"
+#else
+#include "itkTransformToDeformationFieldSource.h"
+#endif
+
+//-------------------------------------------------------------------
+// Convert Coefficient image to DVF
+//-------------------------------------------------------------------
+template<class DisplacementFieldType>
+typename DisplacementFieldType::Pointer
+CoeffsToDVF(std::string fileName, std::string likeFileName, bool verbose = false)
+{
+  typedef clitk::BSplineDeformableTransform<double, DisplacementFieldType::ImageDimension, DisplacementFieldType::ImageDimension> TransformType;
+  typedef typename TransformType::CoefficientImageType CoefficientImageType;
+
+  typedef itk::ImageFileReader<CoefficientImageType> CoeffReaderType;
+  typename CoeffReaderType::Pointer reader = CoeffReaderType::New();
+  reader->SetFileName(fileName);
+  reader->Update();
+
+  typename TransformType::Pointer transform = TransformType::New();
+  transform->SetCoefficientImage(reader->GetOutput());
+  
+#if ITK_VERSION_MAJOR >= 4
+      typedef itk::TransformToDisplacementFieldSource<DisplacementFieldType, double> ConvertorType;
+#else
+      typedef itk::TransformToDeformationFieldSource<DisplacementFieldType, double> ConvertorType;
+#endif
+
+  typedef itk::ImageIOBase ImageIOType;
+  typename ImageIOType::Pointer imageIO = itk::ImageIOFactory::CreateImageIO(likeFileName.c_str(), itk::ImageIOFactory::ReadMode);
+  imageIO->SetFileName(likeFileName);
+  imageIO->ReadImageInformation();
+
+  typename ConvertorType::Pointer convertor= ConvertorType::New();
+  typename ConvertorType::SizeType output_size;
+  typename ConvertorType::SpacingType output_spacing;
+  typename ConvertorType::OriginType output_origin;
+  typename ConvertorType::DirectionType output_direction;
+  for (unsigned int i = 0; i < DisplacementFieldType::ImageDimension; i++) {
+    output_size[i] = imageIO->GetDimensions(i);
+    output_spacing[i] = imageIO->GetSpacing(i);
+    output_origin[i] = imageIO->GetOrigin(i);
+    for (unsigned int j = 0; j < DisplacementFieldType::ImageDimension; j++)
+      output_direction[i][j] = imageIO->GetDirection(i)[j];
+  }
+  
+  if (verbose) {
+    std::cout << "Interpolating coefficients with grid:" << std::endl;
+    std::cout << output_size << output_spacing << std::endl;
+  }
+  
+  convertor->SetNumberOfThreads(1);
+  convertor->SetTransform(transform);
+  convertor->SetOutputOrigin(output_origin);
+  convertor->SetOutputSpacing(output_spacing);
+  convertor->SetOutputSize(output_size);
+  convertor->SetOutputDirection(output_direction);
+  convertor->Update();
+
+  return convertor->GetOutput();
+}
+
+#endif
index b86f22b92974d89bcc8c1f8fc0d30a0de172cde5..4d9d9f5d71873380e403827c1f3c9e507fa80968 100644 (file)
@@ -67,7 +67,6 @@ namespace clitk
     //Templated members
     template<unsigned int Dimension> void UpdateWithDim(std::string PixelType);
     template<unsigned int Dimension, class PixelType> void UpdateWithDimAndPixelType(); 
-    template<class DisplacementFieldType> typename DisplacementFieldType::Pointer CoeffsToDVF(std::string fileName, std::string likeFileName);
     
     std::string m_InputName1;
     std::string m_InputName2;
index 07df9882fe573cafd268b4b78a83d9e914cc80f7..8b0ecdf7ad5408c99eb5adbd4a8d0227da83e029 100644 (file)
 #define __clitkComposeVFGenericFilter_txx
 #include "clitkComposeVFGenericFilter.h"
 
-#include "clitkBSplineDeformableTransform.h"
-#include "clitkBSplineDeformableTransformInitializer.h"
-#if ITK_VERSION_MAJOR >= 4
-#include "itkTransformToDisplacementFieldSource.h"
-#else
-#include "itkTransformToDeformationFieldSource.h"
-#endif
+#include "clitkCoeffsToDVF.h"
 
 namespace clitk
 {
@@ -43,59 +37,6 @@ namespace clitk
       }
   }
 
-  template<class DisplacementFieldType>
-  typename DisplacementFieldType::Pointer ComposeVFGenericFilter::CoeffsToDVF(std::string fileName, std::string likeFileName)
-  {
-    typedef clitk::BSplineDeformableTransform<double, DisplacementFieldType::ImageDimension, DisplacementFieldType::ImageDimension> TransformType;
-    typedef typename TransformType::CoefficientImageType CoefficientImageType;
-
-    typedef itk::ImageFileReader<CoefficientImageType> CoeffReaderType;
-    typename CoeffReaderType::Pointer reader = CoeffReaderType::New();
-    reader->SetFileName(fileName);
-    reader->Update();
-
-    typename TransformType::Pointer transform = TransformType::New();
-    transform->SetCoefficientImage(reader->GetOutput());
-    
-#if ITK_VERSION_MAJOR >= 4
-        typedef itk::TransformToDisplacementFieldSource<DisplacementFieldType, double> ConvertorType;
-#else
-        typedef itk::TransformToDeformationFieldSource<DisplacementFieldType, double> ConvertorType;
-#endif
-
-    typedef itk::ImageIOBase ImageIOType;
-    typename ImageIOType::Pointer imageIO = itk::ImageIOFactory::CreateImageIO(likeFileName.c_str(), itk::ImageIOFactory::ReadMode);
-    imageIO->SetFileName(likeFileName);
-    imageIO->ReadImageInformation();
-
-    typename ConvertorType::Pointer convertor= ConvertorType::New();
-    typename ConvertorType::SizeType output_size;
-    typename ConvertorType::SpacingType output_spacing;
-    typename ConvertorType::OriginType output_origin;
-    typename ConvertorType::DirectionType output_direction;
-    for (unsigned int i = 0; i < DisplacementFieldType::ImageDimension; i++) {
-      output_size[i] = imageIO->GetDimensions(i);
-      output_spacing[i] = imageIO->GetSpacing(i);
-      output_origin[i] = imageIO->GetOrigin(i);
-      for (unsigned int j = 0; j < DisplacementFieldType::ImageDimension; j++)
-        output_direction[i][j] = imageIO->GetDirection(i)[j];
-    }
-    
-    if (m_Verbose) {
-      std::cout << "Interpolating coefficients with grid:" << std::endl;
-      std::cout << output_size << output_spacing << std::endl;
-    }
-    
-    convertor->SetNumberOfThreads(1);
-    convertor->SetTransform(transform);
-    convertor->SetOutputOrigin(output_origin);
-    convertor->SetOutputSpacing(output_spacing);
-    convertor->SetOutputSize(output_size);
-    convertor->SetOutputDirection(output_direction);
-    convertor->Update();
-
-    return convertor->GetOutput();
-  }
   
   template<unsigned int Dimension, class PixelType>
   void ComposeVFGenericFilter::UpdateWithDimAndPixelType()
@@ -106,8 +47,8 @@ namespace clitk
 
     //Define the image type
     if (m_Type == 1) {
-      input1 = this->CoeffsToDVF<ImageType>(m_InputName1, m_LikeImage);
-      input2 = this->CoeffsToDVF<ImageType>(m_InputName2, m_LikeImage);
+      input1 = CoeffsToDVF<ImageType>(m_InputName1, m_LikeImage, m_Verbose);
+      input2 = CoeffsToDVF<ImageType>(m_InputName2, m_LikeImage, m_Verbose);
     }
     else {
       //Read the input1
index cdaf86fcd900e726783a8dafbd933ec3bd7197b5..012524c8fc266d817b9cda408cde4259def6d5d9 100644 (file)
@@ -99,7 +99,6 @@ namespace clitk
     //----------------------------------------  
     template <unsigned int Dimension>  void UpdateWithDim(std::string PixelType);
     template <unsigned int Dimension, class PixelType>  void UpdateWithDimAndPixelType();
-    template<class DisplacementFieldType> typename DisplacementFieldType::Pointer CoeffsToDVF(std::string fileName, std::string likeFileName);
 
 
     //----------------------------------------  
index 5cd12c4ce202129c38ac67d59d5629fa590b58e1..e4b590e060b757c4d298dbabe7d0e33a88e7ec10 100644 (file)
 #define clitkInvertVFGenericFilter_txx
 
 #include "itkVectorResampleImageFilter.h"
-#include "clitkBSplineDeformableTransform.h"
-#if ITK_VERSION_MAJOR >= 4
-#include "itkTransformToDisplacementFieldSource.h"
-#else
-#include "itkTransformToDeformationFieldSource.h"
-#endif
+#include "clitkCoeffsToDVF.h"
 
 /* =================================================
  * @file   clitkInvertVFGenericFilter.txx
@@ -106,65 +101,6 @@ InvertVFGenericFilter<args_info_type>::UpdateWithDim(std::string PixelType)
   // }
 }
 
-//-------------------------------------------------------------------
-// Convert Coefficient image to DVF
-//-------------------------------------------------------------------
-template<class args_info_type>
-template<class DisplacementFieldType>
-typename DisplacementFieldType::Pointer
-InvertVFGenericFilter<args_info_type>::CoeffsToDVF(std::string fileName, std::string likeFileName)
-{
-  typedef clitk::BSplineDeformableTransform<double, DisplacementFieldType::ImageDimension, DisplacementFieldType::ImageDimension> TransformType;
-  typedef typename TransformType::CoefficientImageType CoefficientImageType;
-
-  typedef itk::ImageFileReader<CoefficientImageType> CoeffReaderType;
-  typename CoeffReaderType::Pointer reader = CoeffReaderType::New();
-  reader->SetFileName(fileName);
-  reader->Update();
-
-  typename TransformType::Pointer transform = TransformType::New();
-  transform->SetCoefficientImage(reader->GetOutput());
-  
-#if ITK_VERSION_MAJOR >= 4
-      typedef itk::TransformToDisplacementFieldSource<DisplacementFieldType, double> ConvertorType;
-#else
-      typedef itk::TransformToDeformationFieldSource<DisplacementFieldType, double> ConvertorType;
-#endif
-
-  typedef itk::ImageIOBase ImageIOType;
-  typename ImageIOType::Pointer imageIO = itk::ImageIOFactory::CreateImageIO(likeFileName.c_str(), itk::ImageIOFactory::ReadMode);
-  imageIO->SetFileName(likeFileName);
-  imageIO->ReadImageInformation();
-
-  typename ConvertorType::Pointer convertor= ConvertorType::New();
-  typename ConvertorType::SizeType output_size;
-  typename ConvertorType::SpacingType output_spacing;
-  typename ConvertorType::OriginType output_origin;
-  typename ConvertorType::DirectionType output_direction;
-  for (unsigned int i = 0; i < DisplacementFieldType::ImageDimension; i++) {
-    output_size[i] = imageIO->GetDimensions(i);
-    output_spacing[i] = imageIO->GetSpacing(i);
-    output_origin[i] = imageIO->GetOrigin(i);
-    for (unsigned int j = 0; j < DisplacementFieldType::ImageDimension; j++)
-      output_direction[i][j] = imageIO->GetDirection(i)[j];
-  }
-  
-  if (m_Verbose) {
-    std::cout << "Interpolating coefficients with grid:" << std::endl;
-    std::cout << output_size << output_spacing << std::endl;
-  }
-  
-  convertor->SetNumberOfThreads(1);
-  convertor->SetTransform(transform);
-  convertor->SetOutputOrigin(output_origin);
-  convertor->SetOutputSpacing(output_spacing);
-  convertor->SetOutputSize(output_size);
-  convertor->SetOutputDirection(output_direction);
-  convertor->Update();
-
-  return convertor->GetOutput();
-}
-
 
 //-------------------------------------------------------------------
 // Update with the number of dimensions and the pixeltype
index 9394d8d5798d06ea2af9489822e37edf1716610f..0ce874deb74c26f8773f7f0438c7ce1efecd5fa0 100644 (file)
@@ -98,7 +98,6 @@ namespace clitk
     //----------------------------------------  
     template <unsigned int Dimension>  void UpdateWithDim(std::string PixelType);
     template <unsigned int Dimension, class PixelType>  void UpdateWithDimAndPixelType();
-    template<class DisplacementFieldType> typename DisplacementFieldType::Pointer CoeffsToDVF(std::string fileName, std::string likeFileName);
 
     //----------------------------------------  
     // Data members
index ce979aa330690a1e08375cfca90024cfefbf208e..c64080e2fb23a2748505e7962ffaf02afce1e8a1 100644 (file)
  ===================================================*/
 
 #include "itkVectorResampleImageFilter.h"
-#include "clitkBSplineDeformableTransform.h"
-#if ITK_VERSION_MAJOR >= 4
-#include "itkTransformToDisplacementFieldSource.h"
-#else
-#include "itkTransformToDeformationFieldSource.h"
-#endif
+#include "clitkCoeffsToDVF.h"
 
 namespace clitk
 {
@@ -71,63 +66,6 @@ WarpImageGenericFilter::UpdateWithDim(std::string PixelType)
   }
 }
 
-//-------------------------------------------------------------------
-// Convert Coefficient image to DVF
-//-------------------------------------------------------------------
-template<class DisplacementFieldType>
-typename DisplacementFieldType::Pointer
-WarpImageGenericFilter::CoeffsToDVF(std::string fileName, std::string likeFileName)
-{
-  typedef clitk::BSplineDeformableTransform<double, DisplacementFieldType::ImageDimension, DisplacementFieldType::ImageDimension> TransformType;
-  typedef typename TransformType::CoefficientImageType CoefficientImageType;
-
-  typedef itk::ImageFileReader<CoefficientImageType> CoeffReaderType;
-  typename CoeffReaderType::Pointer reader = CoeffReaderType::New();
-  reader->SetFileName(fileName);
-  reader->Update();
-
-  typename TransformType::Pointer transform = TransformType::New();
-  transform->SetCoefficientImage(reader->GetOutput());
-  
-#if ITK_VERSION_MAJOR >= 4
-      typedef itk::TransformToDisplacementFieldSource<DisplacementFieldType, double> ConvertorType;
-#else
-      typedef itk::TransformToDeformationFieldSource<DisplacementFieldType, double> ConvertorType;
-#endif
-
-  typedef itk::ImageIOBase ImageIOType;
-  typename ImageIOType::Pointer imageIO = itk::ImageIOFactory::CreateImageIO(likeFileName.c_str(), itk::ImageIOFactory::ReadMode);
-  imageIO->SetFileName(likeFileName);
-  imageIO->ReadImageInformation();
-
-  typename ConvertorType::Pointer convertor= ConvertorType::New();
-  typename ConvertorType::SizeType output_size;
-  typename ConvertorType::SpacingType output_spacing;
-  typename ConvertorType::OriginType output_origin;
-  typename ConvertorType::DirectionType output_direction;
-  for (unsigned int i = 0; i < DisplacementFieldType::ImageDimension; i++) {
-    output_size[i] = imageIO->GetDimensions(i);
-    output_spacing[i] = imageIO->GetSpacing(i);
-    output_origin[i] = imageIO->GetOrigin(i);
-    for (unsigned int j = 0; j < DisplacementFieldType::ImageDimension; j++)
-      output_direction[i][j] = imageIO->GetDirection(i)[j];
-  }
-  
-  if (m_Verbose) {
-    std::cout << "Interpolating coefficients with grid:" << std::endl;
-    std::cout << output_size << output_spacing << std::endl;
-  }
-  
-  convertor->SetNumberOfThreads(1);
-  convertor->SetTransform(transform);
-  convertor->SetOutputOrigin(output_origin);
-  convertor->SetOutputSpacing(output_spacing);
-  convertor->SetOutputSize(output_size);
-  convertor->SetOutputDirection(output_direction);
-  convertor->Update();
-
-  return convertor->GetOutput();
-}
 
 //-------------------------------------------------------------------
 // Update with the number of dimensions and the pixeltype
@@ -152,7 +90,7 @@ WarpImageGenericFilter::UpdateWithDimAndPixelType()
   
   typename DeformationFieldType::Pointer deformationField;
   if (m_ArgsInfo.coeff_given) {
-    deformationField = CoeffsToDVF<DeformationFieldType>(m_ArgsInfo.coeff_arg, m_InputFileName);
+    deformationField = CoeffsToDVF<DeformationFieldType>(m_ArgsInfo.coeff_arg, m_InputFileName, m_Verbose);
   }
   else {
     //Read the deformation field