]> Creatis software - clitk.git/blobdiff - common/clitkImage2DicomRTStructFilter.txx
Revert previsous commit with Roi name
[clitk.git] / common / clitkImage2DicomRTStructFilter.txx
index 1b241a99824a233a094058f1b8196d09c4e99772..eba1251cfa0377927923016446135b22be4fda3f 100644 (file)
 #include "vtkCamera.h"
 #include "vtkProperty.h"
 #include "vtkProperty2D.h"
+#include <vtksys/SystemTools.hxx>
 
 // itk
 #include <itkImage.h>
 #include <itkVTKImageToImageFilter.h>
 
 // gdcm
+#include <vtkRTStructSetProperties.h>
+#include <vtkGDCMPolyDataReader.h>
 #include <vtkGDCMPolyDataWriter.h>
 
 //--------------------------------------------------------------------
@@ -61,6 +64,7 @@ clitk::Image2DicomRTStructFilter<PixelType>::Image2DicomRTStructFilter()
   m_DicomFolder = "";
   m_OutputFilename = "default-output.dcm";
   m_ThresholdValue = 0.5;
+  m_SkipInitialStructuresFlag = false;
 }
 //--------------------------------------------------------------------
 
@@ -76,9 +80,8 @@ clitk::Image2DicomRTStructFilter<PixelType>::~Image2DicomRTStructFilter()
 //--------------------------------------------------------------------
 template<class PixelType>
 void
-clitk::Image2DicomRTStructFilter<PixelType>::SetROIName(std::string name, std::string type)
+clitk::Image2DicomRTStructFilter<PixelType>::SetROIType(std::string type)
 {
-  m_ROIName = name;
   m_ROIType = type;
 }
 //--------------------------------------------------------------------
@@ -111,10 +114,18 @@ void clitk::Image2DicomRTStructFilter<PixelType>::Update()
     std::cout << "Number of structures in the dicom-rt-struct : " 
               << p->GetNumberOfStructureSetROIs() << std::endl;
   }
-  
+
+  // number of additional contours
+  int m = m_InputFilenames.size();
+
   // Init writer
   vtkGDCMPolyDataWriter * writer = vtkGDCMPolyDataWriter::New();
-  int numMasks = reader->GetNumberOfOutputPorts() + 1;//add one more
+  int numMasks = reader->GetNumberOfOutputPorts() + m;
+
+  if (m_SkipInitialStructuresFlag) {
+    numMasks = m;
+  }
+
   writer->SetNumberOfInputPorts(numMasks);    
   writer->SetFileName(m_OutputFilename.c_str());
   writer->SetMedicalImageProperties(reader->GetMedicalImageProperties());
@@ -126,35 +137,72 @@ void clitk::Image2DicomRTStructFilter<PixelType>::Update()
   roiNames->SetNumberOfValues(numMasks);
   roiAlgorithms->SetNumberOfValues(numMasks);
   roiTypes->SetNumberOfValues(numMasks);
-  
+
   // Convert the image into a mesh
-  typedef clitk::BinaryImageToMeshFilter<ImageType> BinaryImageToMeshFilterType;
-  typename BinaryImageToMeshFilterType::Pointer convert = BinaryImageToMeshFilterType::New();
-  convert->SetThresholdValue(m_ThresholdValue);
-  convert->SetInput(m_Input);
-  convert->Update();
-  vtkPolyData* mesh = convert->GetOutputMesh();
-  if (GetVerboseFlag()) {
-    std::cout << "Mesh has " << mesh->GetNumberOfLines() << " lines." << std::endl;
+  std::vector<vtkSmartPointer<vtkPolyData> > meshes;
+  std::vector<std::string> m_ROINames;
+  meshes.resize(m);
+  m_ROINames.resize(m);
+  for(unsigned int i=0; i<m; i++) {
+
+    // read image
+    //    typedef float PixelType;
+    //typedef itk::Image<PixelType, 3> ImageType;
+    ImagePointer input = clitk::readImage<ImageType>(m_InputFilenames[i], false);
+
+    std::ostringstream oss;
+    oss << vtksys::SystemTools::
+      GetFilenameName(vtksys::SystemTools::GetFilenameWithoutLastExtension(m_InputFilenames[i]));
+    std::string name = oss.str();
+    m_ROINames[i] = name;
+
+    // convert to mesh
+    typedef clitk::BinaryImageToMeshFilter<ImageType> BinaryImageToMeshFilterType;
+    typename BinaryImageToMeshFilterType::Pointer convert = BinaryImageToMeshFilterType::New();
+    convert->SetThresholdValue(m_ThresholdValue);
+    convert->SetInput(input);
+    convert->Update();
+    meshes[i] = convert->GetOutputMesh();
+    if (GetVerboseFlag()) {
+      std::cout << "Mesh has " << meshes[i]->GetNumberOfLines() << " lines." << std::endl;
+    }
+
+    /*
+    // debug mesh write  FIXME
+    vtkSmartPointer<vtkPolyDataWriter> wr = vtkSmartPointer<vtkPolyDataWriter>::New();
+    wr->SetInputConnection(convert->GetOutputPort()); //psurface->GetOutputPort()
+    wr->SetFileName("bidon.obj");
+    wr->Update();
+    wr->Write();
+    */
   }
 
   // Copy previous contours
-  for (unsigned int i = 0; i < numMasks-1; ++i) {
+  for (unsigned int i = 0; i < numMasks-m; ++i) {
+#if VTK_MAJOR_VERSION <= 5
     writer->SetInput(i, reader->GetOutput(i));
+#else
+    writer->SetInputData(i, reader->GetOutput(i));
+#endif
     std::string theString = reader->GetRTStructSetProperties()->GetStructureSetROIName(i);
     roiNames->InsertValue(i, theString);
     theString = reader->GetRTStructSetProperties()->GetStructureSetROIGenerationAlgorithm(i);
     roiAlgorithms->InsertValue(i, theString);
     theString = reader->GetRTStructSetProperties()->GetStructureSetRTROIInterpretedType(i);
     roiTypes->InsertValue(i, theString);
-  }  
+  }
 
-  // Add new one
-  int last = numMasks-1;
-  writer->SetInput(last, mesh);
-  roiNames->InsertValue(last, m_ROIName);
-  roiAlgorithms->InsertValue(last, "CLITK_CREATED");
-  roiTypes->InsertValue(last, m_ROIType);
+  // Add new ones
+  for (unsigned int i = numMasks-m; i < numMasks; ++i) {
+#if VTK_MAJOR_VERSION <= 5
+    writer->SetInput(i, meshes[i-numMasks+m]);
+#else
+    writer->SetInputData(i, meshes[i-numMasks+m]);
+#endif
+    roiNames->InsertValue(i, m_ROIType);
+    roiAlgorithms->InsertValue(i, "CLITK_CREATED");
+    roiTypes->InsertValue(i, m_ROIType);
+  }
 
   /*
   //  Visu DEBUG
@@ -194,7 +242,7 @@ void clitk::Image2DicomRTStructFilter<PixelType>::Update()
   writer->InitializeRTStructSet(m_DicomFolder,
                                 reader->GetRTStructSetProperties()->GetStructureSetLabel(),
                                 reader->GetRTStructSetProperties()->GetStructureSetName(),
-                                roiNames, roiAlgorithms, roiTypes);  
+                                roiNames, roiAlgorithms, roiTypes);
   writer->Write();
   reader->Delete();
   roiNames->Delete();