]> Creatis software - clitk.git/commitdiff
Merge branch 'master' of /home/dsarrut/clitk3.server
authorRomulo Pinho <pinho@lyon.fnclcc.fr>
Thu, 20 Oct 2011 09:07:24 +0000 (11:07 +0200)
committerRomulo Pinho <pinho@lyon.fnclcc.fr>
Thu, 20 Oct 2011 09:07:24 +0000 (11:07 +0200)
19 files changed:
itk/clitkMeshToBinaryImageFilter.h
itk/clitkMeshToBinaryImageFilter.txx
itk/clitkSegmentationUtils.txx
segmentation/clitkBool.ggo
segmentation/clitkExtractLymphStation_1RL.txx
segmentation/clitkExtractLymphStation_2RL.txx
segmentation/clitkExtractLymphStation_3A.txx
segmentation/clitkExtractLymphStation_3P.txx
segmentation/clitkExtractLymphStation_Supports.txx
segmentation/clitkExtractLymphStationsFilter.txx
segmentation/clitkExtractMediastinalVessels.ggo
segmentation/clitkExtractMediastinalVesselsFilter.h
segmentation/clitkExtractMediastinalVesselsFilter.txx
segmentation/clitkExtractMediastinalVesselsGenericFilter.txx
segmentation/clitkExtractMediastinumFilter.txx
segmentation/clitkFilterWithAnatomicalFeatureDatabaseManagement.cxx
tools/CMakeLists.txt
tools/clitkMeshToBinaryImage.cxx [new file with mode: 0644]
tools/clitkMeshToBinaryImage.ggo [new file with mode: 0644]

index 67607b8c05846446ceaf193649054ed6eaaac471..acecc385558081696b5b5db039770c2edccbf313 100644 (file)
 // clitk
 #include "clitkCommon.h"
 
