]> Creatis software - clitk.git/blob - segmentation/clitkExtractLungFilter.h
option to auto-detect motion mask's initial ellipse
[clitk.git] / segmentation / clitkExtractLungFilter.h
1 /*=========================================================================
2   Program:   vv                     http://www.creatis.insa-lyon.fr/rio/vv
3
4   Authors belong to: 
5   - University of LYON              http://www.universite-lyon.fr/
6   - Léon Bérard cancer center       http://www.centreleonberard.fr
7   - CREATIS CNRS laboratory         http://www.creatis.insa-lyon.fr
8
9   This software is distributed WITHOUT ANY WARRANTY; without even
10   the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
11   PURPOSE.  See the copyright notices for more information.
12
13   It is distributed under dual licence
14
15   - BSD        See included LICENSE.txt file
16   - CeCILL-B   http://www.cecill.info/licences/Licence_CeCILL-B_V1-en.html
17   ===========================================================================**/
18
19 #ifndef CLITKEXTRACTLUNGSFILTER_H
20 #define CLITKEXTRACTLUNGSFILTER_H
21
22 // clitk 
23 #include "clitkFilterBase.h"
24 #include "clitkDecomposeAndReconstructImageFilter.h"
25 #include "clitkExplosionControlledThresholdConnectedImageFilter.h"
26 #include "clitkSegmentationUtils.h"
27 #include "clitkFilterWithAnatomicalFeatureDatabaseManagement.h"
28
29 // itk
30 #include "itkStatisticsImageFilter.h"
31
32 namespace clitk {
33   
34   //--------------------------------------------------------------------
35   /*
36     Try to extract the Lung part of a thorax CT. Inspired by
37     Rikxoort2009, Section IIA, MedPhys.
38
39     - First, all air besides lungs and thrachea is removed, by
40     removing the second largest label of the firstLabelImage and
41     setting the remainder to 0HU . This modified input is optimally
42     thresholded (Otsu1979).
43
44     - Trachea and bronchi are grown from seeds in the top of the image
45     by explosion controlled region growing, slightly dilated and
46     removed from the second label image.
47
48     - Left and right lung are separated (if necessary) by erosion and
49     reconstructed by conditional dilation. 
50
51     - TRACHEA is available at the end
52
53     TODO ********** Remaining holes can be       filled afterwards (clitkFillMask).
54
55   */
56   //--------------------------------------------------------------------
57   
58   //--------------------------------------------------------------------
59   template <class TImageType>
60   class ITK_EXPORT ExtractLungFilter: 
61     public virtual clitk::FilterBase, 
62     public clitk::FilterWithAnatomicalFeatureDatabaseManagement,
63     public itk::ImageToImageFilter<TImageType, itk::Image<uchar, TImageType::ImageDimension> > 
64   {
65     
66   public:
67     /** Standard class typedefs. */
68     typedef itk::Image<uchar, TImageType::ImageDimension> MaskImageType;
69     typedef itk::ImageToImageFilter<TImageType, MaskImageType> Superclass;
70     typedef ExtractLungFilter              Self;
71     typedef itk::SmartPointer<Self>        Pointer;
72     typedef itk::SmartPointer<const Self>  ConstPointer;
73     
74     /** Method for creation through the object factory. */
75     itkNewMacro(Self);  
76     
77     /** Run-time type information (and related methods). */
78     itkTypeMacro(ExtractLungFilter, ImageToImageFilter);
79     FILTERBASE_INIT;
80
81     /** Some convenient typedefs */
82     typedef TImageType                       ImageType;
83     typedef typename ImageType::ConstPointer InputImageConstPointer;
84     typedef typename ImageType::Pointer      InputImagePointer;
85     typedef typename ImageType::RegionType   InputImageRegionType; 
86     typedef typename ImageType::PixelType    InputImagePixelType; 
87     typedef typename ImageType::SizeType     InputImageSizeType; 
88     typedef typename ImageType::IndexType    InputImageIndexType; 
89     typedef typename ImageType::PointType    InputImagePointType; 
90         
91     typedef typename MaskImageType::ConstPointer MaskImageConstPointer;
92     typedef typename MaskImageType::Pointer      MaskImagePointer;
93     typedef typename MaskImageType::RegionType   MaskImageRegionType; 
94     typedef typename MaskImageType::PixelType    MaskImagePixelType; 
95     typedef typename MaskImageType::SizeType     MaskImageSizeType; 
96     typedef typename MaskImageType::IndexType    MaskImageIndexType; 
97     typedef typename MaskImageType::PointType    MaskImagePointType; 
98
99     itkStaticConstMacro(ImageDimension, unsigned int, ImageType::ImageDimension);
100     //    typedef int InternalPixelType;
101     typedef uchar InternalPixelType;
102     typedef itk::Image<InternalPixelType, ImageType::ImageDimension> InternalImageType;
103     typedef typename InternalImageType::Pointer                      InternalImagePointer;
104     typedef typename InternalImageType::IndexType                    InternalIndexType;
105     typedef LabelizeParameters<InternalPixelType>                    LabelParamType;
106     
107     /** Connect inputs */
108     void SetInput(const ImageType * image);
109     itkSetMacro(PatientMaskBackgroundValue, MaskImagePixelType);
110     itkGetConstMacro(PatientMaskBackgroundValue, MaskImagePixelType);
111
112     // Output filename  (for AFBD)
113     itkSetMacro(OutputLungFilename, std::string);
114     itkGetMacro(OutputLungFilename, std::string);
115
116     itkSetMacro(OutputTracheaFilename, std::string);
117     itkGetMacro(OutputTracheaFilename, std::string);
118
119     // Get output (only availabe after update !)
120     typename MaskImageType::Pointer GetTracheaImage() { return trachea; }
121
122     // Background / Foreground
123     itkGetConstMacro(BackgroundValue, MaskImagePixelType);
124     itkGetConstMacro(ForegroundValue, MaskImagePixelType);
125
126     // For common segmentation processes
127     itkSetMacro(MinimalComponentSize, int);
128     itkGetConstMacro(MinimalComponentSize, int);
129
130     // Step 1 options RemoveAir
131     itkSetMacro(UpperThreshold, InputImagePixelType);
132     itkGetConstMacro(UpperThreshold, InputImagePixelType);
133
134     itkSetMacro(NumberOfSlicesToSkipBeforeSearchingSeed, int);
135     itkGetConstMacro(NumberOfSlicesToSkipBeforeSearchingSeed, int);
136     
137     itkSetMacro(LowerThreshold, InputImagePixelType);
138     itkGetConstMacro(LowerThreshold, InputImagePixelType);
139     itkSetMacro(UseLowerThreshold, bool);
140     itkGetConstMacro(UseLowerThreshold, bool);
141     itkBooleanMacro(UseLowerThreshold);
142
143     void SetLabelizeParameters1(LabelParamType * a) { m_LabelizeParameters1 = a; }
144     itkGetConstMacro(LabelizeParameters1, LabelParamType*);
145     
146     itkSetMacro(TracheaSeedAlgorithm, int);
147     itkGetConstMacro(TracheaSeedAlgorithm, int);
148
149     // Step 2 options FindTrachea
150     itkSetMacro(UpperThresholdForTrachea, InputImagePixelType);
151     itkGetConstMacro(UpperThresholdForTrachea, InputImagePixelType);
152
153     itkSetMacro(MultiplierForTrachea, double);
154     itkGetConstMacro(MultiplierForTrachea, double);
155
156     itkSetMacro(ThresholdStepSizeForTrachea, InputImagePixelType);
157     itkGetConstMacro(ThresholdStepSizeForTrachea, InputImagePixelType);
158
159     // options FindTrachea2
160     itkSetMacro(NumSlices, int);
161     itkGetConstMacro(NumSlices, int);
162     itkSetMacro(MaxElongation, double);
163     itkGetConstMacro(MaxElongation, double);
164     itkSetMacro(SeedPreProcessingThreshold, int);
165     itkGetConstMacro(SeedPreProcessingThreshold, int);
166
167     void AddSeed(InternalIndexType s);
168     std::vector<InternalIndexType> & GetSeeds() { return  m_Seeds; }
169
170     itkSetMacro(TracheaVolumeMustBeCheckedFlag, bool);
171     itkGetConstMacro(TracheaVolumeMustBeCheckedFlag, bool);
172     itkBooleanMacro(TracheaVolumeMustBeCheckedFlag);    
173
174     itkSetMacro(VerboseRegionGrowingFlag, bool);
175     itkGetConstMacro(VerboseRegionGrowingFlag, bool);
176     itkBooleanMacro(VerboseRegionGrowingFlag);    
177
178     // Step 3 options ExtractLung
179     itkSetMacro(NumberOfHistogramBins, int);
180     itkGetConstMacro(NumberOfHistogramBins, int);
181
182     void SetLabelizeParameters2(LabelParamType* a) { m_LabelizeParameters2 = a; }
183     itkGetConstMacro(LabelizeParameters2, LabelParamType*);
184
185     // Step 4 options RemoveTrachea
186     itkSetMacro(RadiusForTrachea, int);
187     itkGetConstMacro(RadiusForTrachea, int);
188     
189     void SetLabelizeParameters3(LabelParamType * a) { m_LabelizeParameters3 = a; }
190     itkGetConstMacro(LabelizeParameters3, LabelParamType*);
191
192     // Step 5 final openclose
193     itkSetMacro(OpenCloseFlag, bool);
194     itkGetConstMacro(OpenCloseFlag, bool);
195     itkBooleanMacro(OpenCloseFlag);
196
197     itkSetMacro(OpenCloseRadius, int);
198     itkGetConstMacro(OpenCloseRadius, int);
199     
200     // Step 6 fill holes
201     itkSetMacro(FillHolesFlag, bool);
202     itkGetConstMacro(FillHolesFlag, bool);
203     itkBooleanMacro(FillHolesFlag);
204
205     // Separate lungs
206     itkSetMacro(SeparateLungsFlag, bool);
207     itkGetConstMacro(SeparateLungsFlag, bool);
208     itkBooleanMacro(SeparateLungsFlag);
209
210     // Step Auto Crop
211     itkSetMacro(AutoCrop, bool);
212     itkGetConstMacro(AutoCrop, bool);
213     itkBooleanMacro(AutoCrop);
214     
215   protected:
216     ExtractLungFilter();
217     virtual ~ExtractLungFilter() {}
218
219     // Main members
220     InputImageConstPointer input;
221     MaskImagePointer patient;
222     InputImagePointer working_input;
223     std::string m_OutputLungFilename;
224     std::string m_OutputTracheaFilename;
225     MaskImagePointer working_mask;  
226     MaskImagePointer trachea;
227     unsigned int m_MaxSeedNumber;
228
229     // Global options
230     itkSetMacro(BackgroundValue, MaskImagePixelType);
231     itkSetMacro(ForegroundValue, MaskImagePixelType);
232     MaskImagePixelType m_PatientMaskBackgroundValue;
233     MaskImagePixelType m_BackgroundValue;
234     MaskImagePixelType m_ForegroundValue;
235     int m_MinimalComponentSize;
236     bool m_AutoCrop;
237
238     // Step 1
239     InputImagePixelType m_UpperThreshold;
240     InputImagePixelType m_LowerThreshold;
241     bool m_UseLowerThreshold;
242     LabelParamType* m_LabelizeParameters1;
243
244     // Step 2
245     int m_TracheaSeedAlgorithm;
246     InputImagePixelType m_UpperThresholdForTrachea;
247     InputImagePixelType m_ThresholdStepSizeForTrachea;
248     double m_MultiplierForTrachea;
249     std::vector<InternalIndexType> m_Seeds;
250     int m_NumberOfSlicesToSkipBeforeSearchingSeed;
251     bool m_TracheaVolumeMustBeCheckedFlag;
252     bool m_VerboseRegionGrowingFlag;
253     int m_NumSlices;
254     double m_MaxElongation;
255     int m_SeedPreProcessingThreshold;
256
257     // Step 3
258     int m_NumberOfHistogramBins;
259     LabelParamType* m_LabelizeParameters2;
260
261     // Step 4
262     int m_RadiusForTrachea;
263     LabelParamType* m_LabelizeParameters3;
264
265     // Step 5
266     bool m_OpenCloseFlag;    
267     int m_OpenCloseRadius;
268
269     // Step 6
270     bool m_FillHolesFlag;    
271     InputImageSizeType m_FillHolesDirections;
272
273     bool m_SeparateLungsFlag;
274     
275     // Main functions
276     virtual void GenerateOutputInformation();
277     virtual void GenerateInputRequestedRegion();
278     virtual void GenerateData();
279     
280     // Functions for trachea extraction
281     bool SearchForTracheaSeed(int skip);
282     bool SearchForTracheaSeed2(int numberOfSlices);
283     void SearchForTrachea();
284     void TracheaRegionGrowing();
285     double ComputeTracheaVolume();
286
287   private:
288     ExtractLungFilter(const Self&); //purposely not implemented
289     void operator=(const Self&); //purposely not implemented
290     
291   }; // end class
292   //--------------------------------------------------------------------
293
294 } // end namespace clitk
295 //--------------------------------------------------------------------
296
297 #ifndef ITK_MANUAL_INSTANTIATION
298 #include "clitkExtractLungFilter.txx"
299 #endif
300
301 #endif