]> Creatis software - clitk.git/blob - segmentation/clitkExtractLymphStation_1RL.txx
motion masks with and without bands
[clitk.git] / segmentation / clitkExtractLymphStation_1RL.txx
1
2
3 // vtk
4 #include <vtkAppendPolyData.h>
5 #include <vtkPolyDataWriter.h>
6 #include <vtkCellArray.h>
7
8 // clitk
9 #include "clitkMeshToBinaryImageFilter.h"
10
11 // itk
12 #include <itkImageDuplicator.h>
13
14 //--------------------------------------------------------------------
15 template <class ImageType>
16 void 
17 clitk::ExtractLymphStationsFilter<ImageType>::
18 ExtractStation_1RL_SetDefaultValues()
19 {
20   
21 }
22 //--------------------------------------------------------------------
23
24
25 //--------------------------------------------------------------------
26 template <class TImageType>
27 void 
28 clitk::ExtractLymphStationsFilter<TImageType>::
29 ExtractStation_1RL()
30 {
31   if ((!CheckForStation("1R")) && (!CheckForStation("1L"))) return;
32
33   StartNewStep("Stations 1RL");
34   StartSubStep(); 
35
36   // Get the current support 
37   StartNewStep("[Station 1RL] Get the current 1RL suppport");
38   m_ListOfStations["1R"] = m_ListOfSupports["S1R"];
39   m_ListOfStations["1L"] = m_ListOfSupports["S1L"];
40   StopCurrentStep<MaskImageType>(m_ListOfStations["1R"]);
41     
42   // Specific processes
43   ExtractStation_1RL_Ant_Limits();
44   ExtractStation_1RL_Post_Limits();
45   m_Working_Support = m_ListOfStations["1R"];
46   Remove_Structures(" 1R", "ScaleneMuscleAnt");
47   Remove_Structures(" 1R", "RightCommonCarotidArtery");  
48   m_Working_Support = m_ListOfStations["1L"];
49   Remove_Structures(" 1L", "ScaleneMuscleAnt");
50   Remove_Structures(" 1L", "LeftCommonCarotidArtery");  
51  
52   // Generic RelativePosition processes
53   m_ListOfStations["1R"] = this->ApplyRelativePositionList("Station_1R", m_ListOfStations["1R"]);
54   m_ListOfStations["1L"] = this->ApplyRelativePositionList("Station_1L", m_ListOfStations["1L"]);
55
56   // Store image filenames into AFDB 
57   writeImage<MaskImageType>(m_ListOfStations["1R"], "seg/Station1R.mhd");
58   writeImage<MaskImageType>(m_ListOfStations["1L"], "seg/Station1L.mhd");
59   GetAFDB()->SetImageFilename("Station1R", "seg/Station1R.mhd"); 
60   GetAFDB()->SetImageFilename("Station1L", "seg/Station1L.mhd"); 
61   WriteAFDB(); 
62   StopSubStep();
63 }
64 //--------------------------------------------------------------------
65
66
67 //--------------------------------------------------------------------
68 template <class ImageType>
69 void 
70 clitk::ExtractLymphStationsFilter<ImageType>::
71 ExtractStation_1RL_Ant_Limits()
72 {
73   // -----------------------------------------------------
74   StartNewStep("[Station 1RL] anterior limits with Trachea and Thyroid");
75   
76   /*
77     The idea here it to consider the most anterior points int the
78     Trachea or the Thyroid and cut all parts anterior to it. This is
79     an heuristic, not written explicitely in the articles.
80   */
81
82   MaskImagePointer Trachea = GetAFDB()->template GetImage<MaskImageType>("Trachea");
83   MaskImagePointer Thyroid = GetAFDB()->template GetImage<MaskImageType>("Thyroid");
84   MaskImagePointer S1R = m_ListOfStations["1R"];
85   MaskImagePointer S1L = m_ListOfStations["1L"];
86   MaskImagePointer support = m_ListOfSupports["S1RL"];
87  
88   // Resize like S1R. Warning: for Thyroid, only Right part is thus
89   // taken into account
90   Trachea = clitk::ResizeImageLike<MaskImageType>(Trachea, S1R, GetBackgroundValue());
91   Thyroid = clitk::ResizeImageLike<MaskImageType>(Thyroid, S1R, GetBackgroundValue());
92
93   // Search for most Ant point, slice by slice, between Trachea and Thyroid
94   std::vector<MaskSlicePointer> Trachea_slices;
95   clitk::ExtractSlices<MaskImageType>(Trachea, 2, Trachea_slices);
96   std::vector<MaskSlicePointer> Thyroid_slices;
97   clitk::ExtractSlices<MaskImageType>(Thyroid, 2, Thyroid_slices);
98   std::vector<typename ImageType::PointType> A;
99   std::vector<typename ImageType::PointType> B;
100   for(uint i=0; i<Trachea_slices.size(); i++) {
101     MaskSlicePointType p;
102     MaskSlicePointType q;
103     FindExtremaPointInAGivenDirection<MaskSliceType>(Trachea_slices[i], 
104                                                      GetBackgroundValue(), 
105                                                      1, true, p);
106     FindExtremaPointInAGivenDirection<MaskSliceType>(Thyroid_slices[i], 
107                                                      GetBackgroundValue(), 
108                                                      1, true, q);
109     if (q[1] < p[1]) p = q; // Now p is the most ant.
110     // Add a little margin, 3mm
111     p[1] -= 3;
112     // Convert in 3D
113     ImagePointType x; //dummy
114     A.push_back(x);
115     B.push_back(x);
116     clitk::PointsUtils<MaskImageType>::Convert2DTo3D(p, Trachea, i, A[i]);
117     p[0] += 10;    
118     clitk::PointsUtils<MaskImageType>::Convert2DTo3D(p, Trachea, i, B[i]);
119   }  
120
121   // Remove anterior to this line (keep +10 offset)
122   clitk::SliceBySliceSetBackgroundFromLineSeparation<MaskImageType>(S1R, A, B, 
123                                                                     GetBackgroundValue(), 1, +10);
124   m_ListOfStations["1R"] = S1R;
125
126   clitk::SliceBySliceSetBackgroundFromLineSeparation<MaskImageType>(S1L, A, B, 
127                                                                     GetBackgroundValue(), 1, +10);
128   StopCurrentStep<MaskImageType>(S1L);
129   m_ListOfStations["1L"] = S1L;
130 }
131 //--------------------------------------------------------------------
132
133
134 //--------------------------------------------------------------------
135 template <class ImageType>
136 void 
137 clitk::ExtractLymphStationsFilter<ImageType>::
138 ExtractStation_1RL_Post_Limits()
139 {
140   // -----------------------------------------------------
141   StartNewStep("[Station 1RL] Posterior limits with VertebralArtery");
142   
143   MaskImagePointer VertebralArtery = GetAFDB()->template GetImage<MaskImageType>("VertebralArtery");
144   MaskImagePointer S1R = m_ListOfStations["1R"];
145   MaskImagePointer S1L = m_ListOfStations["1L"];
146
147   // Resize like S1R.
148   VertebralArtery = clitk::ResizeImageLike<MaskImageType>(VertebralArtery, S1R, GetBackgroundValue());
149
150   // Search for most Ant point
151   std::vector<MaskSlicePointer> VertebralArtery_slices;
152   clitk::ExtractSlices<MaskImageType>(VertebralArtery, 2, VertebralArtery_slices);
153   std::vector<typename ImageType::PointType> A;
154   std::vector<typename ImageType::PointType> B;
155   for(uint i=0; i<VertebralArtery_slices.size(); i++) {
156     MaskSlicePointType p;
157     FindExtremaPointInAGivenDirection<MaskSliceType>(VertebralArtery_slices[i], 
158                                                      GetBackgroundValue(), 
159                                                      1, false, p);
160     // Add a little margin ? No.
161     //p[1] += 0;
162     //DD(p);
163     // Convert in 3D
164     ImagePointType x; //dummy
165     A.push_back(x);
166     B.push_back(x);
167     clitk::PointsUtils<MaskImageType>::Convert2DTo3D(p, VertebralArtery, i, A[i]);
168     p[0] += 10;    
169     clitk::PointsUtils<MaskImageType>::Convert2DTo3D(p, VertebralArtery, i, B[i]);
170   }  
171
172   // Remove anterior to this line (keep -10 offset)
173   clitk::SliceBySliceSetBackgroundFromLineSeparation<MaskImageType>(S1R, A, B, 
174                                                                     GetBackgroundValue(), 1, -10);
175   m_ListOfStations["1R"] = S1R;
176
177   // Remove anterior to this line (keep -10 offset)
178   clitk::SliceBySliceSetBackgroundFromLineSeparation<MaskImageType>(S1L, A, B, 
179                                                                     GetBackgroundValue(), 1, -10);
180   StopCurrentStep<MaskImageType>(S1L);
181   m_ListOfStations["1L"] = S1L;
182 }
183 //--------------------------------------------------------------------
184
185