]> Creatis software - clitk.git/blobdiff - tools/clitkSUVPeakGenericFilter.txx
Add option to define the volume of the filter in clitkSUVPeak in cc
[clitk.git] / tools / clitkSUVPeakGenericFilter.txx
index 8a41726c87ee40d990c6af82065950744a041744..893750dd3961b0d888ad4327346d7abfc66cccbe 100644 (file)
@@ -19,6 +19,9 @@
 #define CLITKSUVPEAKGENERICFILTER_TXX
 
 #include "clitkImageCommon.h"
+#include "clitkCropLikeImageFilter.h"
+#include "clitkResampleImageWithOptionsFilter.h"
+
 #include <itkConvolutionImageFilter.h>
 
 namespace clitk
@@ -54,9 +57,9 @@ void SUVPeakGenericFilter<args_info_type>::SetArgsInfo(const args_info_type & a)
 
   // Set value
   this->SetIOVerbose(mArgsInfo.verbose_flag);
-  
+
   if (mArgsInfo.input_given) this->AddInputFilename(mArgsInfo.input_arg);
-  
+
   if (mArgsInfo.mask_given)  this->AddInputFilename(mArgsInfo.mask_arg);
   }
 //--------------------------------------------------------------------
@@ -75,6 +78,34 @@ void SUVPeakGenericFilter<args_info_type>::UpdateWithInputImageType()
   typename MaskImageType::Pointer mask;
   if(mArgsInfo.mask_given) {
       mask = this->template GetInput<MaskImageType>(1);
+      // Check mask sampling/size
+      if (!HaveSameSizeAndSpacing<MaskImageType, ImageType>(mask, input)) {
+        if (mArgsInfo.allow_resize_flag) {
+          if (mArgsInfo.verbose_flag) {
+            std::cout << "Resize mask image like input" << std::endl;
+          }
+          typedef clitk::ResampleImageWithOptionsFilter<MaskImageType> ResamplerType;
+          typename ResamplerType::Pointer resampler = ResamplerType::New();
+          resampler->SetInput(mask); //By default the interpolation in NN, Ok for mask
+          resampler->SetOutputSpacing(input->GetSpacing());
+          resampler->SetOutputOrigin(mask->GetOrigin());
+          resampler->SetGaussianFilteringEnabled(false);
+          resampler->Update();
+          mask = resampler->GetOutput();
+
+          typedef clitk::CropLikeImageFilter<MaskImageType> FilterType;
+          typename FilterType::Pointer crop = FilterType::New();
+          crop->SetInput(mask);
+          crop->SetCropLikeImage(input);
+          crop->Update();
+          mask = crop->GetOutput();
+
+        }
+        else {
+          std::cerr << "Mask image has a different size/spacing than input. Abort. (Use option --allow_resize)" << std::endl;
+          exit(-1);
+        }
+      }
   }
   else {
       mask = MaskImageType::New();
@@ -86,9 +117,11 @@ void SUVPeakGenericFilter<args_info_type>::UpdateWithInputImageType()
   }
 
   double volume = 1000; //1 cc into mc
+  if (mArgsInfo.volume_given)
+    volume *= mArgsInfo.volume_arg;
   const double PI = 3.141592653589793238463;
   double radius = std::pow(3*volume/(4*PI),1./3);
-  
+
   typename ImageType::Pointer kernel = ComputeMeanFilterKernel<ImageType>(input->GetSpacing(), radius);
 
   // Perform the convolution
@@ -98,8 +131,8 @@ void SUVPeakGenericFilter<args_info_type>::UpdateWithInputImageType()
   filter->SetKernelImage(kernel);
   filter->Update();
   typename ImageType::Pointer output = filter->GetOutput();
-  
-  
+
+
   typedef itk::ImageRegionConstIteratorWithIndex<ImageType> IteratorType;
   typedef itk::ImageRegionConstIteratorWithIndex<MaskImageType> MIteratorType;
   IteratorType iters(output, output->GetLargestPossibleRegion());
@@ -120,7 +153,7 @@ void SUVPeakGenericFilter<args_info_type>::UpdateWithInputImageType()
   }
   typename ImageType::PointType p;
   output->TransformIndexToPhysicalPoint(index, p);
-  std::cout<<"SUV Peak found in "<< p << " with the value " << max << std::endl;
+  std::cout<<"SUV Peak found in "<< p << " mm with the value " << max << std::endl;
 }
 //--------------------------------------------------------------------
 
@@ -128,7 +161,7 @@ void SUVPeakGenericFilter<args_info_type>::UpdateWithInputImageType()
 //--------------------------------------------------------------------
 template<class args_info_type>
 template<class ImageType>
-typename ImageType::Pointer 
+typename ImageType::Pointer
 SUVPeakGenericFilter<args_info_type>::ComputeMeanFilterKernel(const typename ImageType::SpacingType & spacing, double radius)
 {
   // Some kind of cache to speed up a bit