]> Creatis software - clitk.git/blob - tools/clitkNVectorImageTo4DImageGenericFilter.txx
COMP: fix compilation errors for ITKv5
[clitk.git] / tools / clitkNVectorImageTo4DImageGenericFilter.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 clitkNVectorImageTo4DImageGenericFilter_txx
19 #define clitkNVectorImageTo4DImageGenericFilter_txx
20
21 /* =================================================
22  * @file   clitkNVectorImageTo4DImageGenericFilter.txx
23  * @author 
24  * @date   
25  * 
26  * @brief 
27  * 
28  ===================================================*/
29
30 #include "itkVectorImageToImageAdaptor.h"
31
32 namespace clitk
33 {
34
35   //-------------------------------------------------------------------
36   // Update with the number of dimensions
37   //-------------------------------------------------------------------
38   template<unsigned int Dimension>
39   void 
40   NVectorImageTo4DImageGenericFilter::UpdateWithDim(std::string PixelType, unsigned int Components)
41   {
42     if (m_Verbose) std::cout << "Image was detected to be "<<Dimension<<"D and "<< PixelType<<"..."<<std::endl;
43
44     if (PixelType == "short") {
45       UpdateWithDimAndPixelType<Dimension, short>();
46     } else if (PixelType == "unsigned short") {
47       UpdateWithDimAndPixelType<Dimension, unsigned short>();
48     } else if (PixelType == "unsigned_short") {
49       UpdateWithDimAndPixelType<Dimension, unsigned short>();
50     } else if (PixelType == "char") {
51       UpdateWithDimAndPixelType<Dimension, char>();
52     } else if (PixelType == "unsigned_char") {
53       UpdateWithDimAndPixelType<Dimension, unsigned char>();
54     } else if (PixelType == "int") {
55       UpdateWithDimAndPixelType<Dimension, int>();
56     } else if (PixelType == "unsigned_int") {
57       UpdateWithDimAndPixelType<Dimension, unsigned int>();
58     } else if (PixelType == "double") {
59       UpdateWithDimAndPixelType<Dimension, double>();
60     } else if (PixelType == "float") {
61       UpdateWithDimAndPixelType<Dimension, float>();
62     } else {
63       std::cerr << "Error, pixel type : \"" << PixelType << "\" unknown !" << std::endl;
64     }
65   }
66
67
68   //-------------------------------------------------------------------
69   // Update with the number of dimensions and the pixeltype
70   //-------------------------------------------------------------------
71   template <unsigned int Dimension, class  PixelType> 
72   void 
73   NVectorImageTo4DImageGenericFilter::UpdateWithDimAndPixelType()
74   {
75     // ImageTypes
76     typedef itk::VectorImage<PixelType, Dimension> InputImageType;
77     typedef itk::Image<PixelType, Dimension+1> OutputImageType;
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     //Filter
87     typedef itk::VectorImageToImageAdaptor<PixelType, Dimension> ImageAdaptorType;
88     typedef itk::ImageFileWriter<OutputImageType> WriterType;
89     typename ImageAdaptorType::Pointer adaptor = ImageAdaptorType::New();
90     typename OutputImageType::Pointer output = OutputImageType::New();
91     typename WriterType::Pointer writer = WriterType::New();
92
93     adaptor->SetExtractComponentIndex(0);
94     adaptor->SetImage(input);
95     std::string fileName=m_ArgsInfo.output_arg;
96     
97     //Create the output
98     typename OutputImageType::IndexType index;
99     index.Fill(0);
100     typename OutputImageType::SizeType size;
101     size.Fill(input->GetNumberOfComponentsPerPixel());
102     typename OutputImageType::SpacingType spacing;
103     spacing.Fill(1);
104     typename OutputImageType::PointType origin;
105     origin.Fill(0);
106     typename OutputImageType::DirectionType direction;
107     direction.SetIdentity();
108     for (unsigned int pixelDim=0; pixelDim<Dimension; ++pixelDim)
109     {
110       size[pixelDim]=adaptor->GetLargestPossibleRegion().GetSize(pixelDim);
111       spacing[pixelDim]=input->GetSpacing()[pixelDim];
112       origin[pixelDim]=input->GetOrigin()[pixelDim];
113       for (unsigned int pixelDim2=0; pixelDim2<Dimension; ++pixelDim2)
114       {
115         direction[pixelDim][pixelDim2]=input->GetDirection()[pixelDim][pixelDim2];
116       }
117     }
118     typename OutputImageType::RegionType region;
119     region.SetSize(size);
120     region.SetIndex(index);    
121     output->SetRegions(region);
122     output->SetOrigin(origin);
123     output->SetDirection(direction);
124     output->SetSpacing(spacing);
125     output->Allocate();
126     writer->SetInput(output);
127     
128     //Copy each channel
129     for (unsigned int pixelDim=0; pixelDim<input->GetNumberOfComponentsPerPixel(); ++pixelDim)
130     {
131       adaptor->SetExtractComponentIndex(pixelDim);
132       
133       itk::ImageRegionIterator<InputImageType> imageIterator(input,input->GetLargestPossibleRegion());
134
135       while(!imageIterator.IsAtEnd())
136       {
137         typename OutputImageType::IndexType indexVector;
138         indexVector.Fill(0);
139         for (unsigned int indexDim=0; indexDim<Dimension; ++indexDim)
140         {
141           indexVector[indexDim]=imageIterator.GetIndex().GetElement(indexDim);
142         }
143         indexVector[Dimension]=pixelDim;
144         
145         output->SetPixel(indexVector, adaptor->GetPixel(imageIterator.GetIndex()));
146         ++imageIterator;
147       }
148     }
149     // Output
150     writer->SetFileName(fileName);
151     writer->Update();
152   }
153 }//end clitk
154  
155 #endif //#define clitkNVectorImageTo4DImageGenericFilter_txx