+// itk
+#include <itkImageSource.h>
+
+// vtk
+#include <vtkSmartPointer.h>
+#include <vtkPolyData.h>
+
 namespace clitk {
     
   /* --------------------------------------------------------------------     
@@ -67,6 +74,9 @@ namespace clitk {
     itkSetMacro(LikeImage, ImagePointer);
     itkGetConstMacro(LikeImage, ImagePointer);
 
+    itkSetMacro(Extrude, bool);
+    itkGetMacro(Extrude, bool);
+
   protected:
     MeshToBinaryImageFilter();
     virtual ~MeshToBinaryImageFilter() {}
@@ -74,6 +84,7 @@ namespace clitk {
     virtual void GenerateOutputInformation();
     virtual void GenerateData();
 
+    bool m_Extrude;
     ImagePointer m_LikeImage;
     vtkSmartPointer<vtkPolyData> m_Mesh;
 
index 77252584885783336cde934ed433a83d81d3677b..94424b8631ca4e9fa7d7e1f3c4af3fbf90bf1d70 100644 (file)
@@ -30,7 +30,7 @@
 //--------------------------------------------------------------------
 template <class ImageType>
 clitk::MeshToBinaryImageFilter<ImageType>::
-MeshToBinaryImageFilter():itk::ImageSource<ImageType>()
+MeshToBinaryImageFilter() : itk::ImageSource<ImageType>(), m_Extrude(true)
 {
 }
 //--------------------------------------------------------------------
@@ -102,11 +102,16 @@ GenerateData()
   sts->SetInformationInput(binary_image);
     
   // Extrusion
-  vtkSmartPointer<vtkLinearExtrusionFilter> extrude=vtkSmartPointer<vtkLinearExtrusionFilter>::New();
-  extrude->SetInput(m_Mesh);
-  // We extrude in the -slice_spacing direction to respect the FOCAL convention 
-  extrude->SetVector(0, 0, -m_LikeImage->GetSpacing()[2]);
-  sts->SetInput(extrude->GetOutput());
+  if (m_Extrude)
+  {
+    vtkSmartPointer<vtkLinearExtrusionFilter> extrude=vtkSmartPointer<vtkLinearExtrusionFilter>::New();
+    extrude->SetInput(m_Mesh);
+    // We extrude in the -slice_spacing direction to respect the FOCAL convention 
+    extrude->SetVector(0, 0, -m_LikeImage->GetSpacing()[2]);
+    sts->SetInput(extrude->GetOutput());
+  }
+  else
+    sts->SetInput(m_Mesh);
 
   // Stencil
   vtkSmartPointer<vtkImageStencil> stencil=vtkSmartPointer<vtkImageStencil>::New();
@@ -134,7 +139,7 @@ GenerateData()
   m_Exporter->SetInput( stencil->GetOutput() );
   m_Importer->Update();
 
-  writeImage<ImageType>(m_Importer->GetOutput(), "f.mhd");
+  // writeImage<ImageType>(m_Importer->GetOutput(), "f.mhd");
 
   this->SetNthOutput(0, m_Importer->GetOutput());
 
index adca988098ad8c78fee5cfb132b890fb0bef0c09..2f0d935a7cbbc88d038000768d6bed04bfd2451f 100644 (file)
@@ -464,10 +464,12 @@ namespace clitk {
     typename ImageType::RegionType region;
     typename ImageType::SizeType size = image->GetLargestPossibleRegion().GetSize();
     typename ImageType::PointType p = image->GetOrigin();
-    p[dim] = min;
+    if (min > p[dim]) p[dim] = min; // Check if not outside the image
     typename ImageType::IndexType start;
     image->TransformPhysicalPointToIndex(p, start);
-    p[dim] = max;
+    double m = image->GetOrigin()[dim] + size[dim]*image->GetSpacing()[dim];
+    if (max > m) p[dim] = m; // Check if not outside the image
+    else p[dim] = max;
     typename ImageType::IndexType end;
     image->TransformPhysicalPointToIndex(p, end);
     size[dim] = abs(end[dim]-start[dim]);
index ede9646053273eed8fa4fc020407b6b42beaf423..f8856553031e95e1d1ad54d71e5818f18caaedff 100644 (file)
@@ -11,4 +11,4 @@ option "input"         i      "Input mask 1"               string     yes
 option "input2"        j       "Input mask 2"               string     yes
 option "output"        o       "Output filename"            string     yes
 
-option "type"          t       "Bool type "               int default="1" no
+option "type"          t       "Bool type (0=AndNot, 1=And, 2=Or)"               int default="1" no
index 7fa72ef4b92e8d3bc2bf471f8b8e6b97f1ab3730..6758ab85220455679a2354e3550b2d95a929ac0d 100644 (file)
@@ -44,10 +44,10 @@ ExtractStation_1RL()
   ExtractStation_1RL_Post_Limits();
   m_Working_Support = m_ListOfStations["1R"];
   Remove_Structures(" 1R", "ScaleneMuscleAnt");
-  Remove_Structures(" 1R", "CommonCarotidArteryRight");  
+  Remove_Structures(" 1R", "RightCommonCarotidArtery");  
   m_Working_Support = m_ListOfStations["1L"];
   Remove_Structures(" 1L", "ScaleneMuscleAnt");
-  Remove_Structures(" 1L", "CommonCarotidArteryLeft");  
+  Remove_Structures(" 1L", "LeftCommonCarotidArtery");  
  
   // Generic RelativePosition processes
   m_ListOfStations["1R"] = this->ApplyRelativePositionList("Station_1R", m_ListOfStations["1R"]);
index e31f2d519c016ed31a4dba8497efd419c1c98b6b..7332dee29921c804c621c161e44127f17bbdc66c 100644 (file)
@@ -184,10 +184,10 @@ ExtractStation_2RL_Remove_Structures(std::string s)
   // m_Working_Support must be set
   Remove_Structures(s, "BrachioCephalicVein");
   Remove_Structures(s, "BrachioCephalicArtery");
-  Remove_Structures(s, "CommonCarotidArteryLeft");
-  Remove_Structures(s, "CommonCarotidArteryRight");
-  Remove_Structures(s, "SubclavianArteryLeft");
-  Remove_Structures(s, "SubclavianArteryRight");
+  Remove_Structures(s, "LeftCommonCarotidArtery");
+  Remove_Structures(s, "RightCommonCarotidArtery");
+  Remove_Structures(s, "LeftSubclavianArtery");
+  Remove_Structures(s, "RightSubclavianArtery");
   Remove_Structures(s, "Thyroid");
   Remove_Structures(s, "Aorta");
 }
index 26e478b81b763d190f318432da353ae01aef3c69..e973ed02d07655caac505bbbf2c78a24403320dd 100644 (file)
@@ -207,11 +207,11 @@ clitk::ExtractLymphStationsFilter<ImageType>::
 ExtractStation_3A_Remove_Structures() 
 {
   Remove_Structures(" 3A", "Aorta");
-  Remove_Structures(" 3A", "SubclavianArteryLeft");
-  Remove_Structures(" 3A", "SubclavianArteryRight");
+  Remove_Structures(" 3A", "LeftSubclavianArtery");
+  Remove_Structures(" 3A", "RightSubclavianArtery");
   Remove_Structures(" 3A", "Thyroid");
-  Remove_Structures(" 3A", "CommonCarotidArteryLeft");
-  Remove_Structures(" 3A", "CommonCarotidArteryRight");
+  Remove_Structures(" 3A", "LeftCommonCarotidArtery");
+  Remove_Structures(" 3A", "RightCommonCarotidArtery");
   Remove_Structures(" 3A", "BrachioCephalicArtery");
   //  Remove_Structures("3A", "Bones"); --> should be in extractmediastinum
   //  Remove_Structures("3A", "BrachioCephalicVein"); ?
index 712049809f464fc84c7f9db448c6669cfc59f03f..274fc867c740596227b8c2c92f84b31b4a5d5a02 100644 (file)
@@ -26,8 +26,6 @@ ExtractStation_3P()
   m_ListOfStations["3P"] = m_Working_Support;
   StopCurrentStep<MaskImageType>(m_Working_Support);
 
-  // ExtractStation_3P_LR_inf_Limits();  // <-- done in RelPosList
-  
   // Generic RelativePosition processes
   m_ListOfStations["3P"] = this->ApplyRelativePositionList("Station_3P", m_ListOfStations["3P"]);
   m_Working_Support = m_ListOfStations["3P"];
@@ -39,6 +37,7 @@ ExtractStation_3P()
   // LR limits superiorly => not here for the moment because not clear in the def
   // ExtractStation_3P_LR_sup_Limits_2(); //TODO
   // ExtractStation_3P_LR_sup_Limits();   // old version to change
+  // ExtractStation_3P_LR_inf_Limits();  // <-- done in RelPosList 
 
   // Store image filenames into AFDB 
   writeImage<MaskImageType>(m_ListOfStations["3P"], "seg/Station3P.mhd");
@@ -94,8 +93,8 @@ ExtractStation_3P_LR_sup_Limits()
   Trachea = clitk::ResizeImageLike<MaskImageType>(Trachea, m_Working_Support, GetBackgroundValue());
   SubclavianArtery = clitk::ResizeImageLike<MaskImageType>(SubclavianArtery, m_Working_Support, GetBackgroundValue());
   
-  writeImage<MaskImageType>(Trachea, "tr.mhd");
-  writeImage<MaskImageType>(SubclavianArtery, "sca.mhd");
+  // writeImage<MaskImageType>(Trachea, "tr.mhd");
+  // writeImage<MaskImageType>(SubclavianArtery, "sca.mhd");
   
   // Get list of slices
   std::vector<MaskSlicePointer> slices_support;
index 1fb82f86b5a6e555e12c8db81e8eff8a83eaa9c8..2cde97dd98b7465ae1addea1e6bcfc11965a4281 100644 (file)
@@ -21,10 +21,10 @@ ExtractStationSupports()
                                                      m_CarinaZ, true, GetBackgroundValue());
   m_ListOfSupports["Support_Superior_to_Carina"] = m_Support_Superior_to_Carina;
   m_ListOfSupports["Support_Inferior_to_Carina"] = m_Support_Inferior_to_Carina;
-  writeImage<MaskImageType>(m_Support_Inferior_to_Carina, "seg/Support_Inf_Carina.mhd");
-  this->GetAFDB()->SetImageFilename("Support_Inf_Carina", "seg/Support_Inf_Carina.mhd");
-  writeImage<MaskImageType>(m_Support_Superior_to_Carina, "seg/Support_Sup_Carina.mhd");
-  this->GetAFDB()->SetImageFilename("Support_Sup_Carina", "seg/Support_Sup_Carina.mhd");
+  writeImage<MaskImageType>(m_Support_Inferior_to_Carina, "seg/Support_Inf_Carina.mha");
+  this->GetAFDB()->SetImageFilename("Support_Inf_Carina", "seg/Support_Inf_Carina.mha");
+  writeImage<MaskImageType>(m_Support_Superior_to_Carina, "seg/Support_Sup_Carina.mha");
+  this->GetAFDB()->SetImageFilename("Support_Sup_Carina", "seg/Support_Sup_Carina.mha");
 
   // S1RL
   Support_SupInf_S1RL();
@@ -57,45 +57,45 @@ ExtractStationSupports()
   m_ListOfSupports["S11"] = clitk::Clone<MaskImageType>(m_Support_Inferior_to_Carina);
 
   // Store image filenames into AFDB 
-  writeImage<MaskImageType>(m_ListOfSupports["S1R"], "seg/Support_S1R.mhd");
-  this->GetAFDB()->SetImageFilename("Support_S1R", "seg/Support_S1R.mhd");
-  writeImage<MaskImageType>(m_ListOfSupports["S1L"], "seg/Support_S1L.mhd");
-  this->GetAFDB()->SetImageFilename("Support_S1L", "seg/Support_S1L.mhd");
+  writeImage<MaskImageType>(m_ListOfSupports["S1R"], "seg/Support_S1R.mha");
+  this->GetAFDB()->SetImageFilename("Support_S1R", "seg/Support_S1R.mha");
+  writeImage<MaskImageType>(m_ListOfSupports["S1L"], "seg/Support_S1L.mha");
+  this->GetAFDB()->SetImageFilename("Support_S1L", "seg/Support_S1L.mha");
 
-  writeImage<MaskImageType>(m_ListOfSupports["S2L"], "seg/Support_S2L.mhd");
-  this->GetAFDB()->SetImageFilename("Support_S2L", "seg/Support_S2L.mhd");
-  writeImage<MaskImageType>(m_ListOfSupports["S2R"], "seg/Support_S2R.mhd");
-  this->GetAFDB()->SetImageFilename("Support_S2R", "seg/Support_S2R.mhd");
+  writeImage<MaskImageType>(m_ListOfSupports["S2L"], "seg/Support_S2L.mha");
+  this->GetAFDB()->SetImageFilename("Support_S2L", "seg/Support_S2L.mha");
+  writeImage<MaskImageType>(m_ListOfSupports["S2R"], "seg/Support_S2R.mha");
+  this->GetAFDB()->SetImageFilename("Support_S2R", "seg/Support_S2R.mha");
 
-  writeImage<MaskImageType>(m_ListOfSupports["S3P"], "seg/Support_S3P.mhd");
-  this->GetAFDB()->SetImageFilename("Support_S3P", "seg/Support_S3P.mhd");
-  writeImage<MaskImageType>(m_ListOfSupports["S3A"], "seg/Support_S3A.mhd");
-  this->GetAFDB()->SetImageFilename("Support_S3A", "seg/Support_S3A.mhd");
+  writeImage<MaskImageType>(m_ListOfSupports["S3P"], "seg/Support_S3P.mha");
+  this->GetAFDB()->SetImageFilename("Support_S3P", "seg/Support_S3P.mha");
+  writeImage<MaskImageType>(m_ListOfSupports["S3A"], "seg/Support_S3A.mha");
+  this->GetAFDB()->SetImageFilename("Support_S3A", "seg/Support_S3A.mha");
 
-  writeImage<MaskImageType>(m_ListOfSupports["S4L"], "seg/Support_S4L.mhd");
-  this->GetAFDB()->SetImageFilename("Support_S4L", "seg/Support_S4L.mhd");
-  writeImage<MaskImageType>(m_ListOfSupports["S4R"], "seg/Support_S4R.mhd");
-  this->GetAFDB()->SetImageFilename("Support_S4R", "seg/Support_S4R.mhd");
+  writeImage<MaskImageType>(m_ListOfSupports["S4L"], "seg/Support_S4L.mha");
+  this->GetAFDB()->SetImageFilename("Support_S4L", "seg/Support_S4L.mha");
+  writeImage<MaskImageType>(m_ListOfSupports["S4R"], "seg/Support_S4R.mha");
+  this->GetAFDB()->SetImageFilename("Support_S4R", "seg/Support_S4R.mha");
 
-  writeImage<MaskImageType>(m_ListOfSupports["S5"], "seg/Support_S5.mhd");
-  this->GetAFDB()->SetImageFilename("Support_S5", "seg/Support_S5.mhd");
-  writeImage<MaskImageType>(m_ListOfSupports["S6"], "seg/Support_S6.mhd");
-  this->GetAFDB()->SetImageFilename("Support_S6", "seg/Support_S6.mhd");
+  writeImage<MaskImageType>(m_ListOfSupports["S5"], "seg/Support_S5.mha");
+  this->GetAFDB()->SetImageFilename("Support_S5", "seg/Support_S5.mha");
+  writeImage<MaskImageType>(m_ListOfSupports["S6"], "seg/Support_S6.mha");
+  this->GetAFDB()->SetImageFilename("Support_S6", "seg/Support_S6.mha");
 
-  writeImage<MaskImageType>(m_ListOfSupports["S7"], "seg/Support_S7.mhd");
-  this->GetAFDB()->SetImageFilename("Support_S7", "seg/Support_S7.mhd");
+  writeImage<MaskImageType>(m_ListOfSupports["S7"], "seg/Support_S7.mha");
+  this->GetAFDB()->SetImageFilename("Support_S7", "seg/Support_S7.mha");
 
-  writeImage<MaskImageType>(m_ListOfSupports["S8"], "seg/Support_S8.mhd");
-  this->GetAFDB()->SetImageFilename("Support_S8", "seg/Support_S8.mhd");
+  writeImage<MaskImageType>(m_ListOfSupports["S8"], "seg/Support_S8.mha");
+  this->GetAFDB()->SetImageFilename("Support_S8", "seg/Support_S8.mha");
 
-  writeImage<MaskImageType>(m_ListOfSupports["S9"], "seg/Support_S9.mhd");
-  this->GetAFDB()->SetImageFilename("Support_S9", "seg/Support_S9.mhd");
+  writeImage<MaskImageType>(m_ListOfSupports["S9"], "seg/Support_S9.mha");
+  this->GetAFDB()->SetImageFilename("Support_S9", "seg/Support_S9.mha");
 
-  writeImage<MaskImageType>(m_ListOfSupports["S10"], "seg/Support_S10.mhd");
-  this->GetAFDB()->SetImageFilename("Support_S10", "seg/Support_S10.mhd");  
+  writeImage<MaskImageType>(m_ListOfSupports["S10"], "seg/Support_S10.mha");
+  this->GetAFDB()->SetImageFilename("Support_S10", "seg/Support_S10.mha");  
 
-  writeImage<MaskImageType>(m_ListOfSupports["S11"], "seg/Support_S11.mhd");
-  this->GetAFDB()->SetImageFilename("Support_S11", "seg/Support_S11.mhd");  
+  writeImage<MaskImageType>(m_ListOfSupports["S11"], "seg/Support_S11.mha");
+  this->GetAFDB()->SetImageFilename("Support_S11", "seg/Support_S11.mha");  
   WriteAFDB();
 }
 //--------------------------------------------------------------------
@@ -544,7 +544,7 @@ void
 clitk::ExtractLymphStationsFilter<ImageType>::
 Support_S5()
 {
-  StartNewStep("[Support] Sup-Inf limits S5 with aorta");
+  StartNewStep("[Support] Sup-Inf limits S5 with Aorta and MainPulmonaryArtery");
 
   // Initial S5 support
   MaskImagePointer S5 = 
@@ -555,11 +555,11 @@ Support_S5()
   
   // Inf limits with "upper rim of the left main pulmonary artery"
   // For the moment only, it will change.
-  MaskImagePointer PulmonaryTrunk = this->GetAFDB()->template GetImage<MaskImageType>("PulmonaryTrunk");
+  MaskImagePointer MainPulmonaryArtery = this->GetAFDB()->template GetImage<MaskImageType>("MainPulmonaryArtery");
   MaskImagePointType p;
   p[0] = p[1] = p[2] =  0.0; // to avoid warning
-  clitk::FindExtremaPointInAGivenDirection<MaskImageType>(PulmonaryTrunk, GetBackgroundValue(), 2, false, p);
-  p[2] += PulmonaryTrunk->GetSpacing()[2];
+  clitk::FindExtremaPointInAGivenDirection<MaskImageType>(MainPulmonaryArtery, GetBackgroundValue(), 2, false, p);
+  p[2] += MainPulmonaryArtery->GetSpacing()[2];
   
   // Cut Sup/Inf
   S5 = clitk::CropImageAlongOneAxis<MaskImageType>(S5, 2, p[2], sup, true, GetBackgroundValue());
index 9051ec154df9cf5c8a8cc94cb640556a27f73819..ec93c098b3421ea3d8a9105f5677c1212771ed45 100644 (file)
@@ -126,10 +126,10 @@ GenerateOutputInformation() {
   }
 
   // Extract Stations
+  ExtractStation_1RL();
+  ExtractStation_2RL();
   ExtractStation_3P();
   ExtractStation_3A();
-  ExtractStation_2RL();
-  ExtractStation_1RL();
   ExtractStation_4RL();
   ExtractStation_5();
   ExtractStation_6();
@@ -137,13 +137,13 @@ GenerateOutputInformation() {
   // ---------- TODO -----------------------
 
   // Extract Station8
-  ExtractStation_8();
+  //  ExtractStation_8();
 
   // Extract Station7
-  this->StartNewStep("Station 7");
-  this->StartSubStep();
-  ExtractStation_7();
-  this->StopSubStep();
+  //this->StartNewStep("Station 7");
+  //this->StartSubStep();
+  //ExtractStation_7();
+  //this->StopSubStep();
   
 }
 //--------------------------------------------------------------------
@@ -988,8 +988,8 @@ FindAntPostVesselsOLD()
   binarizedContour = clitk::CropImageAlongOneAxis<MaskImageType>(binarizedContour, 2, inf, sup, 
                                                         false, GetBackgroundValue());
   // Update the AFDB
-  writeImage<MaskImageType>(binarizedContour, "seg/AntPostVesselsSeparation.mhd");
-  this->GetAFDB()->SetImageFilename("AntPostVesselsSeparation", "seg/AntPostVesselsSeparation.mhd");  
+  writeImage<MaskImageType>(binarizedContour, "seg/AntPostVesselsSeparation.mha");
+  this->GetAFDB()->SetImageFilename("AntPostVesselsSeparation", "seg/AntPostVesselsSeparation.mha");
   this->WriteAFDB();
   return binarizedContour;
 
@@ -1076,10 +1076,10 @@ FindAntPostVessels2()
   typedef std::map<std::string, MaskImagePointer>::iterator MapIter;
   MapOfStructures["BrachioCephalicArtery"] = this->GetAFDB()->template GetImage<MaskImageType>("BrachioCephalicArtery");
   MapOfStructures["BrachioCephalicVein"] = this->GetAFDB()->template GetImage<MaskImageType>("BrachioCephalicVein");
-  MapOfStructures["CommonCarotidArteryLeft"] = this->GetAFDB()->template GetImage<MaskImageType>("CommonCarotidArteryLeft");
-  MapOfStructures["CommonCarotidArteryRight"] = this->GetAFDB()->template GetImage<MaskImageType>("CommonCarotidArteryRight");
-  MapOfStructures["SubclavianArteryLeft"] = this->GetAFDB()->template GetImage<MaskImageType>("SubclavianArteryLeft");
-  MapOfStructures["SubclavianArteryRight"] = this->GetAFDB()->template GetImage<MaskImageType>("SubclavianArteryRight");
+  MapOfStructures["CommonCarotidArteryLeft"] = this->GetAFDB()->template GetImage<MaskImageType>("LeftCommonCarotidArtery");
+  MapOfStructures["CommonCarotidArteryRight"] = this->GetAFDB()->template GetImage<MaskImageType>("RightCommonCarotidArtery");
+  MapOfStructures["SubclavianArteryLeft"] = this->GetAFDB()->template GetImage<MaskImageType>("LeftSubclavianArtery");
+  MapOfStructures["SubclavianArteryRight"] = this->GetAFDB()->template GetImage<MaskImageType>("RightSubclavianArtery");
   MapOfStructures["Thyroid"] = this->GetAFDB()->template GetImage<MaskImageType>("Thyroid");
   MapOfStructures["Aorta"] = this->GetAFDB()->template GetImage<MaskImageType>("Aorta");
   MapOfStructures["Trachea"] = this->GetAFDB()->template GetImage<MaskImageType>("Trachea");
@@ -1347,8 +1347,8 @@ FindAntPostVessels2()
                                                         false, GetBackgroundValue());
 
   // Update the AFDB
-  writeImage<MaskImageType>(binarizedContour, "seg/AntPostVesselsSeparation.mhd");
-  this->GetAFDB()->SetImageFilename("AntPostVesselsSeparation", "seg/AntPostVesselsSeparation.mhd");  
+  writeImage<MaskImageType>(binarizedContour, "seg/AntPostVesselsSeparation.mha");
+  this->GetAFDB()->SetImageFilename("AntPostVesselsSeparation", "seg/AntPostVesselsSeparation.mha");
   this->WriteAFDB();
   return binarizedContour;
 
index 77bd309988d39eac91c130421ccc2f463e5d0821..2e537971f0840a2a080e44737bbb9e64fdae5eef 100644 (file)
@@ -10,7 +10,8 @@ option "verboseStep"    -  "Verbose each step"                  flag          off
 option "writeStep"      w  "Write image at each step"    flag          off
 option "verboseOption"  -  "Display options values"       flag          off
 option "verboseWarningOff" -  "Do not display warning"    flag          off
-option "verboseMemory"  -  "Display memory usage"         flag          off
+option "verboseMemory"     -  "Display memory usage"      flag          off
+option "verboseTracking"   -  "Verbose tracking"          flag          off
 
 section "I/O"
 
index 44bb53af2fe267e2d9a6b3db6a62efe81b477616..0f3df6cbf09c41bb942b5d87c622fd2ea90d507a 100644 (file)
@@ -113,6 +113,10 @@ namespace clitk {
     itkGetConstMacro(DebugFlag, bool);
     itkBooleanMacro(DebugFlag);
 
+    itkSetMacro(VerboseTrackingFlag, bool);
+    itkGetConstMacro(VerboseTrackingFlag, bool);
+    itkBooleanMacro(VerboseTrackingFlag);
+
     itkSetMacro(SoughtVesselSeedName, std::string);
     itkGetConstMacro(SoughtVesselSeedName, std::string);
 
@@ -137,6 +141,7 @@ namespace clitk {
     virtual void GenerateData();
     
     bool               m_DebugFlag;
+    bool               m_VerboseTrackingFlag;
     ImagePointer       m_Input;
     MaskImagePointer   m_Working_Support;
     MaskImagePointer   m_Mediastinum;
index 29709f48df1d8fc6bb95f85b8f7af80c02afb8ce..537b14867b1d3cfa74692b8e4326f4d65873f04f 100644 (file)
@@ -58,6 +58,7 @@ ExtractMediastinalVesselsFilter():
   SetMaxDistanceRightToCarina(35);
   SetSoughtVesselSeedName("NoSeedNameGiven");
   SetFinalOpeningRadius(1);
+  VerboseTrackingFlagOff();
 }
 //--------------------------------------------------------------------
 
@@ -213,7 +214,7 @@ GenerateOutputInformation() {
 
   // Build the final 3D image from the previous slice by slice processing
   m_SoughtVessel = clitk::JoinSlices<MaskImageType>(m_slice_recon, m_Mask, 2);
-  writeImage<MaskImageType>(m_SoughtVessel, "recon2.mhd");
+  // writeImage<MaskImageType>(m_SoughtVessel, "recon2.mhd");
   
   // Set binary image, (remove other labels).  
   m_SoughtVessel = 
@@ -294,6 +295,10 @@ CropInputImage() {
   m_Input = 
     clitk::CropImageRemoveLowerThan<ImageType>(m_Input, 2, m_CarinaZ, false, GetBackgroundValue());  
 
+  // // Get seed, crop around
+  // MaskImagePointType SoughtVesselSeedPoint;
+  // GetAFDB()->GetPoint3D(m_SoughtVesselSeedName, SoughtVesselSeedPoint);
+
   // Crop post
   double m_CarinaY = centroids[1][1];
   m_Input = clitk::CropImageRemoveGreaterThan<ImageType>(m_Input, 1,
@@ -359,9 +364,20 @@ TrackBifurcationFromPoint(MaskImagePointer & recon,
   uint currentSlice=index[2]; 
   bool found = false;
   LabelType previous_slice_label=recon->GetPixel(index);
+
+
+  if (GetVerboseTrackingFlag()) {
+    std::cout << "TrackBifurcationFromPoint " << std::endl;
+    std::cout << "\t point3D = " << point3D << std::endl;
+    std::cout << "\t pointMaxSlice = " << pointMaxSlice << std::endl;
+    std::cout << "\t newLabel = " << newLabel << std::endl;
+  }
+
   //  DD(slices_recon.size());
   do {
-    //    DD(currentSlice);
+    if (GetVerboseTrackingFlag()) {
+      std::cout << "currentSlice = " << currentSlice << std::endl;
+    }
     // Consider current reconstructed slice
     MaskSlicePointer s = slices_recon[currentSlice];
     MaskSlicePointer previous;
@@ -413,7 +429,10 @@ TrackBifurcationFromPoint(MaskImagePointer & recon,
     // -------------------------
     // If no centroid were found
     if (centroids.size() == 0) {
-      // Last attempt to find -> check if previous centroid is inside a CCL
+      if (GetVerboseTrackingFlag()) {
+        std::cout << "no centroid" << std::endl;
+      }
+       // Last attempt to find -> check if previous centroid is inside a CCL
       //      if in s -> get value, getcentroid add.
       //      DD(currentSlice);
       //DD("Last change to find");
@@ -434,6 +453,9 @@ TrackBifurcationFromPoint(MaskImagePointer & recon,
     // Some centroid were found
     // If several centroids, we found a bifurcation
     if (centroids.size() > 1) {
+      if (GetVerboseTrackingFlag()) {
+        std::cout << "Centroid" << centroids.size() << std::endl;
+      }
       //      int n = centroids.size();
       // Debug point
       std::vector<ImagePointType> d;
@@ -498,6 +520,10 @@ TrackBifurcationFromPoint(MaskImagePointer & recon,
     
     // if only one centroids, we change the current image with the current label 
     if (centroids.size() == 1) {
+      if (GetVerboseTrackingFlag()) {
+        std::cout << "Centroid" << centroids.size() << std::endl;
+      }
+      
       //DD(centroids_label[0]);
       s = clitk::SetBackground<MaskSliceType, MaskSliceType>(s, s, centroids_label[0], newLabel, true);
       slices_recon[currentSlice] = s;
@@ -534,6 +560,10 @@ TrackBifurcationFromPoint(MaskImagePointer & recon,
       previousCenter = centroids[1];
     }
 
+    if (GetVerboseTrackingFlag()) {
+        std::cout << "End iteration c=" << centroids.size() << std::endl;
+    }
+    
     if (centroids.size() == 0) {
       //      DD("ZERO");
       found = true;
@@ -598,14 +628,25 @@ TrackVesselsFromPoint(MaskImagePointer & recon,
   std::vector<MaskSlicePointer> shapeObjectsSliceList; // keep slice, useful for 'union'
   typename LabelObjectType::Pointer previousShapeObject;
 
+ if (GetVerboseTrackingFlag()) {
+    std::cout << "TrackBifurcationFromPoint " << std::endl;
+    std::cout << "\t seedPoint = " << seedPoint << std::endl;
+    std::cout << "\t pointMaxSlice = " << pointMaxSlice << std::endl;
+    std::cout << "\t newLabel = " << newLabel << std::endl;
+  }
+
   do {
     // Debug
-    std::cout << std::endl;
+    //std::cout << std::endl;
     //DD(currentSlice);
     ImagePointType c;
     clitk::PointsUtils<MaskImageType>::Convert2DTo3D(previousCentroid, m_Mask, currentSlice-1, c);
     //DD(c);
 
+    if (GetVerboseTrackingFlag()) {
+      std::cout << "Loop slice = " << currentSlice << " c = " << c << std::endl;
+    }
+
     // Consider current reconstructed slice
     MaskSlicePointer s = slices[currentSlice];
     MaskSlicePointer previous;
index 71f8c6b7ec0713e6717dc30a3dbad219880c3283..4a681a6d2fad3a599efe58359221119fe5c74244 100644 (file)
@@ -83,6 +83,8 @@ SetOptionsFromArgsInfoToFilter(FilterType * f)
   f->SetMaxNumberOfFoundBifurcation(mArgsInfo.bif_arg);
 
   f->SetFinalOpeningRadius(mArgsInfo.open_arg);
+
+  f->SetVerboseTrackingFlag(mArgsInfo.verboseTracking_flag);
   
   // Output filename
   this->AddOutputFilename(mArgsInfo.output_arg);
index 312fd14ed01d518cddfed3236b5a0ff720526f00..39e5542b5df9fdaf670c5e66bc354f103532ab0f 100644 (file)
@@ -54,7 +54,7 @@ ExtractMediastinumFilter():
   SetForegroundValueLeftLung(1);
   SetForegroundValueRightLung(2);
   SetDistanceMaxToAnteriorPartOfTheVertebralBody(10);  
-  SetOutputMediastinumFilename("mediastinum.mhd");  
+  SetOutputMediastinumFilename("mediastinum.mha");  
   UseBonesOff();
 }
 //--------------------------------------------------------------------
@@ -250,9 +250,9 @@ GenerateOutputInformation() {
   left_lung = clitk::SetBackground<MaskImageType, MaskImageType>(left_lung, left_lung, 2, 1, false);
   right_lung = clitk::ResizeImageLike<MaskImageType>(right_lung, output, this->GetBackgroundValue());
   left_lung = clitk::ResizeImageLike<MaskImageType>(left_lung, output, this->GetBackgroundValue());
-  this->GetAFDB()->template SetImage<MaskImageType>("RightLung", "seg/RightLung.mhd", 
+  this->GetAFDB()->template SetImage<MaskImageType>("RightLung", "seg/RightLung.mha", 
                                                     right_lung, true);
-  this->GetAFDB()->template SetImage<MaskImageType>("LeftLung", "seg/LeftLung.mhd", 
+  this->GetAFDB()->template SetImage<MaskImageType>("LeftLung", "seg/LeftLung.mha", 
                                                     left_lung, true);
   this->GetAFDB()->Write();
   this->template StopCurrentStep<MaskImageType>(output);
@@ -319,12 +319,12 @@ GenerateOutputInformation() {
     bones_post = roiFilter->GetOutput();
     // writeImage<MaskImageType>(bones_post, "b_post.mhd");
 
-     // Now, insert this image in the AFDB
-     this->GetAFDB()->template SetImage<MaskImageType>("Bones_Post", "seg/Bones_Post.mhd", 
+    // Now, insert this image in the AFDB ==> (needed because used in the RelativePosition config file)
+    this->GetAFDB()->template SetImage<MaskImageType>("Bones_Post", "seg/Bones_Post.mha", 
                                                  bones_post, true);
-     this->GetAFDB()->template SetImage<MaskImageType>("Bones_Ant", "seg/Bones_Ant.mhd", 
+    this->GetAFDB()->template SetImage<MaskImageType>("Bones_Ant", "seg/Bones_Ant.mha", 
                                                  bones_ant, true);
-     this->GetAFDB()->Write();
+    this->GetAFDB()->Write();
 
     this->template StopCurrentStep<MaskImageType>(output);
   }
@@ -345,11 +345,16 @@ GenerateOutputInformation() {
   // RelativePosition to avoid some issue due to superior boundaries.
   this->StartNewStep("[Mediastinum] Keep inferior to CricoidCartilag");
   // load Cricoid, get centroid, cut above (or below), lower bound
-  MaskImagePointer CricoidCartilag = this->GetAFDB()->template GetImage <MaskImageType>("CricoidCartilag");
   MaskImagePointType p;
-  p[0] = p[1] = p[2] =  0.0; // to avoid warning
-  clitk::FindExtremaPointInAGivenDirection<MaskImageType>(CricoidCartilag, 
-                                                          this->GetBackgroundValue(), 2, true, p);
+  try {
+    MaskImagePointer CricoidCartilag = this->GetAFDB()->template GetImage <MaskImageType>("CricoidCartilag");
+    p[0] = p[1] = p[2] =  0.0; // to avoid warning
+    clitk::FindExtremaPointInAGivenDirection<MaskImageType>(CricoidCartilag, 
+                                                            this->GetBackgroundValue(), 2, true, p);
+  } catch (clitk::ExceptionObject e) {
+    //DD("CricoidCartilag image not found, try CricoidCartilagZ");
+    this->GetAFDB()->GetPoint3D("CricoidCartilagPoint", p);
+  }
   output = clitk::CropImageRemoveGreaterThan<MaskImageType>(output, 2, p[2], true, this->GetBackgroundValue());
   this->template StopCurrentStep<MaskImageType>(output);
 
index 6e3e5946c7228f0f3b5be1fc73c43d5e47aa6893..bec6b9e59a8a24d380b2b4f45c086e2765ae64b5 100644 (file)
@@ -45,7 +45,7 @@ void clitk::FilterWithAnatomicalFeatureDatabaseManagement::LoadAFDB()
   try {
     GetAFDB()->Load();
   } catch (clitk::ExceptionObject e) {
-    std::cout << "Could not read '" << GetAFDBFilename() << "', create one AFDB." << std::endl;
+    std::cout << "******* Could not read '" << GetAFDBFilename() << "', create one AFDB. ********" << std::endl;
     GetAFDB();
   }
 }
index 8bf0cbc26cf5ebf2a0a443040fac459a03bc1f1e..5222f12b9650d8d0ba9e42fca6b5f3530d942d34 100644 (file)
@@ -346,6 +346,11 @@ IF (CLITK_BUILD_TOOLS)
     TARGET_LINK_LIBRARIES(clitkBinaryImageToMesh ${ITK_LIBRARIES} vtkCommon vtkRendering)
     SET(TOOLS_INSTALL ${TOOLS_INSTALL} clitkBinaryImageToMesh)
 
+    WRAP_GGO(clitkMeshToBinaryImage_GGO_C clitkMeshToBinaryImage.ggo)
+    ADD_EXECUTABLE(clitkMeshToBinaryImage clitkMeshToBinaryImage.cxx ${clitkMeshToBinaryImage_GGO_C})
+    TARGET_LINK_LIBRARIES(clitkMeshToBinaryImage clitkCommon ${ITK_LIBRARIES} ${VTK_LIBRARIES})
+    SET(TOOLS_INSTALL ${TOOLS_INSTALL} clitkMeshToBinaryImage)
+
     WRAP_GGO(clitkMeshViewer_GGO_C clitkMeshViewer.ggo)
     ADD_EXECUTABLE(clitkMeshViewer clitkMeshViewer.cxx ${clitkMeshViewer_GGO_C})
     TARGET_LINK_LIBRARIES(clitkMeshViewer vtkCommon vtkRendering)
diff --git a/tools/clitkMeshToBinaryImage.cxx b/tools/clitkMeshToBinaryImage.cxx
new file mode 100644 (file)
index 0000000..2c3d8cc
--- /dev/null
@@ -0,0 +1,94 @@
+/*=========================================================================
+  Program:   vv                     http://www.creatis.insa-lyon.fr/rio/vv
+
+  Authors belong to:
+  - University of LYON              http://www.universite-lyon.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
+  the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
+  PURPOSE.  See the copyright notices for more information.
+
+  It is distributed under dual licence
+
+  - BSD        See included LICENSE.txt file
+  - CeCILL-B   http://www.cecill.info/licences/Licence_CeCILL-B_V1-en.html
+===========================================================================*/
+#include "clitkMeshToBinaryImage_ggo.h"
+#include "clitkMeshToBinaryImageFilter.h"
+#include "clitkCommon.h"
+#include "clitkImageCommon.h"
+
+#include <itkImageIOBase.h>
+#include <itkInvertIntensityImageFilter.h>
+
+#include <vtkOBJReader.h>
+
+#include <string>
+
+template <unsigned dim>
+void run(const args_info_clitkMeshToBinaryImage& argsInfo,
+         const itk::ImageIOBase * header)
+{
+  typedef itk::Image<unsigned char, dim> ImageType;
+  typename ImageType::Pointer like = ImageType::New();
+
+  typename ImageType::SizeType size;
+  for (unsigned i = 0; i < dim; ++i)
+    size[i] = header->GetDimensions(i);
+  typename ImageType::RegionType region;
+  region.SetSize(size);
+  like->SetRegions(region);
+
+  typename ImageType::SpacingType spacing;
+  for (unsigned i = 0; i < dim; ++i)
+    spacing[i] = header->GetSpacing(i);
+  like->SetSpacing(spacing);
+
+
+  typename ImageType::PointType origin;
+  for (unsigned i = 0; i < dim; ++i)
+    origin[i] = header->GetOrigin(i);
+  like->SetOrigin(origin);
+
+  vtkSmartPointer<vtkOBJReader> reader = vtkOBJReader::New();
+  reader->SetFileName(argsInfo.input_arg);
+  reader->Update();
+
+  typename clitk::MeshToBinaryImageFilter<ImageType>::Pointer filter =
+    clitk::MeshToBinaryImageFilter<ImageType>::New();
+  filter->SetExtrude(false);
+  filter->SetMesh(reader->GetOutput());
+  filter->SetLikeImage(like);
+  filter->Update();
+
+  typedef itk::InvertIntensityImageFilter<ImageType> InvertFilterType;
+  typename InvertFilterType::Pointer ifilter = InvertFilterType::New();
+  ifilter->SetInput(filter->GetOutput());
+  ifilter->SetMaximum(1);
+
+  clitk::writeImage(ifilter->GetOutput(), argsInfo.output_arg);
+}
+
+int main(int argc, char** argv)
+{
+  GGO(clitkMeshToBinaryImage, args_info);
+
+  itk::ImageIOBase::Pointer header =
+    clitk::readImageHeader(args_info.like_arg);
+  switch (header->GetNumberOfDimensions())
+  {
+    case 2:
+      run<2>(args_info, header);
+      break;
+    case 3:
+      run<3>(args_info, header);
+      break;
+    case 4:
+      run<4>(args_info, header);
+      break;
+  }
+
+  return EXIT_SUCCESS;
+}
diff --git a/tools/clitkMeshToBinaryImage.ggo b/tools/clitkMeshToBinaryImage.ggo
new file mode 100644 (file)
index 0000000..64577b3
--- /dev/null
@@ -0,0 +1,11 @@
+#File clitkBinaryImageToMesh.ggo
+package "clitkMeshToBinaryImage"
+version "1.0"
+purpose "Transform a vtk obj file into a binary mask"
+
+option "config"   - "Config file"   string  no
+option "verbose"  v     "Verbose"   flag  off
+
+option "input"      i "Input obj file"   string  yes
+option "output"     o "Output binary image"  string yes 
+option "like"       l "Output image will be like this one (size, spacing)"  string yes