]> Creatis software - clitk.git/blob - tools/clitkAverageTemporalDimensionGenericFilter.txx
473e36fee73ab9fdcb66a38205f71b2740be4cb0
[clitk.git] / tools / clitkAverageTemporalDimensionGenericFilter.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://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 #ifndef clitkAverageTemporalDimensionGenericFilter_txx
19 #define clitkAverageTemporalDimensionGenericFilter_txx
20
21 /* =================================================
22  * @file   clitkAverageTemporalDimensionGenericFilter.txx
23  * @author 
24  * @date   
25  * 
26  * @brief 
27  * 
28  ===================================================*/
29
30
31 namespace clitk
32 {
33
34   //-----------------------------------------------------------
35   // Constructor
36   //-----------------------------------------------------------
37   template<class args_info_type>
38   AverageTemporalDimensionGenericFilter<args_info_type>::AverageTemporalDimensionGenericFilter()
39   {
40     m_Verbose=false;
41     m_InputFileName="";
42   }
43
44
45   //-----------------------------------------------------------
46   // Update
47   //-----------------------------------------------------------
48   template<class args_info_type>
49   void AverageTemporalDimensionGenericFilter<args_info_type>::Update()
50   {
51     // Read the Dimension and PixelType
52     int Dimension, Components;
53     std::string PixelType;
54     ReadImageDimensionAndPixelType(m_InputFileName, Dimension, PixelType, Components);
55
56     
57     // Call UpdateWithDim
58     if (m_ArgsInfo.input_given>1) Dimension+=1;
59     if(Dimension==2) UpdateWithDim<2>(PixelType, Components);
60     else if(Dimension==3)UpdateWithDim<3>(PixelType, Components);
61     else if (Dimension==4)UpdateWithDim<4>(PixelType, Components); 
62     else 
63       {
64         std::cout<<"Error, Only for 2, 3 or 4D!!!"<<std::endl ;
65         return;
66       }
67   }
68
69   //-------------------------------------------------------------------
70   // Update with the number of dimensions
71   //-------------------------------------------------------------------
72   template<class args_info_type>
73   template<unsigned int Dimension>
74   void 
75   AverageTemporalDimensionGenericFilter<args_info_type>::UpdateWithDim(const std::string PixelType, const int Components)
76   {
77     if (m_Verbose) std::cout << "Image was detected to be "<<Dimension<<"D with "<<Components<<" component(s) of "<<  PixelType<<"..."<<std::endl;
78
79     if (Components==1)
80       {
81         if(PixelType == "short"){  
82           if (m_Verbose) std::cout << "Launching filter in "<< Dimension <<"D and signed short..." << std::endl;
83           UpdateWithDimAndPixelType<Dimension, signed short>(); 
84         }
85         //    else if(PixelType == "unsigned_short"){  
86         //       if (m_Verbose) std::cout  << "Launching filter in "<< Dimension <<"D and unsigned_short..." << std::endl;
87         //       UpdateWithDimAndPixelType<Dimension, unsigned short>(); 
88         //     }
89         
90         else if (PixelType == "unsigned_char"){ 
91           if (m_Verbose) std::cout  << "Launching filter in "<< Dimension <<"D and unsigned_char..." << std::endl;
92           UpdateWithDimAndPixelType<Dimension, unsigned char>();
93         }
94         
95         //     else if (PixelType == "char"){ 
96         //       if (m_Verbose) std::cout  << "Launching filter in "<< Dimension <<"D and signed_char..." << std::endl;
97         //       UpdateWithDimAndPixelType<Dimension, signed char>();
98         //     }
99         else {
100           if (m_Verbose) std::cout  << "Launching filter in "<< Dimension <<"D and float..." << std::endl;
101           UpdateWithDimAndPixelType<Dimension, float>();
102         }
103       }
104
105     //     else if (Components==2)
106     //       {
107     //  if (m_Verbose) std::cout  << "Launching transform in "<< Dimension <<"D and 3D float (DVF)" << std::endl;
108     //  UpdateWithDimAndPixelType<Dimension, itk::Vector<float, 2> >();
109     //       }
110     
111     else if (Components==3)
112       {
113         if (m_Verbose) std::cout  << "Launching transform in "<< Dimension <<"D and 3D float (DVF)" << std::endl;
114         UpdateWithDimAndPixelType<Dimension, itk::Vector<float, 3> >();
115       }
116     
117     //     else if (Components==4)
118     //       {
119     //  if (m_Verbose) std::cout  << "Launching transform in "<< Dimension <<"D and 3D float (DVF)" << std::endl;
120     //  UpdateWithDimAndPixelType<Dimension, itk::Vector<float, 3> >();
121     //       }
122     else std::cerr<<"Number of components is "<<Components<<", not supported!"<<std::endl;
123   }
124
125
126   //-------------------------------------------------------------------
127   // Update with the number of dimensions and the pixeltype
128   //-------------------------------------------------------------------
129   template<class args_info_type>
130   template <unsigned int Dimension, class  PixelType> 
131   void 
132   AverageTemporalDimensionGenericFilter<args_info_type>::UpdateWithDimAndPixelType()
133   {
134
135     // ImageTypes
136     typedef itk::Image<PixelType, Dimension> InputImageType;
137     typedef itk::Image<PixelType, Dimension-1> OutputImageType;
138     typename InputImageType::Pointer input;
139
140     if (m_ArgsInfo.input_given ==1 )
141       {
142         // Read the input
143         typedef itk::ImageFileReader<InputImageType> InputReaderType;
144         typename InputReaderType::Pointer reader = InputReaderType::New();
145         reader->SetFileName( m_InputFileName);
146         reader->Update();
147         input= reader->GetOutput();
148       }
149
150     else 
151       {
152         // Read and join multiple inputs
153         if (m_Verbose) std::cout<<m_ArgsInfo.input_given<<" inputs given..."<<std::endl;
154         std::vector<std::string> filenames;
155         for(unsigned int i=0; i<m_ArgsInfo.input_given;i++)
156           {
157             if (m_Verbose) std::cout<<m_ArgsInfo.input_arg[i]<<std::endl;
158             filenames.push_back(m_ArgsInfo.input_arg[i]);
159           }
160         
161         typedef itk::ImageSeriesReader<InputImageType> ImageReaderType;
162         typename  ImageReaderType::Pointer reader= ImageReaderType::New();
163         reader->SetFileNames(filenames);
164         reader->Update();
165         input =reader->GetOutput();
166       }
167
168     
169     // Output properties
170     typename OutputImageType::RegionType region;
171     typename OutputImageType::RegionType::SizeType size;
172     typename OutputImageType::IndexType index;
173     typename OutputImageType::SpacingType spacing;
174     typename OutputImageType::PointType origin;
175     typename InputImageType::RegionType region4D=input->GetLargestPossibleRegion();
176     typename InputImageType::RegionType::SizeType size4D=region4D.GetSize();
177     typename InputImageType::IndexType index4D=region4D.GetIndex();
178     typename InputImageType::SpacingType spacing4D=input->GetSpacing();
179     typename InputImageType::PointType origin4D=input->GetOrigin();
180
181     for (unsigned int i=0; i< Dimension-1; i++)
182       {
183         size[i]=size4D[i];
184         index[i]=index4D[i];
185         spacing[i]=spacing4D[i];
186         origin[i]=origin4D[i];
187       }
188     region.SetSize(size);
189     region.SetIndex(index);
190     typename OutputImageType::Pointer output= OutputImageType::New();
191     output->SetRegions(region);
192     output->SetSpacing(spacing);
193     output->SetOrigin(origin);
194     output->Allocate();
195
196
197     // Region iterators
198     typedef itk::ImageRegionIterator<InputImageType> IteratorType;
199     std::vector<IteratorType> iterators(size4D[Dimension-1]);
200     for (unsigned int i=0; i< size4D[Dimension-1]; i++)
201       {
202         typename InputImageType::RegionType regionIt=region4D;
203         typename InputImageType::RegionType::SizeType sizeIt=regionIt.GetSize();
204         sizeIt[Dimension-1]=1;
205         regionIt.SetSize(sizeIt);
206         typename InputImageType::IndexType indexIt=regionIt.GetIndex();
207         indexIt[Dimension-1]=i;
208         regionIt.SetIndex(indexIt);
209         iterators[i]=IteratorType(input, regionIt);
210        }
211
212     typedef itk::ImageRegionIterator<OutputImageType> OutputIteratorType;
213     OutputIteratorType avIt(output, output->GetLargestPossibleRegion());
214
215     // Average
216     PixelType vector;
217     PixelType zeroVector=itk::NumericTraits<PixelType>::Zero;
218     //zeroVector.Fill(0.0);
219     while (!(iterators[0]).IsAtEnd())
220       {
221         vector=zeroVector;
222         for (unsigned int i=0; i<size4D[Dimension-1]; i++)
223           {
224             vector+=iterators[i].Get();
225             ++(iterators[i]);
226           }
227         vector/=size4D[Dimension-1];
228         avIt.Set(vector);
229         ++avIt;
230       }
231     
232     // Output
233     typedef itk::ImageFileWriter<OutputImageType> WriterType;
234     typename WriterType::Pointer writer = WriterType::New();
235     writer->SetFileName(m_ArgsInfo.output_arg);
236     writer->SetInput(output);
237     writer->Update();
238
239   }
240
241
242 }//end clitk
243  
244 #endif //#define clitkAverageTemporalDimensionGenericFilter_txx