]> Creatis software - clitk.git/commitdiff
Add pixel coordinate instead of mm coordinates into Xml2DicomRTStruct
authortbaudier <thomas.baudier@creatis.insa-lyon.fr>
Fri, 8 Mar 2019 08:26:09 +0000 (09:26 +0100)
committertbaudier <thomas.baudier@creatis.insa-lyon.fr>
Fri, 8 Mar 2019 08:26:09 +0000 (09:26 +0100)
For some CT with transformation matrix, it's better to convert Xml struct using pixel coordinate. In such case, I need the path to the dicom CT

common/clitkXml2DicomRTStructFilter.h
common/clitkXml2DicomRTStructFilter.txx
tools/clitkXml2DicomRTStruct.cxx
tools/clitkXml2DicomRTStruct.ggo

index 32027169110e12027a2edfb9e3b5d6e6c166dfbb..cd13a4aaab4eb9735e1a5959f453dc6df03e8fe8 100644 (file)
@@ -40,6 +40,8 @@ namespace clitk {
     typedef typename ImageType::Pointer ImagePointer;
     typedef typename clitk::DicomRT_StructureSet::Pointer DicomRTStructPointer;
 
+    itkSetMacro(UsePixel, bool);
+    itkSetMacro(ImageMHD, std::string);
     itkSetMacro(InputFilename, std::string);
     itkSetMacro(StructureSetFilename, std::string);
     itkSetMacro(DicomFolder, std::string);
@@ -49,6 +51,8 @@ namespace clitk {
     void Update();
 
   protected:
+    bool m_UsePixel;
+    std::string m_ImageMHD;
     std::string m_StructureSetFilename;
     std::string m_DicomFolder;
     std::string m_OutputFilename;
index 9c106ec4ac8a49f76c613a1cf73a03b1e260f391..4945a55675d1803a402b87f7b77ea12719e467d0 100644 (file)
 // xml parser
 #include "../utilities/pugixml/pugixml.hpp"
 
+// itk
+#include "itkImage.h"
+#include "itkImageFileReader.h"
+
 // vtk
 #include <vtkSmartPointer.h>
 #include <vtkMetaImageWriter.h>
 #include <vtkImageData.h>
 #include <vtkPolygon.h>
 #include <vtkAppendPolyData.h>
+#include <vtkPointData.h>
 
 // FIXME to remove
 #include "vtkPolyDataMapper.h"
@@ -51,6 +56,7 @@
 template<class PixelType>
 clitk::Xml2DicomRTStructFilter<PixelType>::Xml2DicomRTStructFilter()
 {
+  m_ImageMHD = "";
   m_StructureSetFilename = "";
   m_DicomFolder = "";
   m_OutputFilename = "default-output.dcm";
@@ -87,6 +93,20 @@ void clitk::Xml2DicomRTStructFilter<PixelType>::Update()
   reader->SetFileName(m_StructureSetFilename.c_str());
   reader->Update();
 
+  //Check id use pixel is set and if imageMHD is present in such case
+  typedef itk::Image<double, 3> Input3DType;
+  typedef itk::ImageFileReader<Input3DType> InputReader3DType;
+  InputReader3DType::Pointer mhdReader = InputReader3DType::New();
+  if (m_UsePixel) {
+    if (m_ImageMHD != "") {
+      mhdReader->SetFileName(m_ImageMHD);
+      mhdReader->Update();
+    } else {
+      std::cout << "Set MHD image (option -j) when pixel is set" << std::endl;
+      return;
+    }
+  }
+
   // Get properties
   vtkRTStructSetProperties * p = reader->GetRTStructSetProperties();
   if (GetVerboseFlag()) {
@@ -107,6 +127,14 @@ void clitk::Xml2DicomRTStructFilter<PixelType>::Update()
   pugi::xml_node mainDict = doc.child("plist").child("dict").child("array");
   for (pugi::xml_node_iterator mainDictIt = mainDict.begin(); mainDictIt != mainDict.end(); ++mainDictIt) { //Look for all slice (one slice per dict)
     if (!std::strcmp(mainDictIt->name(), "dict")) {
+      pugi::xml_node_iterator sliceDescriptionIt = mainDictIt->begin();
+      while (sliceDescriptionIt != mainDictIt->end() && std::abs(std::strcmp(sliceDescriptionIt->child_value(), "ImageIndex"))) //Look for the ImageIndex for the current slice
+        ++sliceDescriptionIt;
+      int sliceNb = 0;
+      if (sliceDescriptionIt != mainDictIt->end()) { //take the following node to have the values
+        ++sliceDescriptionIt;
+        sliceNb = std::atoi(sliceDescriptionIt->child_value());
+      }
       for (pugi::xml_node_iterator sliceIt = mainDictIt->child("array").begin(); sliceIt != mainDictIt->child("array").end(); ++sliceIt) { //Look for all struct in the current slice (one struct per dict)
         if (!std::strcmp(sliceIt->name(), "dict")) {
           pugi::xml_node_iterator structureIt = sliceIt->begin();
@@ -115,15 +143,20 @@ void clitk::Xml2DicomRTStructFilter<PixelType>::Update()
           if (structureIt != sliceIt->end()) { //take the following node to have the values
             ++structureIt;
             std::string name(structureIt->child_value());
-            while (structureIt != sliceIt->end() && std::abs(std::strcmp(structureIt->child_value(), "Point_mm"))) //Look for the Point_mm for the current struct
-              ++structureIt;
+            if (m_UsePixel) {
+              while (structureIt != sliceIt->end() && std::abs(std::strcmp(structureIt->child_value(), "Point_px"))) //Look for the Point_px for the current struct
+                ++structureIt;
+            } else {
+              while (structureIt != sliceIt->end() && std::abs(std::strcmp(structureIt->child_value(), "Point_mm"))) //Look for the Point_mm for the current struct
+                ++structureIt;
+            }
             if (structureIt != sliceIt->end()) { //take the following node to have the values
               ++structureIt;
               // Insert all points for 1 struct into vtkPoints
               vtkSmartPointer<vtkPoints> points = vtkSmartPointer<vtkPoints>::New();
               for (pugi::xml_node_iterator pointIt = structureIt->begin(); pointIt != structureIt->end(); ++pointIt) { //Look for all points in mm in the current struct (one point per string)
                 if (!std::strcmp(pointIt->name(), "string")) {
-                  // Split the string to save the 3 values (remove the useless char) as double
+                  // Split the string to save the 2 or 3 values (remove the useless char) as double
                   char* copy = (char*)pointIt->child_value();
                   std::vector<char*> v;
                   char* chars_array = strtok(copy, ", ");
@@ -131,12 +164,28 @@ void clitk::Xml2DicomRTStructFilter<PixelType>::Update()
                       v.push_back(chars_array);
                       chars_array = strtok(NULL, ", ");
                   }
-                  v[0] = v[0] + 1;
-                  v[2][strlen(v[2])-1] = '\0';
                   double x, y, z;
-                  x = std::atof(v[0]);
-                  y = std::atof(v[1]);
-                  z = std::atof(v[2]);
+                  v[0] = v[0] + 1;
+                  if (m_UsePixel) {
+                    v[1][strlen(v[1])-1] = '\0';
+                    x = std::atof(v[0]);
+                    y = std::atof(v[1]);
+                    z = sliceNb;
+                    Input3DType::PointType position;
+                    Input3DType::IndexType index;
+                    index[0] = x;
+                    index[1] = y;
+                    index[2] = z;
+                    mhdReader->GetOutput()->TransformIndexToPhysicalPoint(index, position);
+                    x = position[0];
+                    y = position[1];
+                    z = position[2];
+                  } else {
+                    v[2][strlen(v[2])-1] = '\0';
+                    x = std::atof(v[0]);
+                    y = std::atof(v[1]);
+                    z = std::atof(v[2]);
+                  }
                   points->InsertNextPoint(x, y, z);
                 }
               }
@@ -207,8 +256,12 @@ void clitk::Xml2DicomRTStructFilter<PixelType>::Update()
   /*
   //  Visu DEBUG
   vtkPolyDataMapper *cubeMapper = vtkPolyDataMapper::New();
-  cubeMapper->SetInputData( mapName2Data[0]->GetOutput() );
-  cubeMapper->SetScalarRange(0,7);
+  cubeMapper->SetInputData( mapName2Data.begin()->second->GetOutput() );
+  double range[2];
+  mapName2Data.begin()->second->GetOutput()->GetPointData()->GetScalars()->GetRange(range);
+  std::cout << numMasks << std::endl;
+  cubeMapper->SetScalarRange(range);
+  std::cout << range[0] << " " << range[1] << std::endl;
   vtkActor *cubeActor = vtkActor::New();
   cubeActor->SetMapper(cubeMapper);
   vtkProperty * property = cubeActor->GetProperty();
index 96c160fd0b660daebec9327ecadb9aecba582bee..34fe5107f62053c80ee9c38ecd6af5566ef55b9b 100644 (file)
@@ -31,6 +31,15 @@ int main(int argc, char * argv[]) {
   typedef float PixelType;
   clitk::Xml2DicomRTStructFilter<PixelType> filter;
   filter.SetVerboseFlag(args_info.verbose_flag);
+  filter.SetUsePixel(args_info.pixel_flag);
+  if (args_info.pixel_flag) {
+    if (args_info.image_given)
+      filter.SetImageMHD(args_info.image_arg);
+    else {
+      std::cout << "Set MHD image (option -j) when pixel is set" << std::endl;
+      return -1;
+    }
+  }
   filter.SetInputFilename(args_info.input_arg);
   filter.SetDicomFolder(args_info.dicom_arg);
   filter.SetStructureSetFilename(args_info.rtstruct_arg);
index ab49cce6bd49560ef4b7727d614117d83296a75a..0b053833b4653699fe5ab82fb5d600dd19397456 100644 (file)
@@ -10,5 +10,7 @@ option "verbose"   v "Verbose"
 
 option "input"     i "Input xml file with struct to be converted into DicomRTStruct"  string  yes
 option "rtstruct"  r "Input rt struct"                                                string  yes
-option "dicom"     d "Input folder where the initial dicom ct is"                     string  yes
+option "dicom"     d "Input folder where the initial dicom used for the input is"     string  yes
 option "output"    o "Output DicomRT filename"                                        string  yes
+option "pixel"     p "Use pixel coordinates instead of mm coordinates (Usefull with transformation matrix)"  flag    off
+option "image"     j "If pixel is set, path to the mhd image converted from dicom"    string  no