]> Creatis software - clitk.git/blob - itk/clitkReconstructThroughDilationImageFilter.txx
d2d3d31f5eae738967a25afef41635be3df0ec62
[clitk.git] / itk / clitkReconstructThroughDilationImageFilter.txx
1 #ifndef clitkReconstructThroughDilationImageFilter_txx
2 #define clitkReconstructThroughDilationImageFilter_txx
3
4 /* =================================================
5  * @file   clitkReconstructThroughDilationImageFilter.txx
6  * @author 
7  * @date   
8  * 
9  * @brief 
10  * 
11  ===================================================*/
12
13
14 namespace clitk
15 {
16
17   //-------------------------------------------------------------------
18   // Update with the number of dimensions
19   //-------------------------------------------------------------------
20   template<class InputImageType, class OutputImageType>
21   ReconstructThroughDilationImageFilter<InputImageType, OutputImageType>::ReconstructThroughDilationImageFilter()
22   {
23    m_Verbose=false;
24    m_BackgroundValue=0;
25    m_ForegroundValue=1;
26    m_ErosionPaddingValue=static_cast<InputPixelType>(-1);
27    for (unsigned int i=0; i<InputImageDimension; i++)
28      m_Radius[i]=1;
29    m_MaximumNumberOfLabels=10;
30   }
31
32
33   //-------------------------------------------------------------------
34   // Update with the number of dimensions and the pixeltype
35   //-------------------------------------------------------------------
36   template<class InputImageType, class  OutputImageType> 
37   void 
38   ReconstructThroughDilationImageFilter<InputImageType, OutputImageType>::GenerateData()
39   {
40
41     //---------------------------------
42     // Typedefs 
43     //--------------------------------- 
44     
45     // Internal type
46     typedef itk::Image<InternalPixelType, InputImageDimension> InternalImageType;
47
48     // Filters used
49     typedef itk::CastImageFilter<InputImageType, InternalImageType> InputCastImageFilterType;   
50     typedef itk::ThresholdImageFilter<InternalImageType> InputThresholdImageFilterType;
51     typedef itk::StatisticsImageFilter<InternalImageType> StatisticsImageFilterType;
52     typedef itk::BinaryBallStructuringElement<InternalPixelType,InputImageDimension > KernelType;
53     typedef clitk::ConditionalBinaryDilateImageFilter<InternalImageType, InternalImageType , KernelType> ConditionalBinaryDilateImageFilterType;
54     typedef itk::DifferenceImageFilter<InternalImageType, InternalImageType> DifferenceImageFilterType;
55     typedef itk::CastImageFilter<InternalImageType, OutputImageType> OutputCastImageFilterType;
56     typedef clitk::SetBackgroundImageFilter<InternalImageType, InternalImageType, InternalImageType> SetBackgroundImageFilterType;
57
58     //---------------------------------
59     // Cast
60     //---------------------------------
61     typename InputCastImageFilterType::Pointer castImageFilter=InputCastImageFilterType::New();
62     castImageFilter->SetInput(this->GetInput());
63     castImageFilter->Update();
64
65     //---------------------------------
66     // Threshold
67     //---------------------------------
68     typename InputThresholdImageFilterType::Pointer thresholdImageFilter=InputThresholdImageFilterType::New();
69     thresholdImageFilter->SetInput(castImageFilter->GetOutput());
70     thresholdImageFilter->ThresholdAbove(m_MaximumNumberOfLabels);
71     thresholdImageFilter->SetOutsideValue(m_ForegroundValue);
72     if(m_Verbose) std::cout<<"Thresholding the input to "<<m_MaximumNumberOfLabels<<" labels ..."<<std::endl;
73     thresholdImageFilter->Update();
74
75     //---------------------------------
76     // Set -1 to padding value
77     //---------------------------------
78     typename SetBackgroundImageFilterType::Pointer setBackgroundFilter =SetBackgroundImageFilterType::New();
79     setBackgroundFilter->SetInput(thresholdImageFilter->GetOutput());
80     setBackgroundFilter->SetInput2(castImageFilter->GetOutput());
81     setBackgroundFilter->SetMaskValue(m_ErosionPaddingValue);
82     setBackgroundFilter->SetOutsideValue(-1);
83     if(m_Verbose) std::cout<<"Setting the eroded region from "<<m_ErosionPaddingValue<<" to -1..."<<std::endl;
84
85     
86     //---------------------------------
87     // Count the initial labels
88     //---------------------------------
89     typename StatisticsImageFilterType::Pointer inputStatisticsImageFilter=StatisticsImageFilterType::New();
90     inputStatisticsImageFilter->SetInput(setBackgroundFilter->GetOutput());
91     if(m_Verbose) std::cout<<"Counting the initial labels..."<<std::endl;
92     inputStatisticsImageFilter->Update();
93     unsigned int initialNumberOfLabels= inputStatisticsImageFilter->GetMaximum();
94     if(m_Verbose) std::cout<<"The input contained "<<initialNumberOfLabels<<" disctictive label(s)..."<<std::endl;
95     unsigned int numberOfConsideredLabels=std::min(initialNumberOfLabels, m_MaximumNumberOfLabels);
96     if(m_Verbose) std::cout<<"Performing dilation the first "<<numberOfConsideredLabels<<" disctictive labels..."<<std::endl;
97
98     //---------------------------------
99     // Dilate while change
100     //---------------------------------
101     typename itk::NumericTraits<InputPixelType>::AccumulateType difference=1;
102     typename InternalImageType::Pointer labelImage=inputStatisticsImageFilter->GetOutput();
103     typename InternalImageType::Pointer oldLabelImage=inputStatisticsImageFilter->GetOutput();
104
105     // element
106     KernelType structuringElement;
107     structuringElement.SetRadius(m_Radius);
108     structuringElement.CreateStructuringElement();
109
110     while( difference)
111       {
112         // Dilate all labels once
113         for ( int label=0; label<(int)numberOfConsideredLabels+1;label++)  
114           if ( m_BackgroundValue != label)
115             {
116               typename ConditionalBinaryDilateImageFilterType::Pointer dilateFilter=ConditionalBinaryDilateImageFilterType::New();
117               dilateFilter->SetBoundaryToForeground(false);
118               dilateFilter->SetKernel(structuringElement);
119               dilateFilter->SetBackgroundValue (-1);
120               dilateFilter->SetInput (labelImage);
121               dilateFilter->SetForegroundValue (label);
122               if(m_Verbose) std::cout<<"Dilating the label "<<label<<"..."<<std::endl;
123               dilateFilter->Update();
124               labelImage=dilateFilter->GetOutput();
125             }
126   
127         // Difference with previous labelImage
128         typename DifferenceImageFilterType::Pointer differenceFilter=DifferenceImageFilterType::New();
129         differenceFilter->SetValidInput(oldLabelImage);
130         differenceFilter->SetTestInput(labelImage);
131         differenceFilter->Update();
132         difference =differenceFilter->GetTotalDifference();
133         if(m_Verbose) std::cout<<"The change in this iteration was "<<difference<<"..."<<std::endl;     
134         oldLabelImage=labelImage;
135       }
136       
137     //---------------------------------
138     // Set -1 to padding value
139     //---------------------------------
140     typename SetBackgroundImageFilterType::Pointer setBackgroundFilter2 =SetBackgroundImageFilterType::New();
141     setBackgroundFilter2->SetInput(labelImage);
142     setBackgroundFilter2->SetInput2(labelImage);
143     setBackgroundFilter2->SetMaskValue(-1);
144     setBackgroundFilter2->SetOutsideValue(m_ErosionPaddingValue);
145     if(m_Verbose) std::cout<<"Setting the eroded region to "<<m_ErosionPaddingValue<<"..."<<std::endl;
146
147     //---------------------------------
148     // Cast
149     //---------------------------------
150     typename OutputCastImageFilterType::Pointer outputCastImageFilter=OutputCastImageFilterType::New();
151     outputCastImageFilter->SetInput(setBackgroundFilter2->GetOutput());
152     if(m_Verbose) std::cout<<"Casting the output..."<<std::endl;
153     outputCastImageFilter->Update();
154
155     //---------------------------------
156     // SetOutput
157     //---------------------------------
158     this->SetNthOutput(0, outputCastImageFilter->GetOutput());
159
160
161   }
162
163
164 }//end clitk
165  
166 #endif //#define clitkReconstructThroughDilationImageFilter_txx