]> Creatis software - clitk.git/blob - segmentation/clitkExtractPatientFilter.h
changes in license header
[clitk.git] / segmentation / clitkExtractPatientFilter.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 CLITKEXTRACTPATIENTFILTER_H
20 #define CLITKEXTRACTPATIENTFILTER_H
21
22 #include "clitkFilterBase.h"
23 #include "clitkFilterWithAnatomicalFeatureDatabaseManagement.h"
24
25 namespace clitk {
26   
27   //--------------------------------------------------------------------
28   /*
29     Try to extract the Patient part of a thorax CT.  
30
31     Prefer high resolution input and resample (NN) output at the end
32     (like). Input is binarized using initial thresholds, connected
33     components are labeled (firstLabel). The air label (1) is
34     removed. The remaining is binarized and relabeled, patient should
35     now be the principal label (secondLabel). Two mechanismes are
36     provided to influence the label images. Crop to reduce
37     connectivity (image is restored to original size), eg for
38     SBF. Decomposition through ersion and reconstruction through
39     dilation (slow), eg for Pulmo bellows. Choose which labels to keep
40     from second Label image. Final mask is cleaned by opening and
41     closing.
42
43   */
44   //--------------------------------------------------------------------
45   
46   template <class TInputImageType>
47   class ITK_EXPORT ExtractPatientFilter: 
48     public virtual clitk::FilterBase, 
49     public clitk::FilterWithAnatomicalFeatureDatabaseManagement,
50     public itk::ImageToImageFilter<TInputImageType, 
51                                    itk::Image<uchar, TInputImageType::ImageDimension> > 
52   {
53   public:
54     /** Standard class typedefs. */
55     typedef itk::Image<uchar, TInputImageType::ImageDimension>      MaskImageType;
56     typedef ExtractPatientFilter                                    Self;
57     typedef itk::ImageToImageFilter<TInputImageType, MaskImageType> Superclass;
58     typedef itk::SmartPointer<Self>                                 Pointer;
59     typedef itk::SmartPointer<const Self>                           ConstPointer;
60     
61     /** Method for creation through the object factory. */
62     itkNewMacro(Self);
63     
64     /** Run-time type information (and related methods). */
65     itkTypeMacro(ExtractPatientFilter, ImageToImageFilter);
66     FILTERBASE_INIT;
67
68     /** Some convenient typedefs. */
69     typedef TInputImageType                       InputImageType;
70     typedef typename InputImageType::ConstPointer InputImageConstPointer;
71     typedef typename InputImageType::Pointer      InputImagePointer;
72     typedef typename InputImageType::RegionType   InputImageRegionType; 
73     typedef typename InputImageType::PixelType    InputImagePixelType; 
74     typedef typename InputImageType::SizeType     InputImageSizeType; 
75     typedef typename InputImageType::IndexType    InputImageIndexType; 
76         
77     typedef typename MaskImageType::ConstPointer MaskImageConstPointer;
78     typedef typename MaskImageType::Pointer      MaskImagePointer;
79     typedef typename MaskImageType::RegionType   MaskImageRegionType; 
80     typedef typename MaskImageType::PixelType    MaskImagePixelType; 
81     typedef typename MaskImageType::SizeType     MaskImageSizeType; 
82     typedef typename MaskImageType::IndexType    MaskImageIndexType; 
83
84     itkStaticConstMacro(ImageDimension, unsigned int, InputImageType::ImageDimension);
85     typedef int InternalPixelType;
86     typedef itk::Image<InternalPixelType, ImageDimension> InternalImageType;
87     typedef typename InternalImageType::SizeType          InternalImageSizeType;
88     
89     /** Connect inputs */
90     void SetInput(const TInputImageType * image);
91     itkSetMacro(OutputPatientFilename, std::string);
92     itkGetMacro(OutputPatientFilename, std::string);
93
94     // Step 1
95     itkSetMacro(UpperThreshold, InputImagePixelType);
96     itkGetMacro(UpperThreshold, InputImagePixelType);
97
98     itkSetMacro(LowerThreshold, InputImagePixelType);
99     itkGetMacro(LowerThreshold, InputImagePixelType);
100     itkSetMacro(UseLowerThreshold, bool);    
101     itkGetConstMacro(UseLowerThreshold, bool);    
102     itkBooleanMacro(UseLowerThreshold);
103
104     // Step 2
105     itkSetMacro(DecomposeAndReconstructDuringFirstStep, bool);
106     itkGetConstMacro(DecomposeAndReconstructDuringFirstStep, bool);
107     itkBooleanMacro(DecomposeAndReconstructDuringFirstStep);
108
109     itkSetMacro(Radius1, InternalImageSizeType);
110     itkGetConstMacro(Radius1, InternalImageSizeType);
111
112     itkSetMacro(MaximumNumberOfLabels1, int);
113     itkGetConstMacro(MaximumNumberOfLabels1, int);
114
115     itkSetMacro(NumberOfNewLabels1, int);
116     itkGetConstMacro(NumberOfNewLabels1, int);
117
118     // Step 2
119     itkSetMacro(DecomposeAndReconstructDuringSecondStep, bool);
120     itkGetConstMacro(DecomposeAndReconstructDuringSecondStep, bool);
121     itkBooleanMacro(DecomposeAndReconstructDuringSecondStep);
122
123     itkSetMacro(Radius2, InternalImageSizeType);
124     itkGetConstMacro(Radius2, InternalImageSizeType);
125
126     itkSetMacro(MaximumNumberOfLabels2, int);
127     itkGetConstMacro(MaximumNumberOfLabels2, int);
128
129     itkSetMacro(NumberOfNewLabels2, int);
130     itkGetConstMacro(NumberOfNewLabels2, int);
131
132     // Step 3
133     itkSetMacro(FirstKeep, int);
134     itkGetConstMacro(FirstKeep, int);
135
136     itkSetMacro(LastKeep, int);
137     itkGetConstMacro(LastKeep, int);
138
139     // Step 4
140     itkSetMacro(FinalOpenClose, bool);
141     itkGetConstMacro(FinalOpenClose, bool);
142     itkBooleanMacro(FinalOpenClose);
143
144     // Step 4
145     itkSetMacro(AutoCrop, bool);
146     itkGetConstMacro(AutoCrop, bool);
147     itkBooleanMacro(AutoCrop);
148
149   protected:
150     ExtractPatientFilter();
151     virtual ~ExtractPatientFilter() {}
152     
153     itkSetMacro(BackgroundValue, MaskImagePixelType);
154     itkSetMacro(ForegroundValue, MaskImagePixelType);
155     itkGetConstMacro(BackgroundValue, MaskImagePixelType);
156     itkGetConstMacro(ForegroundValue, MaskImagePixelType);
157     MaskImagePixelType m_BackgroundValue;
158     MaskImagePixelType m_ForegroundValue;
159
160     std::string m_OutputPatientFilename;
161     InputImagePixelType m_UpperThreshold;
162     InputImagePixelType m_LowerThreshold;
163     bool m_UseLowerThreshold;
164     bool m_DecomposeAndReconstructDuringFirstStep;
165     bool m_DecomposeAndReconstructDuringSecondStep;
166     bool m_FinalOpenClose;
167     InternalImageSizeType m_Radius1;
168     InternalImageSizeType m_Radius2;
169     int m_MaximumNumberOfLabels1;
170     int m_MaximumNumberOfLabels2;
171     int m_NumberOfNewLabels1;
172     int m_NumberOfNewLabels2;
173     int m_FirstKeep;
174     int m_LastKeep;
175     bool m_AutoCrop;
176
177     virtual void GenerateOutputInformation();
178     virtual void GenerateData();
179     
180     InputImageConstPointer input;
181     MaskImagePointer output;
182     typename InternalImageType::Pointer working_image;
183         
184   private:
185     ExtractPatientFilter(const Self&); //purposely not implemented
186     void operator=(const Self&); //purposely not implemented
187     
188   }; // end class
189   //--------------------------------------------------------------------
190
191 } // end namespace clitk
192 //--------------------------------------------------------------------
193
194 #ifndef ITK_MANUAL_INSTANTIATION
195 #include "clitkExtractPatientFilter.txx"
196 #endif
197
198 #endif