]> Creatis software - clitk.git/blobdiff - common/clitkDicomRT_ROI_ConvertToImageFilter.cxx
Fix for bug Bug #670: closing and opening do the same as successive
[clitk.git] / common / clitkDicomRT_ROI_ConvertToImageFilter.cxx
index 969463990b6fb71cb2534590299a9d254613cbb5..e3b34723fdc264865ee925ddac9cd549b836bd67 100644 (file)
@@ -4,7 +4,7 @@
 
   Authors belongs to:
   - University of LYON           http://www.universite-lyon.fr/
-  - Léon Bérard cancer center    http://oncora1.lyon.fnclcc.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
 
   =========================================================================*/
 
+#include <iterator>
+#include <algorithm>
+
+// clitk
 #include "clitkDicomRT_ROI_ConvertToImageFilter.h"
+#include "clitkImageCommon.h"
+
+// vtk
 #include <vtkPolyDataToImageStencil.h>
 #include <vtkSmartPointer.h>
 #include <vtkImageStencil.h>
 #include <vtkLinearExtrusionFilter.h>
-#include <itkVTKImageToImageFilter.h>
-#include "clitkImageCommon.h"
+#include <vtkMetaImageWriter.h>
+
 
 //--------------------------------------------------------------------
 clitk::DicomRT_ROI_ConvertToImageFilter::DicomRT_ROI_ConvertToImageFilter()
 {
   mROI = NULL;
-  mImageInfoIsSet = false;
   mWriteOutput = false;
   mCropMask = true;
 }
@@ -43,6 +49,10 @@ clitk::DicomRT_ROI_ConvertToImageFilter::~DicomRT_ROI_ConvertToImageFilter()
 }
 //--------------------------------------------------------------------
 
+bool clitk::DicomRT_ROI_ConvertToImageFilter::ImageInfoIsSet() const
+{
+  return mSize.size() && mSpacing.size() && mOrigin.size();
+}
 
 //--------------------------------------------------------------------
 void clitk::DicomRT_ROI_ConvertToImageFilter::SetROI(clitk::DicomRT_ROI * roi)
@@ -88,10 +98,23 @@ void clitk::DicomRT_ROI_ConvertToImageFilter::SetImageFilename(std::string f)
     mOrigin[i] = header->GetOrigin(i);
     mSize[i] = header->GetDimensions(i);
   }
-  mImageInfoIsSet = true;
 }
 //--------------------------------------------------------------------
 
+void clitk::DicomRT_ROI_ConvertToImageFilter::SetOutputOrigin(const double* origin)
+{
+  std::copy(origin,origin+3,std::back_inserter(mOrigin));
+}
+//--------------------------------------------------------------------
+void clitk::DicomRT_ROI_ConvertToImageFilter::SetOutputSpacing(const double* spacing)
+{
+  std::copy(spacing,spacing+3,std::back_inserter(mSpacing));
+}
+//--------------------------------------------------------------------
+void clitk::DicomRT_ROI_ConvertToImageFilter::SetOutputSize(const unsigned long* size)
+{
+  std::copy(size,size+3,std::back_inserter(mSize));
+}
 
 //--------------------------------------------------------------------
 void clitk::DicomRT_ROI_ConvertToImageFilter::Update()
@@ -100,15 +123,13 @@ void clitk::DicomRT_ROI_ConvertToImageFilter::Update()
     std::cerr << "Error. No ROI set, please use SetROI." << std::endl;
     exit(0);
   }
-  if (!mImageInfoIsSet) {
+  if (!ImageInfoIsSet()) {
     std::cerr << "Error. Please provide image info (spacing/origin) with SetImageFilename" << std::endl;
     exit(0);
   }
-  // DD("Update");
 
   // Get Mesh
   vtkPolyData * mesh = mROI->GetMesh();
-  DD(mesh->GetNumberOfCells());
 
   // Get bounds
   double *bounds=mesh->GetBounds();
@@ -139,7 +160,7 @@ void clitk::DicomRT_ROI_ConvertToImageFilter::Update()
   }
 
   // Create new output image
-  mBinaryImage = vtkImageData::New();
+  mBinaryImage = vtkSmartPointer<vtkImageData>::New();
   mBinaryImage->SetScalarTypeToUnsignedChar();
   mBinaryImage->SetOrigin(&origin[0]);
   mBinaryImage->SetSpacing(&mSpacing[0]);
@@ -159,7 +180,7 @@ void clitk::DicomRT_ROI_ConvertToImageFilter::Update()
   // Extrude
   vtkSmartPointer<vtkLinearExtrusionFilter> extrude=vtkSmartPointer<vtkLinearExtrusionFilter>::New();
   extrude->SetInput(mesh);
-  ///We extrude in the -slice_spacing direction to respect the FOCAL convention // ?????????????
+  ///We extrude in the -slice_spacing direction to respect the FOCAL convention (NEEDED !)
   extrude->SetVector(0, 0, -mSpacing[2]);
 
   // Binarization
@@ -176,6 +197,14 @@ void clitk::DicomRT_ROI_ConvertToImageFilter::Update()
   stencil->SetInput(mBinaryImage);
   stencil->ReverseStencilOn();
   stencil->Update();
+
+  /*
+  vtkSmartPointer<vtkMetaImageWriter> w = vtkSmartPointer<vtkMetaImageWriter>::New();
+  w->SetInput(stencil->GetOutput());
+  w->SetFileName("binary2.mhd");
+  w->Write();
+  */
+
   mBinaryImage->ShallowCopy(stencil->GetOutput());
 
   if (mWriteOutput) {
@@ -194,6 +223,8 @@ void clitk::DicomRT_ROI_ConvertToImageFilter::Update()
 //--------------------------------------------------------------------
 vtkImageData * clitk::DicomRT_ROI_ConvertToImageFilter::GetOutput()
 {
+  assert(mBinaryImage);
   return mBinaryImage;
 }
 //--------------------------------------------------------------------
+