]> Creatis software - clitk.git/blob - tools/clitkDicom2Image.cxx
Change dicom file sorting procedure in clitkDicom2Image
[clitk.git] / tools / clitkDicom2Image.cxx
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 // clitk includes
20 #include "clitkIO.h"
21 #include "clitkDicom2Image_ggo.h"
22 #include "clitkCommon.h"
23 #include "clitkImageCommon.h"
24 #include "vvImageReader.h"
25 #include "vvImageWriter.h"
26 #include <itkGDCMImageIO.h>
27 #include <itkGDCMSeriesFileNames.h>
28 #include <gdcmFile.h>
29 #include <vtkVersion.h>
30 #include <vtkImageChangeInformation.h>
31 #if GDCM_MAJOR_VERSION >= 2
32   #include <gdcmImageHelper.h>
33   #include <gdcmAttribute.h>
34   #include <gdcmReader.h>
35 #endif
36
37 #include <set>
38
39 //====================================================================
40 int main(int argc, char * argv[])
41 {
42   // init command line
43   GGO(clitkDicom2Image, args_info);
44   CLITK_INIT;
45
46   std::vector<std::string> input_files;
47   ///if std_input is given, read the input files from stdin
48   if (args_info.std_input_given) {
49     while (true) {
50       std::string tmp;
51       std::cin >> tmp;
52       if(std::cin.good())
53         input_files.push_back(tmp);
54       else break;
55     }
56     args_info.inputs_num=input_files.size();
57   } else for (unsigned int i=0; i<args_info.inputs_num; i++)
58       input_files.push_back(args_info.inputs[i]);
59
60   //Get GDCMSeriesFileNames order to sort filenames
61   typedef itk::GDCMSeriesFileNames NamesGeneratorType;
62   NamesGeneratorType::Pointer nameGenerator = NamesGeneratorType::New();
63   nameGenerator->SetUseSeriesDetails(true);
64   std::string folderName=".";
65   const size_t last_slash_idx = input_files[0].rfind('/');
66   if (std::string::npos != last_slash_idx)
67     folderName = input_files[0].substr(0, last_slash_idx);
68   nameGenerator->SetInputDirectory(folderName);
69
70   //===========================================
71   /// Get slices locations ...
72   std::string series_UID = "";
73   std::set<std::string> series_UIDs;
74   std::map< std::string, std::vector<double> > theorigin;
75   std::map< std::string, std::vector<double> > theorientation;
76   std::map< std::string, std::vector<double> > sliceLocations;
77   std::map< std::string, std::vector<std::string> > seriesFiles;
78 #if GDCM_MAJOR_VERSION >= 2
79   if (args_info.verbose_flag)
80     std::cout << "Using GDCM-2.x" << std::endl;
81 #else
82   if (args_info.verbose_flag) {
83     std::cout << "Not using GDCM-2.x" << std::endl;
84     std::cout<< "The image orientation is not supported with this version of GDCM" <<std::endl;
85   }
86 #endif
87   for(unsigned int i=0; i<args_info.inputs_num; i++) {
88     if (args_info.verbose_flag)
89         std::cout << "Reading <" << input_files[i] << std::endl;
90 #if GDCM_MAJOR_VERSION >= 2
91     gdcm::Reader hreader;
92     hreader.SetFileName(input_files[i].c_str());
93     hreader.Read();
94     gdcm::DataSet& ds = hreader.GetFile().GetDataSet();
95
96     if (args_info.extract_series_flag) {
97       gdcm::Attribute<0x20,0x000e> series_UID_att;
98       series_UID_att.SetFromDataSet(ds);
99       series_UID = series_UID_att.GetValue();
100     }
101
102     series_UIDs.insert(series_UID);
103     theorigin[series_UID] = gdcm::ImageHelper::GetOriginValue(hreader.GetFile());
104     theorientation[series_UID] = gdcm::ImageHelper::GetDirectionCosinesValue(hreader.GetFile());
105     if (args_info.patientSystem_flag) {
106       double n1 = theorientation[series_UID][1]*theorientation[series_UID][5]-
107                   theorientation[series_UID][2]*theorientation[series_UID][4];
108       double n2 = theorientation[series_UID][3]*theorientation[series_UID][2]-
109                   theorientation[series_UID][5]*theorientation[series_UID][0];
110       double n3 = theorientation[series_UID][0]*theorientation[series_UID][4]-
111                   theorientation[series_UID][1]*theorientation[series_UID][3];
112       double sloc = theorigin[series_UID][0]*n1+
113                     theorigin[series_UID][1]*n2+
114                     theorigin[series_UID][2]*n3;
115       sliceLocations[series_UID].push_back(sloc);
116     } else
117       sliceLocations[series_UID].push_back(theorigin[series_UID][2]);
118     seriesFiles[series_UID].push_back(input_files[i]);
119
120     gdcm::Attribute<0x28, 0x100> pixel_size;
121     pixel_size.SetFromDataSet(ds);
122     /* if (pixel_size.GetValue() != 16)
123        {
124        std::cerr << "Pixel type not 2 bytes ! " << std::endl;
125        std::cerr << "In file " << input_files[i] << std::endl;
126        exit(0);
127        }
128     */
129 #else
130   gdcm::File *header = new gdcm::File();
131   header->SetFileName(input_files[i]);
132   header->SetMaxSizeLoadEntry(16384); // required ?
133   header->Load();
134
135   if (args_info.extract_series_flag) {
136     series_UID = header->GetEntryValue(0x20,0x000e).c_str();
137   }
138
139   series_UIDs.insert(series_UID);
140   theorigin[series_UID].resize(3);
141   theorigin[series_UID][0] = header->GetXOrigin();
142   theorigin[series_UID][1] = header->GetYOrigin();
143   theorigin[series_UID][2] = header->GetZOrigin();
144   sliceLocations[series_UID].push_back(theorigin[series_UID][2]);
145   seriesFiles[series_UID].push_back(input_files[i]);
146   /*if (header->GetPixelSize() != 2) {
147     std::cerr << "Pixel type 2 bytes ! " << std::endl;
148     std::cerr << "In file " << input_files[i] << std::endl;
149     exit(0);
150   }
151   */
152 #endif
153   }
154
155   //===========================================
156   // Sort slices locations ...
157   std::set<std::string>::iterator sn = series_UIDs.begin();
158   while ( sn != series_UIDs.end() ) {
159     std::vector<double> locs = sliceLocations[*sn];
160     std::vector<double> origin = theorigin[*sn];
161     std::vector<std::string> files = seriesFiles[*sn];
162     std::vector<int> sliceIndex;
163     //clitk::GetSortedIndex(locs, sliceIndex);
164     //Look for files into GDCMSeriesFileNames, because it sorts files correctly and take the order
165     const std::vector<std::string> & temp = nameGenerator->GetFileNames(*sn);
166     for(unsigned int i=0; i<files.size(); i++) {
167       int j(0);
168       while (sliceIndex.size()==i && j<temp.size()) {
169         const size_t last_slash_idx2 = temp[j].rfind('/');
170         std::string tempFilename(temp[j]);
171         if (std::string::npos != last_slash_idx2)
172           tempFilename = temp[j].substr(last_slash_idx2+1, temp[j].size()-1);
173         if (tempFilename == files[i])
174           sliceIndex.push_back(j);
175         ++j;
176       }
177     }
178
179     if (args_info.verbose_flag) {
180       std::cout << locs[sliceIndex[0]] << " -> "
181                 << sliceIndex[0] << " / " << 0 << " => "
182                 << "0 mm "
183                 << files[sliceIndex[0]]
184                 << std::endl;
185       for(unsigned int i=1; i<sliceIndex.size(); i++) {
186         std::cout << locs[sliceIndex[i]] << " -> "
187                   << sliceIndex[i] << " / " << i << " => "
188                   << locs[sliceIndex[i]] - locs[sliceIndex[i-1]]
189                   << "mm "
190                   << files[sliceIndex[i]]
191                   << std::endl;
192       }
193     }
194
195     //===========================================
196     // Analyze slices locations ...
197     double currentDist;
198     double dist=0;
199     double tolerance = args_info.tolerance_arg;
200     double previous = locs[sliceIndex[0]];
201     for(unsigned int i=1; i<sliceIndex.size(); i++) {
202       currentDist = locs[sliceIndex[i]]-previous;
203       if (i!=1) {
204         if (fabs(dist-currentDist) > tolerance) {
205           std::cout << "ERROR : " << std::endl
206                     << "Current slice pos is  = " << locs[sliceIndex[i]] << std::endl
207                     << "Previous slice pos is = " << previous << std::endl
208                     << "Current file is       = " << files[sliceIndex[i]] << std::endl
209                     << "Current index is      = " << i << std::endl
210                     << "Current sortindex is  = " << sliceIndex[i] << std::endl
211                     << "Current slice diff    = " << dist << std::endl
212                     << "Current error         = " << fabs(dist-currentDist) << std::endl;
213           exit(1);
214         }
215       } else dist = currentDist;
216       previous = locs[sliceIndex[i]];
217     }
218
219     //===========================================
220     // Create ordered vector of filenames
221     std::vector<std::string> sorted_files;
222     sorted_files.resize(sliceIndex.size());
223     for(unsigned int i=0; i<sliceIndex.size(); i++)
224       sorted_files[i] = files[ sliceIndex[i] ];
225
226     //===========================================
227     // Read write serie
228     vvImageReader::Pointer reader = vvImageReader::New();
229     reader->SetInputFilenames(sorted_files);
230     reader->SetPatientCoordinateSystem(args_info.patientSystem_flag);
231     reader->Update(vvImageReader::DICOM);
232     if (reader->GetLastError().size() != 0) {
233       std::cerr << reader->GetLastError() << std::endl;
234       return 1;
235     }
236
237     vvImage::Pointer image = reader->GetOutput();
238     vtkImageData* vtk_image = image->GetFirstVTKImageData();
239     vtkImageChangeInformation* modifier = vtkImageChangeInformation::New();
240     if  (args_info.focal_origin_given) {
241       std::vector<double> spacing = image->GetSpacing();
242       std::vector<int> size = image->GetSize();
243       origin[0] = -spacing[0]*size[0]/2.0;
244       origin[1] = -spacing[1]*size[1]/2.0;
245 #if VTK_MAJOR_VERSION <= 5
246       modifier->SetInput(vtk_image);
247 #else
248       modifier->SetInputData(vtk_image);
249 #endif
250       modifier->SetOutputOrigin(origin[0], origin[1], locs[sliceIndex[0]]);
251       modifier->Update();
252       vvImage::Pointer focal_image = vvImage::New();
253       focal_image->AddVtkImage(modifier->GetOutput());
254       image = focal_image;
255     }
256
257     std::string outfile;
258     if (series_UIDs.size() == 1)
259       outfile = args_info.output_arg;
260     else {
261       std::ostringstream name;
262       std::vector<std::string> directory = clitk::SplitFilename(args_info.output_arg);
263       if (directory.size() == 2)
264         name << directory[0] << "/" << *sn << "_" << directory[1];
265       else
266         name << *sn << "_" << args_info.output_arg;
267       outfile = name.str();
268     }
269     vvImageWriter::Pointer writer = vvImageWriter::New();
270     writer->SetInput(image);
271     if (args_info.patientSystem_flag && !image->GetTransform().empty())
272       writer->SetSaveTransform(true);
273     writer->SetOutputFileName(outfile);
274     writer->Update();
275
276     modifier->Delete();
277
278     sn++;
279   }
280
281   return 0;
282 }