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