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