]> Creatis software - clitk.git/blob - tools/clitkImageStatisticsGenericFilter.txx
Merge branch 'master' into GammaIndex3D
[clitk.git] / tools / clitkImageStatisticsGenericFilter.txx
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 #ifndef clitkImageStatisticsGenericFilter_txx
19 #define clitkImageStatisticsGenericFilter_txx
20
21 #include "itkNthElementImageAdaptor.h"
22 #include "itkJoinSeriesImageFilter.h"
23
24 #include "clitkImageStatisticsGenericFilter.h"
25 #include "clitkCropLikeImageFilter.h"
26 #include "clitkResampleImageWithOptionsFilter.h"
27
28 namespace clitk
29 {
30
31   //-------------------------------------------------------------------
32   // Update with the number of dimensions
33   //-------------------------------------------------------------------
34   template<unsigned int Dimension, unsigned int Components>
35   void
36   ImageStatisticsGenericFilter::UpdateWithDim(std::string PixelType)
37   {
38     if (m_Verbose) std::cout << "Image was detected to be "<<Dimension<<"D and "<< PixelType<<" with " << Components << " channel(s)..."<<std::endl;
39
40     if(PixelType == "short"){
41       if (m_Verbose) std::cout << "Launching filter in "<< Dimension <<"D and signed short..." << std::endl;
42       UpdateWithDimAndPixelType<Dimension, signed short, Components>();
43     }
44     else if(PixelType == "unsigned_short"){
45       if (m_Verbose) std::cout  << "Launching filter in "<< Dimension <<"D and unsigned_short..." << std::endl;
46       UpdateWithDimAndPixelType<Dimension, unsigned short, Components>();
47     }
48
49     else if (PixelType == "unsigned_char"){
50       if (m_Verbose) std::cout  << "Launching filter in "<< Dimension <<"D and unsigned_char..." << std::endl;
51       UpdateWithDimAndPixelType<Dimension, unsigned char, Components>();
52     }
53
54     else if(PixelType == "double"){
55       if (m_Verbose) std::cout << "Launching filter in "<< Dimension <<"D and double..." << std::endl;
56       UpdateWithDimAndPixelType<Dimension, double, Components>();
57     }
58     else {
59       if (m_Verbose) std::cout  << "Launching filter in "<< Dimension <<"D and float..." << std::endl;
60       UpdateWithDimAndPixelType<Dimension, float, Components>();
61     }
62   }
63
64
65   //-------------------------------------------------------------------
66   // Update with the number of dimensions and the pixeltype
67   //-------------------------------------------------------------------
68   template <unsigned int Dimension, class  PixelType, unsigned int Components>
69   void
70   ImageStatisticsGenericFilter::UpdateWithDimAndPixelType()
71   {
72
73     // ImageTypes
74     typedef unsigned char LabelPixelType;
75     typedef itk::Image<itk::Vector<PixelType, Components>, Dimension> InputImageType;
76     typedef itk::Image<LabelPixelType, Dimension> LabelImageType;
77
78     // Read the input
79     typedef itk::ImageFileReader<InputImageType> InputReaderType;
80     typename InputReaderType::Pointer reader = InputReaderType::New();
81     reader->SetFileName( m_InputFileName);
82     reader->Update();
83     typename InputImageType::Pointer input= reader->GetOutput();
84
85     typedef itk::NthElementImageAdaptor<InputImageType, PixelType> InputImageAdaptorType;
86     typedef itk::Image<PixelType, Dimension> OutputImageType;
87
88     typename InputImageAdaptorType::Pointer input_adaptor = InputImageAdaptorType::New();
89     input_adaptor->SetImage(input);
90
91     // Filter
92     typedef itk::LabelStatisticsImageFilter<InputImageAdaptorType, LabelImageType> StatisticsImageFilterType;
93     typename StatisticsImageFilterType::Pointer statisticsFilter=StatisticsImageFilterType::New();
94     statisticsFilter->SetInput(input_adaptor);
95
96     // Label image
97     typename LabelImageType::Pointer labelImage;
98     if (m_ArgsInfo.mask_given) {
99       int maskDimension, maskComponents;
100       std::string maskPixelType;
101       ReadImageDimensionAndPixelType(m_ArgsInfo.mask_arg, maskDimension, maskPixelType, maskComponents);
102
103       if (maskDimension == Dimension - 1) {
104         // Due to a limitation of filter itk::LabelStatisticsImageFilter, InputImageType and LabelImageType
105         // must have the same image dimension. However, we want to support label images with Dl = Di - 1,
106         // so we need to replicate the label image as many times as the size along dimension Di.
107         if (m_Verbose)
108           std::cout << "Replicating label image to match input image's dimension... " << std::endl;
109
110         typedef itk::Image<LabelPixelType, Dimension - 1> ReducedLabelImageType;
111         typedef itk::ImageFileReader<ReducedLabelImageType> LabelImageReaderType;
112         typedef itk::JoinSeriesImageFilter<ReducedLabelImageType, LabelImageType> JoinImageFilterType;
113
114         typename LabelImageReaderType::Pointer labelImageReader=LabelImageReaderType::New();
115         labelImageReader->SetFileName(m_ArgsInfo.mask_arg);
116         labelImageReader->Update();
117
118         typename JoinImageFilterType::Pointer joinFilter = JoinImageFilterType::New();
119         typename InputImageType::SizeType size = input->GetLargestPossibleRegion().GetSize();
120         for (unsigned int i = 0; i < size[Dimension - 1]; i++)
121           joinFilter->PushBackInput(labelImageReader->GetOutput());
122
123         joinFilter->Update();
124         labelImage = joinFilter->GetOutput();
125       }
126       else {
127         typedef itk::ImageFileReader<LabelImageType> LabelImageReaderType;
128         typename LabelImageReaderType::Pointer labelImageReader=LabelImageReaderType::New();
129         labelImageReader->SetFileName(m_ArgsInfo.mask_arg);
130         labelImageReader->Update();
131         labelImage= labelImageReader->GetOutput();
132
133         // Check mask sampling/size
134         if (!HaveSameSizeAndSpacing<LabelImageType, InputImageType>(labelImage, input)) {
135           if (m_ArgsInfo.allow_resize_flag) {
136             if (m_ArgsInfo.verbose_flag) {
137               std::cout << "Resize mask image like input" << std::endl;
138             }
139             typedef clitk::ResampleImageWithOptionsFilter<LabelImageType> ResamplerType;
140             typename ResamplerType::Pointer resampler = ResamplerType::New();
141             resampler->SetInput(labelImage);
142             resampler->SetOutputSpacing(input->GetSpacing());
143             resampler->SetOutputOrigin(labelImage->GetOrigin());
144             resampler->SetGaussianFilteringEnabled(false);
145             resampler->Update();
146             labelImage = resampler->GetOutput();
147             //writeImage<LabelImageType>(labelImage, "test1.mha");
148
149             typedef clitk::CropLikeImageFilter<LabelImageType> FilterType;
150             typename FilterType::Pointer crop = FilterType::New();
151             crop->SetInput(labelImage);
152             crop->SetCropLikeImage(input);
153             crop->Update();
154             labelImage = crop->GetOutput();
155             //writeImage<LabelImageType>(labelImage, "test2.mha");
156
157           }
158           else {
159             std::cerr << "Mask image has a different size/spacing than input. Abort. (Use option to resize)" << std::endl;
160             exit(-1);
161           }
162         }
163
164       }
165
166     }
167     else {
168       labelImage=LabelImageType::New();
169       labelImage->SetDirection(input->GetDirection());
170       labelImage->SetRegions(input->GetLargestPossibleRegion());
171       labelImage->SetOrigin(input->GetOrigin());
172       labelImage->SetSpacing(input->GetSpacing());
173       labelImage->SetDirection(input->GetDirection());
174       labelImage->Allocate();
175       labelImage->FillBuffer(m_ArgsInfo.label_arg[0]);
176     }
177     statisticsFilter->SetLabelInput(labelImage);
178
179     // Check/compute spacing
180     const typename LabelImageType::SpacingType& spacing = input->GetSpacing();
181     double spacing_cc = (spacing[0]*spacing[1]*spacing[2])/1000;
182     // std::cout<<"Spacing x : " << spacing[0]<<std::endl;
183     // std::cout<<"Spacing y :  " << spacing[1]<<std::endl;
184     // std::cout<<"Spacing z :  " << spacing[2]<<std::endl;
185     // std::cout <<"spacing_cc : " << spacing_cc << std::endl;
186
187
188     // For each Label
189     typename LabelImageType::PixelType label;
190     unsigned int numberOfLabels;
191     if (m_ArgsInfo.label_given)
192       numberOfLabels=m_ArgsInfo.label_given;
193     else
194       numberOfLabels=1;
195
196     unsigned int firstComponent = 0, lastComponent = 0;
197     if (m_ArgsInfo.channel_arg == -1) {
198       firstComponent = 0;
199       lastComponent = Components - 1;
200     }
201     else {
202       firstComponent = m_ArgsInfo.channel_arg;
203       lastComponent = m_ArgsInfo.channel_arg;
204     }
205
206     for (unsigned int c=firstComponent; c<=lastComponent; c++) {
207       if (m_Verbose) std::cout << std::endl << "Processing channel " << c << std::endl;
208
209       input_adaptor->SelectNthElement(c);
210       input_adaptor->Update();
211
212       for (unsigned int k=0; k< numberOfLabels; k++) {
213         label=m_ArgsInfo.label_arg[k];
214
215         std::cout<<std::endl;
216         if (m_Verbose) std::cout<<"-------------"<<std::endl;
217         if (m_Verbose) std::cout<<"| Label: "<< (int) label<<"  |"<<std::endl;
218         if (m_Verbose) std::cout<<"-------------"<<std::endl;
219
220         // Histograms
221         if (m_ArgsInfo.histogram_given) {
222           statisticsFilter->SetUseHistograms(true);
223           statisticsFilter->SetHistogramParameters(m_ArgsInfo.bins_arg, m_ArgsInfo.lower_arg, m_ArgsInfo.upper_arg);
224         }
225
226         // DVH
227         if(m_ArgsInfo.dvhistogram_given)
228         {
229           statisticsFilter->SetUseHistograms(true);
230           statisticsFilter->SetHistogramParameters(m_ArgsInfo.bins_arg, m_ArgsInfo.lower_arg, m_ArgsInfo.upper_arg);
231         }
232
233         statisticsFilter->Update();
234
235         // Output
236        if (m_Verbose) std::cout<<"N° of pixels: ";
237         std::cout<<statisticsFilter->GetCount(label)<<std::endl;
238         if (m_Verbose) std::cout<<"Mean: ";
239         std::cout<<statisticsFilter->GetMean(label)<<std::endl;
240         if (m_Verbose) std::cout<<"SD: ";
241         std::cout<<statisticsFilter->GetSigma(label)<<std::endl;
242         if (m_Verbose) std::cout<<"Variance: ";
243         std::cout<<statisticsFilter->GetVariance(label)<<std::endl;
244         if (m_Verbose) std::cout<<"Min: ";
245         std::cout<<statisticsFilter->GetMinimum(label)<<std::endl;
246         if (m_Verbose) std::cout<<"Max: ";
247         std::cout<<statisticsFilter->GetMaximum(label)<<std::endl;
248         if (m_Verbose) std::cout<<"Sum: ";
249         std::cout<<statisticsFilter->GetSum(label)<<std::endl;
250         if (m_Verbose) std::cout<<"Volume (cc): ";
251         std::cout<<statisticsFilter->GetCount(label)*spacing_cc<<std::endl;
252         if (m_Verbose) std::cout<<"Bounding box: ";
253         for(unsigned int i =0; i <statisticsFilter->GetBoundingBox(label).size(); i++)
254             std::cout<<statisticsFilter->GetBoundingBox(label)[i]<<" ";
255         std::cout<<std::endl;
256
257         // Histogram
258         if (m_ArgsInfo.histogram_given)
259         {
260           if (m_Verbose) std::cout<<"Median: ";
261           std::cout<<statisticsFilter->GetMedian(label)<<std::endl;
262
263           typename StatisticsImageFilterType::HistogramPointer histogram =statisticsFilter->GetHistogram(label);
264
265           // Screen
266           if (m_Verbose) std::cout<<"Histogram: "<<std::endl;
267             std::cout<<"# MinBin\tMidBin\tMaxBin\tFrequency"<<std::endl;
268           for( int i =0; i <m_ArgsInfo.bins_arg; i++)
269             std::cout<<histogram->GetBinMin(0,i)<<"\t"<<histogram->GetMeasurement(i,0)<<"\t"<<histogram->GetBinMax(0,i)<<"\t"<<histogram->GetFrequency(i)<<std::endl;
270
271           // Add to the file
272           std::ofstream histogramFile(m_ArgsInfo.histogram_arg);
273           histogramFile<<"#Histogram: "<<std::endl;
274           histogramFile<<"#MinBin\tMidBin\tMaxBin\tFrequency"<<std::endl;
275           for( int i =0; i <m_ArgsInfo.bins_arg; i++)
276             histogramFile<<histogram->GetBinMin(0,i)<<"\t"<<histogram->GetMeasurement(i,0)<<"\t"<<histogram->GetBinMax(0,i)<<"\t"<<histogram->GetFrequency(i)<<std::endl;
277         }
278
279         // DVH
280         if(m_ArgsInfo.dvhistogram_given)
281         {
282             typename StatisticsImageFilterType::HistogramPointer dvhistogram = statisticsFilter->GetHistogram(label);
283             double totalVolumeCC = ((statisticsFilter->GetCount(label))*spacing_cc);
284             double totalVolume = statisticsFilter->GetCount(label);
285
286             // Screen
287             std::cout<<"# Total volume : ";
288             std::cout<<totalVolume<<" [No. of voxels]"<<std::endl;
289             std::cout<<"# Total volume : ";
290             std::cout<<totalVolumeCC<<" [cc]"<<std::endl;
291             std::cout<<"# Dose mean: ";
292             std::cout<<statisticsFilter->GetMean(label)<<" [Gy]"<<std::endl;
293             std::cout<<"# Dose min: ";
294             std::cout<<statisticsFilter->GetMinimum(label)<<" [Gy]"<<std::endl;
295             std::cout<<"# Dose max: ";
296             std::cout<<statisticsFilter->GetMaximum(label)<<" [Gy]"<<std::endl;
297             std::cout<<" "<<std::endl;
298             std::cout<<"#Dose_diff[Gy] Volume_diff[No. of voxels] Volume_diff[%] Volume_diff[cc] Volume_cumul[No. of voxels] Volume_cumul[%] Volume_cumul[cc]"<<std::endl;
299             for( int i =0; i <m_ArgsInfo.bins_arg; i++)
300             {
301               double popCumulativeVolume = 0;
302               for(int j=0; j<i; j++)
303               {
304                  popCumulativeVolume+=(dvhistogram->GetFrequency(j));
305               }
306               double cumulativeVolume = (totalVolume - (popCumulativeVolume + (dvhistogram->GetFrequency(i))));
307               double percentCumulativeVolume =(cumulativeVolume*100)/(statisticsFilter->GetCount(label)) ;
308               double ccCumulativeVolume = (totalVolumeCC -((popCumulativeVolume + (dvhistogram->GetFrequency(i)))*spacing_cc));
309               double percentDiffVolume = dvhistogram->GetFrequency(i)*100/(statisticsFilter->GetCount(label));
310               std::cout<<dvhistogram->GetBinMax(0,i)<<"\t  "<<dvhistogram->GetFrequency(i)<<"\t  "<<percentDiffVolume<<"\t "<<((dvhistogram->GetFrequency(i))*spacing_cc)<<"\t "<<"\t "<<cumulativeVolume<<"\t  "<<percentCumulativeVolume<<"\t  "<<ccCumulativeVolume<<"\t "<<std::endl;
311             }
312
313             // Add to the file
314             std::ofstream dvhistogramFile(m_ArgsInfo.dvhistogram_arg);
315             dvhistogramFile<<"# Total volume : ";
316             dvhistogramFile<<statisticsFilter->GetCount(label)<<" [No. of voxels]"<<std::endl;
317             dvhistogramFile<<"# Total volume : ";
318             dvhistogramFile<<totalVolumeCC<<" [cc]"<<std::endl;
319             dvhistogramFile<<"# Dose mean: ";
320             dvhistogramFile<<statisticsFilter->GetMean(label)<<" [Gy]"<<std::endl;
321             dvhistogramFile<<"# Dose min: ";
322             dvhistogramFile<<statisticsFilter->GetMinimum(label)<<" [Gy]"<<std::endl;
323             dvhistogramFile<<"# Dose max: ";
324             dvhistogramFile<<statisticsFilter->GetMaximum(label)<<" [Gy]"<<std::endl;
325             dvhistogramFile<<"  "<<std::endl;
326             dvhistogramFile<<"#Dose_diff[Gy] Volume_diff[No. of voxels] Volume_diff[%] Volume_diff[cc] Volume_cumul[No. of voxels] Volume_cumul[%] Volume_cumul[cc]"<<std::endl;
327             for( int i =0; i <m_ArgsInfo.bins_arg; i++)
328             {
329                double popCumulativeVolume = 0;
330                for(int j=0; j<i; j++)
331                {
332                  popCumulativeVolume+=(dvhistogram->GetFrequency(j));
333                }
334                double cumulativeVolume = (totalVolume - (popCumulativeVolume + (dvhistogram->GetFrequency(i))));
335                double percentCumulativeVolume =(cumulativeVolume*100)/(statisticsFilter->GetCount(label));
336                double ccCumulativeVolume = (totalVolumeCC -((popCumulativeVolume + (dvhistogram->GetFrequency(i)))*spacing_cc));
337                double percentDiffVolume = ((dvhistogram->GetFrequency(i))*100)/(statisticsFilter->GetCount(label));
338                dvhistogramFile<<dvhistogram->GetBinMax(0,i)<<"\t  "<<dvhistogram->GetFrequency(i)<<"\t  "<<percentDiffVolume<<"\t "<<((dvhistogram->GetFrequency(i))*spacing_cc)<<"\t "<<cumulativeVolume<<"\t "<<percentCumulativeVolume<<"\t  "<<ccCumulativeVolume<<std::endl;
339            }
340         }
341       }
342     }
343
344     return;
345
346   }
347
348
349 }//end clitk
350
351 #endif //#define clitkImageStatisticsGenericFilter_txx