]> Creatis software - clitk.git/blob - segmentation/clitkExtractLungFilter.h
76e2640f7d43d73b5586457678a9fe4c9e7c3388
[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://oncora1.lyon.fnclcc.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     // Step 2 options FindTrachea
147     itkSetMacro(UpperThresholdForTrachea, InputImagePixelType);
148     itkGetConstMacro(UpperThresholdForTrachea, InputImagePixelType);
149
150     itkSetMacro(MultiplierForTrachea, double);
151     itkGetConstMacro(MultiplierForTrachea, double);
152
153     itkSetMacro(ThresholdStepSizeForTrachea, InputImagePixelType);
154     itkGetConstMacro(ThresholdStepSizeForTrachea, InputImagePixelType);
155
156     void AddSeed(InternalIndexType s);
157     std::vector<InternalIndexType> & GetSeeds() { return  m_Seeds; }
158
159     itkSetMacro(TracheaVolumeMustBeCheckedFlag, bool);
160     itkGetConstMacro(TracheaVolumeMustBeCheckedFlag, bool);
161     itkBooleanMacro(TracheaVolumeMustBeCheckedFlag);    
162
163     itkSetMacro(VerboseRegionGrowingFlag, bool);
164     itkGetConstMacro(VerboseRegionGrowingFlag, bool);
165     itkBooleanMacro(VerboseRegionGrowingFlag);    
166
167     // Step 3 options ExtractLung
168     itkSetMacro(NumberOfHistogramBins, int);
169     itkGetConstMacro(NumberOfHistogramBins, int);
170
171     void SetLabelizeParameters2(LabelParamType* a) { m_LabelizeParameters2 = a; }
172     itkGetConstMacro(LabelizeParameters2, LabelParamType*);
173
174     // Step 4 options RemoveTrachea
175     itkSetMacro(RadiusForTrachea, int);
176     itkGetConstMacro(RadiusForTrachea, int);
177     
178     void SetLabelizeParameters3(LabelParamType * a) { m_LabelizeParameters3 = a; }
179     itkGetConstMacro(LabelizeParameters3, LabelParamType*);
180
181     // Step 5 final openclose
182     itkSetMacro(OpenCloseFlag, bool);
183     itkGetConstMacro(OpenCloseFlag, bool);
184     itkBooleanMacro(OpenCloseFlag);
185
186     itkSetMacro(OpenCloseRadius, int);
187     itkGetConstMacro(OpenCloseRadius, int);
188     
189     // Step 6 fill holes
190     itkSetMacro(FillHolesFlag, bool);
191     itkGetConstMacro(FillHolesFlag, bool);
192     itkBooleanMacro(FillHolesFlag);
193
194   protected:
195     ExtractLungFilter();
196     virtual ~ExtractLungFilter() {}
197
198     // Main members
199     InputImageConstPointer input;
200     MaskImagePointer patient;
201     InputImagePointer working_input;
202     std::string m_OutputLungFilename;
203     std::string m_OutputTracheaFilename;
204     MaskImagePointer working_mask;  
205     MaskImagePointer trachea;
206     unsigned int m_MaxSeedNumber;
207
208     // Global options
209     itkSetMacro(BackgroundValue, MaskImagePixelType);
210     itkSetMacro(ForegroundValue, MaskImagePixelType);
211     MaskImagePixelType m_PatientMaskBackgroundValue;
212     MaskImagePixelType m_BackgroundValue;
213     MaskImagePixelType m_ForegroundValue;
214     int m_MinimalComponentSize;
215
216     // Step 1
217     InputImagePixelType m_UpperThreshold;
218     InputImagePixelType m_LowerThreshold;
219     bool m_UseLowerThreshold;
220     LabelParamType* m_LabelizeParameters1;
221
222     // Step 2
223     InputImagePixelType m_UpperThresholdForTrachea;
224     InputImagePixelType m_ThresholdStepSizeForTrachea;
225     double m_MultiplierForTrachea;
226     std::vector<InternalIndexType> m_Seeds;
227     int m_NumberOfSlicesToSkipBeforeSearchingSeed;
228     bool m_TracheaVolumeMustBeCheckedFlag;
229     bool m_VerboseRegionGrowingFlag;
230
231     // Step 3
232     int m_NumberOfHistogramBins;
233     LabelParamType* m_LabelizeParameters2;
234
235     // Step 4
236     int m_RadiusForTrachea;
237     LabelParamType* m_LabelizeParameters3;
238
239     // Step 5
240     bool m_OpenCloseFlag;    
241     int m_OpenCloseRadius;
242
243     // Step 6
244     bool m_FillHolesFlag;    
245     InputImageSizeType m_FillHolesDirections;
246
247     // Main functions
248     virtual void GenerateOutputInformation();
249     virtual void GenerateInputRequestedRegion();
250     virtual void GenerateData();
251     
252     // Functions for trachea extraction
253     bool SearchForTracheaSeed(int skip);
254     void SearchForTrachea();
255     void TracheaRegionGrowing();
256     double ComputeTracheaVolume();
257
258   private:
259     ExtractLungFilter(const Self&); //purposely not implemented
260     void operator=(const Self&); //purposely not implemented
261     
262   }; // end class
263   //--------------------------------------------------------------------
264
265 } // end namespace clitk
266 //--------------------------------------------------------------------
267
268 #ifndef ITK_MANUAL_INSTANTIATION
269 #include "clitkExtractLungFilter.txx"
270 #endif
271
272 #endif