]> Creatis software - clitk.git/blobdiff - tools/clitkWarpImageGenericFilter.txx
With ITKv5, change VectorResample and VectorCast Image Filter to Resample and Cast...
[clitk.git] / tools / clitkWarpImageGenericFilter.txx
index c64080e2fb23a2748505e7962ffaf02afce1e8a1..c3f58b37d225dc41181009664831ea7f1f093808 100644 (file)
  *
  ===================================================*/
 
+#if ( ITK_VERSION_MAJOR < 5 )
 #include "itkVectorResampleImageFilter.h"
-#include "clitkCoeffsToDVF.h"
+#else
+#include "itkResampleImageFilter.h"
+#endif
+#include "clitkConvertBLUTCoeffsToVFFilter.h"
 
 namespace clitk
 {
@@ -90,7 +94,13 @@ WarpImageGenericFilter::UpdateWithDimAndPixelType()
   
   typename DeformationFieldType::Pointer deformationField;
   if (m_ArgsInfo.coeff_given) {
-    deformationField = CoeffsToDVF<DeformationFieldType>(m_ArgsInfo.coeff_arg, m_InputFileName, m_Verbose);
+    typedef ConvertBLUTCoeffsToVFFilter<DeformationFieldType> FilterType;
+    typename FilterType::Pointer filter = FilterType::New();
+    filter->SetInputFileName(m_ArgsInfo.coeff_arg);
+    filter->SetLikeFileName(m_InputFileName);
+    filter->SetVerbose(m_Verbose);
+    filter->Update();
+    deformationField = filter->GetOutput();
   }
   else {
     //Read the deformation field
@@ -122,10 +132,16 @@ WarpImageGenericFilter::UpdateWithDimAndPixelType()
     genericInterpolator->SetArgsInfo(m_ArgsInfo);
 
     // Resample to match the extent of the input
+#if ( ITK_VERSION_MAJOR < 5 )
     typename itk::VectorResampleImageFilter<DeformationFieldType, DeformationFieldType >::Pointer
     resampler =itk::VectorResampleImageFilter<DeformationFieldType, DeformationFieldType >::New();
+#else
+    typename itk::ResampleImageFilter<DeformationFieldType, DeformationFieldType >::Pointer
+    resampler =itk::ResampleImageFilter<DeformationFieldType, DeformationFieldType >::New();
+#endif
     resampler->SetInput(deformationField);
     resampler->SetOutputSpacing(deformationField->GetSpacing());
+    resampler->SetOutputDirection(deformationField->GetDirection());
     resampler->SetSize(newSize);
     resampler->SetOutputOrigin(input->GetOrigin());
     resampler->SetInterpolator(genericInterpolator->GetInterpolatorPointer());
@@ -153,10 +169,16 @@ WarpImageGenericFilter::UpdateWithDimAndPixelType()
       newSize[i]=input->GetLargestPossibleRegion().GetSize()[i];
 
     // Resample to match the extent of the input
+#if ( ITK_VERSION_MAJOR < 5 )
     typename itk::VectorResampleImageFilter<DeformationFieldType, DeformationFieldType >::Pointer
     resampler =itk::VectorResampleImageFilter<DeformationFieldType, DeformationFieldType >::New();
+#else
+    typename itk::ResampleImageFilter<DeformationFieldType, DeformationFieldType >::Pointer
+    resampler =itk::ResampleImageFilter<DeformationFieldType, DeformationFieldType >::New();
+#endif
     resampler->SetInput(deformationField);
     resampler->SetOutputSpacing(input->GetSpacing());
+    resampler->SetOutputDirection(deformationField->GetDirection());
     resampler->SetSize(newSize);
     resampler->SetOutputOrigin(input->GetOrigin());
     resampler->SetInterpolator(genericInterpolator->GetInterpolatorPointer());
@@ -199,17 +221,19 @@ WarpImageGenericFilter::UpdateWithDimAndPixelType()
     //Backward mapping
     typedef itk::WarpImageFilter<InputImageType, InputImageType, DeformationFieldType> BackwardWarpFilterType;
     typename BackwardWarpFilterType::Pointer backwardWarpFilter= BackwardWarpFilterType::New();
-#if ITK_VERSION_MAJOR >= 4
     backwardWarpFilter->SetDisplacementField( deformationField );
-#else
-    backwardWarpFilter->SetDeformationField( deformationField );
-#endif
     backwardWarpFilter->SetEdgePaddingValue( static_cast<PixelType>(m_ArgsInfo.pad_arg) );
     backwardWarpFilter->SetOutputSpacing( deformationField->GetSpacing() );
     backwardWarpFilter->SetOutputOrigin( input->GetOrigin() );
     backwardWarpFilter->SetOutputSize( deformationField->GetLargestPossibleRegion().GetSize() );
+    backwardWarpFilter->SetOutputDirection( input->GetDirection() );
+#if ( ITK_VERSION_MAJOR < 5 )
     typename itk::VectorResampleImageFilter<DeformationFieldType, DeformationFieldType >::Pointer
     resampler =itk::VectorResampleImageFilter<DeformationFieldType, DeformationFieldType >::New();
+#else
+    typename itk::ResampleImageFilter<DeformationFieldType, DeformationFieldType >::Pointer
+    resampler =itk::ResampleImageFilter<DeformationFieldType, DeformationFieldType >::New();
+#endif
     backwardWarpFilter->SetInterpolator(genericInterpolator->GetInterpolatorPointer());
     warpFilter=backwardWarpFilter;
   }