]> Creatis software - gdcm.git/blob - vtk/vtkgdcmSerieViewer.cxx
Allow checking reverse order file sorting
[gdcm.git] / vtk / vtkgdcmSerieViewer.cxx
1 /*=========================================================================
2                                                                                 
3   Program:   gdcm
4   Module:    $RCSfile: vtkgdcmSerieViewer.cxx,v $
5   Language:  C++
6   Date:      $Date: 2005/07/29 15:24:04 $
7   Version:   $Revision: 1.3 $
8                                                                                 
9   Copyright (c) CREATIS (Centre de Recherche et d'Applications en Traitement de
10   l'Image). All rights reserved. See Doc/License.txt or
11   http://www.creatis.insa-lyon.fr/Public/Gdcm/License.html for details.
12                                                                                 
13      This software is distributed WITHOUT ANY WARRANTY; without even
14      the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
15      PURPOSE.  See the above copyright notices for more information.
16                                                                                 
17 =========================================================================*/
18 // This example illustrates how the vtkGdcmReader vtk class can be
19 // used in order to:
20 //  * produce a simple (vtk based) Dicom image STACK VIEWER.
21 //  * dump the stack considered as a volume in a vtkStructuredPoints
22 //    vtk file: the vtk gdcm wrappers can be seen as a simple way to convert
23 //    a stack of Dicom images into a native vtk volume.
24 //
25 // Usage:
26 //  * the filenames of the Dicom images constituting the stack should be
27 //    given as command line arguments,
28 //  * you can navigate through the stack by hitting any character key,
29 //  * the produced vtk file is named "foo.vtk" (in the invocation directory).
30 // 
31 //----------------------------------------------------------------------------
32 #include <vtkRenderWindowInteractor.h>
33 #include <vtkImageViewer.h>
34 #include <vtkStructuredPoints.h>
35 #include <vtkStructuredPointsWriter.h>
36 #include <vtkCommand.h>
37 #include <vtkRenderer.h>
38 #include <vtkImageMapToColors.h>
39 #include <vtkLookupTable.h>
40
41 #include "vtkGdcmReader.h"
42 #include "gdcmDocument.h"  // for NO_SHADOWSEQ
43 #include "gdcmSerieHelper.h"
44 #include "gdcmDebug.h"
45 #ifndef vtkFloatingPointType
46 #define vtkFloatingPointType float
47 #endif
48
49 //----------------------------------------------------------------------------
50 // Callback for the interaction
51 class vtkgdcmObserver : public vtkCommand
52 {
53 public:
54    virtual char const *GetClassName() const 
55    { 
56       return "vtkgdcmObserver";
57    }
58    static vtkgdcmObserver *New() 
59    { 
60       return new vtkgdcmObserver; 
61    }
62    vtkgdcmObserver()
63    {
64       this->ImageViewer = NULL;
65    }
66    virtual void Execute(vtkObject *, unsigned long event, void* )
67    {
68       if ( this->ImageViewer )
69       {
70          if ( event == vtkCommand::CharEvent )
71          {
72             int max = ImageViewer->GetWholeZMax();
73             int slice = (ImageViewer->GetZSlice() + 1 ) % ++max;
74             ImageViewer->SetZSlice( slice );
75             ImageViewer->GetRenderer()->ResetCameraClippingRange();
76             ImageViewer->Render();
77          }
78       }
79    }
80    vtkImageViewer *ImageViewer;
81 };
82
83
84 int main(int argc, char *argv[])
85 {
86    if( argc == 1 )
87    {
88       std::cout << "Usage : vtkgdcmSerieViewer directoryName [reverse] [debug]"
89                 << std::endl;
90       return 0;
91    }
92
93   if( argc > 3 )
94     gdcm::Debug::DebugOn();
95   
96    vtkGdcmReader *reader = vtkGdcmReader::New();
97    reader->AllowLookupTableOff();
98
99    // ------------ to check Coherent File List as a parameter
100
101    gdcm::SerieHelper *sh = new gdcm::SerieHelper();
102    sh->SetLoadMode(NO_SHADOWSEQ);
103    if (argc > 2)
104       sh->SetSortOrderToReverse();  // just to check
105    sh->SetDirectory( argv[1], true);
106     
107    // Just to see
108
109    int nbFiles;
110    // For all the Coherent Files lists of the gdcm::Serie
111    gdcm::FileList *l = sh->GetFirstCoherentFileList();
112    if (l == 0 )
113    {
114       std::cout << "Oops! No CoherentFileList found ?!?" << std::endl;
115       return 0;
116    }
117    while (l)
118    { 
119       nbFiles = l->size() ;
120       if ( l->size() > 1 )
121       {
122          std::cout << "Sort list : " << nbFiles << " long" << std::endl;
123          sh->OrderFileList(l);  // sort the list
124          break;  // The first one is OK. user will have to check
125       }
126       else
127       {
128          std::cout << "Oops! Empty CoherentFileList found ?!?" << std::endl;
129       }
130       l = sh->GetNextCoherentFileList();
131    }
132
133    // Only the first FileList is dealt with (just an example)
134    // (The files will not be parsed twice by the reader)
135
136    //---------------------------------------------------------
137    reader->SetCoherentFileList(l);
138    //---------------------------------------------------------
139
140
141 // TODO : allow user to choose Load Mode
142
143  //  reader->SetLoadMode(NO_SHADOWSEQ);  
144    reader->Update();
145
146    //print debug info:
147    reader->GetOutput()->Print( cout );
148
149    vtkRenderWindowInteractor *iren = vtkRenderWindowInteractor::New();
150
151    vtkImageViewer *viewer = vtkImageViewer::New();
152
153    if( reader->GetLookupTable() )
154    {
155       //convert to color:
156       vtkImageMapToColors *map = vtkImageMapToColors::New ();
157       map->SetInput (reader->GetOutput());
158       map->SetLookupTable (reader->GetLookupTable());
159       map->SetOutputFormatToRGB();
160       viewer->SetInput ( map->GetOutput() );
161       map->Delete();
162    }
163    else
164    {
165       vtkFloatingPointType *range = reader->GetOutput()->GetScalarRange();
166       viewer->SetColorLevel (0.5 * (range[1] + range[0]));
167       viewer->SetColorWindow (range[1] - range[0]);
168
169       viewer->SetInput ( reader->GetOutput() );
170    }
171    viewer->SetupInteractor (iren);
172   
173    //vtkFloatingPointType *range = reader->GetOutput()->GetScalarRange();
174    //viewer->SetColorWindow (range[1] - range[0]);
175    //viewer->SetColorLevel (0.5 * (range[1] + range[0]));
176
177    // Here is where we setup the observer, 
178    vtkgdcmObserver *obs = vtkgdcmObserver::New();
179    obs->ImageViewer = viewer;
180    iren->AddObserver(vtkCommand::CharEvent,obs);
181    obs->Delete();
182
183    //viewer->Render();
184    iren->Initialize();
185    iren->Start();
186
187    //if you wish you can export dicom to a vtk file  
188    vtkStructuredPointsWriter *writer = vtkStructuredPointsWriter::New();
189    writer->SetInput( reader->GetOutput());
190    writer->SetFileName( "foo.vtk" );
191    writer->SetFileTypeToBinary();
192    //writer->Write();
193
194    reader->Delete();
195    iren->Delete();
196    viewer->Delete();
197    writer->Delete();
198
199    return 0;
200 }