]> Creatis software - clitk.git/blob - common/vvImageReader.txx
corrected error when opening 4D vector fields
[clitk.git] / common / vvImageReader.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
19 #ifndef VVIMAGEREADER_TXX
20 #define VVIMAGEREADER_TXX
21
22 #include <string>
23 #include <itkImageFileReader.h>
24 #include <itkImageSeriesReader.h>
25 #include <itkImageToVTKImageFilter.h>
26 #include <itkAnalyzeImageIO.h>
27 #include <itkFlexibleVectorCastImageFilter.h>
28
29 #include <vtkTransform.h>
30
31 #include "clitkCommon.h"
32 #include "clitkConfiguration.h"
33 #include "vvFromITK.h"
34 //----------------------------------------------------------------------------
35 template<unsigned int VImageDimension>
36 void vvImageReader::UpdateWithDim(std::string InputPixelType)
37 {
38   if (mType == VECTORFIELD || mType == VECTORFIELDWITHTIME)
39   {
40     if (VImageDimension == 4)
41       UpdateWithDimAndInputVectorPixelType<itk::Vector<float,3>,VImageDimension>();
42     else
43       UpdateWithDimAndInputVectorPixelType<itk::Vector<float,VImageDimension>,VImageDimension>();
44   }
45   else if (InputPixelType == "short")
46     UpdateWithDimAndInputPixelType<short,VImageDimension>();
47   else if (InputPixelType == "unsigned_short")
48     UpdateWithDimAndInputPixelType<unsigned short,VImageDimension>();
49   else if (InputPixelType == "char")
50     UpdateWithDimAndInputPixelType<char,VImageDimension>();
51   else if (InputPixelType == "unsigned_char")
52     UpdateWithDimAndInputPixelType<unsigned char,VImageDimension>();
53   else if (InputPixelType == "int")
54     UpdateWithDimAndInputPixelType<int,VImageDimension>();
55   else if (InputPixelType == "unsigned_int")
56     UpdateWithDimAndInputPixelType<unsigned int,VImageDimension>();
57   else if (InputPixelType == "double")
58     UpdateWithDimAndInputPixelType<double,VImageDimension>();
59   else if (InputPixelType == "float")
60     UpdateWithDimAndInputPixelType<float,VImageDimension>();
61   else
62     std::cerr << "Error, input pixel type : " << InputPixelType << " unknown !" << std::endl;
63
64   if (CLITK_EXPERIMENTAL && mLastError.size()==0) {
65     //ReadNkiImageTransform();
66     ReadMatImageTransform();
67   }
68 }
69 //----------------------------------------------------------------------------
70
71
72 //----------------------------------------------------------------------------
73 template<class InputPixelType, unsigned int VImageDimension>
74 void vvImageReader::UpdateWithDimAndInputPixelType()
75 {
76   itk::AnalyzeImageIO *analyzeImageIO = NULL;
77
78   if (mType == MERGEDWITHTIME)   // In this case we can load the images
79     // one at the time to avoid excessive
80     // memory use
81   {
82     mImage=vvImage::New();
83
84     for (std::vector<std::string>::const_iterator i=mInputFilenames.begin(); i!=mInputFilenames.end(); i++) {
85       typedef itk::Image< InputPixelType, VImageDimension-1 > InputImageType;
86       typedef itk::ImageFileReader<InputImageType> ReaderType;
87       typename ReaderType::Pointer reader = ReaderType::New();
88       reader->ReleaseDataFlagOn();
89       reader->SetFileName(*i);
90       try {
91         mImage->AddItkImage<InputImageType>(reader->GetOutput());
92       } catch ( itk::ExceptionObject & err ) {
93         std::cerr << "Error while reading " << mInputFilenames[0].c_str()
94                   << " " << err << std::endl;
95         std::stringstream error;
96         error << err;
97         mLastError = error.str();
98         return;
99       }
100       analyzeImageIO = dynamic_cast<itk::AnalyzeImageIO*>( reader->GetImageIO() );
101     }
102   } else if (mType == SLICED) {
103     mImage=vvImage::New();
104     typedef itk::Image< InputPixelType, VImageDimension > InputImageType;
105     typedef itk::ImageFileReader<InputImageType> ReaderType;
106     typename ReaderType::Pointer reader = ReaderType::New();
107     reader->SetFileName(mInputFilenames[0]);
108     reader->UpdateOutputInformation();
109
110     typedef itk::Image< InputPixelType, VImageDimension-1 > SlicedImageType;
111     typedef itk::ExtractImageFilter<InputImageType, SlicedImageType> FilterType;
112
113     typename InputImageType::RegionType inputRegion = reader->GetOutput()->GetLargestPossibleRegion();
114     typename InputImageType::SizeType inputSize = inputRegion.GetSize();
115     typename InputImageType::IndexType start = inputRegion.GetIndex();
116     typename InputImageType::SizeType extractedRegionSize = inputSize;
117     typename InputImageType::RegionType extractedRegion;
118     extractedRegionSize[VImageDimension - 1] = 0;
119     extractedRegion.SetSize(extractedRegionSize);
120     start[VImageDimension - 1] = mSlice;
121     extractedRegion.SetIndex(start);
122
123     typename FilterType::Pointer filter = FilterType::New();
124     filter->SetExtractionRegion(extractedRegion);
125     filter->SetInput(reader->GetOutput());
126     filter->ReleaseDataFlagOn();
127 #if ITK_VERSION_MAJOR == 4
128     filter->SetDirectionCollapseToSubmatrix();
129 #endif
130     try {
131       mImage->AddItkImage<SlicedImageType>(filter->GetOutput());
132     } catch ( itk::ExceptionObject & err ) {
133       std::cerr << "Error while slicing " << mInputFilenames[0].c_str()
134                 << "(slice #" << mSlice << ") " << err << std::endl;
135       return;
136     }
137     analyzeImageIO = dynamic_cast<itk::AnalyzeImageIO*>( reader->GetImageIO() );
138   } else {
139     if (mInputFilenames.size() > 1) {
140       typedef itk::Image< InputPixelType, VImageDimension > InputImageType;
141       typedef itk::ImageSeriesReader<InputImageType> ReaderType;
142       typename ReaderType::Pointer reader = ReaderType::New();
143       reader->SetFileNames(mInputFilenames);
144       reader->ReleaseDataFlagOn();
145
146       try {
147         if (mType == IMAGEWITHTIME)
148         {
149           std::cerr << "We should never come here:" << std::endl
150             << "  Calling vvImageReader with multiple images and IMAGEWITHTIME is undefined." << std::endl
151             << "  You are probably looking for MERGEDWITHTIME Type." << std::endl;
152           return;
153         }
154         else
155           mImage=vvImageFromITK<VImageDimension,InputPixelType>(reader->GetOutput());
156       } catch ( itk::ExceptionObject & err ) {
157         std::cerr << "Error while reading image series:" << err << std::endl;
158         std::stringstream error;
159         error << err;
160         mLastError = error.str();
161         return;
162       }
163     } else {
164       typedef itk::Image< InputPixelType, VImageDimension > InputImageType;
165       typedef itk::ImageFileReader<InputImageType> ReaderType;
166       typename ReaderType::Pointer reader = ReaderType::New();
167       reader->SetFileName(mInputFilenames[0]);
168       reader->ReleaseDataFlagOn();
169
170       try {
171         mImage = vvImageFromITK<VImageDimension,InputPixelType>(reader->GetOutput(), mType == IMAGEWITHTIME || mType == VECTORFIELDWITHTIME);
172       } catch ( itk::ExceptionObject & err ) {
173         std::cerr << "Error while reading " << mInputFilenames[0].c_str()
174                   << " " << err << std::endl;
175         std::stringstream error;
176         error << err;
177         mLastError = error.str();
178         return;
179       }
180       analyzeImageIO = dynamic_cast<itk::AnalyzeImageIO*>( reader->GetImageIO() );
181     }
182   }
183
184   // For unknown analyze orientations, we set identity
185   if(analyzeImageIO) {
186     const double m[16] = {1.,0.,0.,0.,
187                           0.,0.,1.,0.,
188                           0.,-1.,0.,0.,
189                           0.,0.,0.,1.};
190     int i;
191     for(i=0; i<16 && m[i]==mImage->GetTransform()->GetMatrix()->GetElement(i%4, i/4); i++);
192     if(i==16) {
193       itkWarningMacro(<< "Analyze image file format detected with unknown orientation. "
194                       << "Forcing identity orientation, use other file format if not ok.");
195       mImage->GetTransform()->Identity();
196     }
197   }
198 }
199 //----------------------------------------------------------------------------
200
201 //----------------------------------------------------------------------------
202 template<class InputPixelType, unsigned int VImageDimension>
203 void vvImageReader::UpdateWithDimAndInputVectorPixelType()
204 {
205   itk::AnalyzeImageIO *analyzeImageIO = NULL;
206
207   typedef itk::Image< InputPixelType, VImageDimension > InputImageType;
208   typename InputImageType::Pointer input;
209
210   if (mInputFilenames.size() > 1) {
211     typedef itk::ImageSeriesReader<InputImageType> ReaderType;
212     typename ReaderType::Pointer reader = ReaderType::New();
213     reader->SetFileNames(mInputFilenames);
214     reader->ReleaseDataFlagOn();
215     try {
216       reader->Update();
217       input = reader->GetOutput();
218     } catch ( itk::ExceptionObject & err ) {
219       std::cerr << "Error while reading image series:" << err << std::endl;
220       std::stringstream error;
221       error << err;
222       mLastError = error.str();
223       return;
224     }
225   } else {
226     typedef itk::ImageFileReader<InputImageType> ReaderType;
227     typename ReaderType::Pointer reader = ReaderType::New();
228     reader->SetFileName(mInputFilenames[0]);
229     reader->ReleaseDataFlagOn();
230     try {
231       reader->Update();
232       input = reader->GetOutput();
233     } catch ( itk::ExceptionObject & err ) {
234       std::cerr << "Error while reading " << mInputFilenames[0].c_str()
235         << " " << err << std::endl;
236       std::stringstream error;
237       error << err;
238       mLastError = error.str();
239       return;
240     }
241     analyzeImageIO = dynamic_cast<itk::AnalyzeImageIO*>( reader->GetImageIO() );
242   }
243   
244   typedef itk::Image< itk::Vector<float , 3>, VImageDimension > VectorImageType;
245   typedef itk::FlexibleVectorCastImageFilter<InputImageType, VectorImageType> CasterType;
246   typename VectorImageType::Pointer casted_input;
247   typename CasterType::Pointer caster = CasterType::New();
248   caster->SetInput(input);
249   caster->Update();
250   casted_input = caster->GetOutput();
251   
252   mImage = vvImageFromITK<VImageDimension, itk::Vector<float, 3> >(casted_input, mType == IMAGEWITHTIME || mType == VECTORFIELDWITHTIME);
253
254   // For unknown analyze orientations, we set identity
255   if (analyzeImageIO)
256   {
257     const double m[16] = {1.,0.,0.,0.,
258                           0.,0.,1.,0.,
259                           0.,-1.,0.,0.,
260                           0.,0.,0.,1.};
261     int i;
262     for (i = 0; i < 16 && m[i] == mImage->GetTransform()->GetMatrix()->GetElement(i % 4, i / 4); i++)
263       ;
264     if (i == 16)
265     {
266       itkWarningMacro(<< "Analyze image file format detected with unknown orientation. "
267                       << "Forcing identity orientation, use other file format if not ok.");
268       mImage->GetTransform()->Identity();
269     }
270   }
271 }
272 //----------------------------------------------------------------------------
273
274 #endif /* end #define vvImageReader_TXX */
275