]> Creatis software - clitk.git/blob - tools/clitkGammaIndex.cxx
clitkGammaIndex switched to vvImage
[clitk.git] / tools / clitkGammaIndex.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 #include <vtkPoints.h>
19 #include <vtkCellArray.h>
20 #include <vtkPointData.h>
21 #include <vtkImageData.h>
22 #include <vtkMetaImageReader.h>
23 #include <vtkMetaImageWriter.h>
24 #include <vtkPNGReader.h>
25 #include <vtkPolyData.h>
26 #include <vtkCellLocator.h>
27 #include <iostream>
28 #include <cmath>
29 #include <cassert>
30 using std::endl;
31 using std::cout;
32
33 #include "clitkGammaIndex_ggo.h"
34
35 #include <vvImage.h>
36 #include <vvImageReader.h>
37
38 void saveImage(vtkImageData *image,const std::string &filename) {
39     cout << "saving " << filename << endl;
40     vtkImageWriter *writer = vtkMetaImageWriter::New();
41     writer->SetFileName(filename.c_str());
42     writer->SetInput(image);
43     writer->Write();
44     writer->Delete();
45 }
46
47
48 void insertTriangles(vtkCellArray *cells, vtkPoints *points, const vtkIdType ids[4]) {
49     double p0[3]; points->GetPoint(ids[0],p0);
50     double p1[3]; points->GetPoint(ids[1],p1);
51     double p2[3]; points->GetPoint(ids[2],p2);
52     double p3[3]; points->GetPoint(ids[3],p3);
53     //cout << "----------------------------------" << endl;
54     //cout << "p0=[" << p0[0] << "," << p0[1] << "," << p0[2] << "]" << endl;
55     //cout << "p1=[" << p1[0] << "," << p1[1] << "," << p1[2] << "]" << endl;
56     //cout << "p2=[" << p2[0] << "," << p2[1] << "," << p2[2] << "]" << endl;
57     //cout << "p3=[" << p3[0] << "," << p3[1] << "," << p3[2] << "]" << endl;
58
59     double center[] = {0,0,0};
60     for (int kk=0; kk<3; kk++) {
61         center[kk] += p0[kk];
62         center[kk] += p1[kk];
63         center[kk] += p2[kk];
64         center[kk] += p3[kk];
65         center[kk] /= 4.;
66     }
67
68     vtkIdType center_id = points->InsertNextPoint(center);
69     //cout << "center=[" << center[0] << "," << center[1] << "," << center[2] << "]" << endl;
70
71     cells->InsertNextCell(3);
72     cells->InsertCellPoint(ids[0]);
73     cells->InsertCellPoint(ids[1]);
74     cells->InsertCellPoint(center_id);
75
76     cells->InsertNextCell(3);
77     cells->InsertCellPoint(ids[1]);
78     cells->InsertCellPoint(ids[3]);
79     cells->InsertCellPoint(center_id);
80
81     cells->InsertNextCell(3);
82     cells->InsertCellPoint(ids[3]);
83     cells->InsertCellPoint(ids[2]);
84     cells->InsertCellPoint(center_id);
85
86     cells->InsertNextCell(3);
87     cells->InsertCellPoint(ids[2]);
88     cells->InsertCellPoint(ids[0]);
89     cells->InsertCellPoint(center_id);
90 }
91
92 double getMaximum(vtkImageData *image) {
93     bool first = true;
94     double maximum = 0;
95
96     for (int kk=0; kk<image->GetNumberOfPoints(); kk++) {
97         double value = image->GetPointData()->GetScalars()->GetTuple1(kk);
98
99         if (first) {
100             maximum = value;
101             first = false;
102             continue;
103         }
104
105         if (maximum<value) maximum = value;
106     }
107
108     return maximum;
109 }
110
111
112 vtkPolyData *buildPlane(vtkImageData *image,double spatial_margin,double dose_margin) {
113     vtkPoints *points = vtkPoints::New();
114     for (int kk=0; kk<image->GetNumberOfPoints(); kk++) {
115         double *point = image->GetPoint(kk);
116         double value = image->GetPointData()->GetScalars()->GetTuple1(kk);
117         assert(value>=0);
118         assert(point[2]==0);
119         point[2] = value;
120
121         point[0] /= spatial_margin;
122         point[1] /= spatial_margin;
123         point[2] /= dose_margin;
124
125 #ifndef NDEBUG
126         vtkIdType point_id = points->InsertNextPoint(point);
127         assert(kk==point_id);
128 #else
129         points->InsertNextPoint(point);
130 #endif
131     }
132
133     vtkCellArray *cells = vtkCellArray::New();
134     for (int kk=0; kk<image->GetNumberOfCells(); kk++) {
135         vtkCell *cell = image->GetCell(kk);
136         assert(cell->GetNumberOfPoints()==4);
137         insertTriangles(cells,points,cell->GetPointIds()->GetPointer(0));
138     }
139
140     vtkPolyData *data = vtkPolyData::New();
141     data->SetPoints(points);
142     data->SetPolys(cells);
143     points->Delete();
144     cells->Delete();
145
146     return data;
147 }
148
149 int main(int argc,char * argv[])
150 {
151     args_info_clitkGammaIndex args_info;
152
153     if (cmdline_parser_clitkGammaIndex(argc, argv, &args_info) != 0)
154         exit(1);
155
156     if (!args_info.absolute_dose_margin_given && !args_info.relative_dose_margin_given) {
157         std::cerr << "Specify either relative or absolute dose margin" << endl;
158         exit(1);
159     }
160
161     bool verbose = args_info.verbose_flag;
162
163     std::string reference_filename(args_info.reference_arg);
164     std::string target_filename(args_info.target_arg);
165     std::string gamma_filename(args_info.output_arg);
166     double space_margin = args_info.spatial_margin_arg;
167     double dose_rel_margin = args_info.relative_dose_margin_arg;
168     double dose_margin = args_info.absolute_dose_margin_arg;
169     bool use_dose_margin = args_info.absolute_dose_margin_given;
170
171     if (verbose) {
172         cout << "reference_filename=" << reference_filename << endl;
173         cout << "target_filename=" << target_filename << endl;
174         cout << "gamma_filename=" << gamma_filename << endl;
175         cout << "space_margin=" << space_margin << endl;
176         if (use_dose_margin) cout << "dose_margin=" << dose_margin << endl;
177         else cout << "dose_rel_margin=" << dose_rel_margin << endl;
178     }
179
180     // load reference
181     vtkImageData* reference = NULL;
182     {
183         vvImageReader::Pointer reader = vvImageReader::New();
184         reader->SetInputFilename(reference_filename);
185         reader->Update();
186         vvImage::Pointer vvimage = reader->GetOutput();
187         reference = vvimage->GetFirstVTKImageData();
188     }
189
190     // intensity normalisation
191     if (!use_dose_margin) {
192         dose_margin = getMaximum(reference)*dose_rel_margin;
193         if (verbose) cout << "dose_margin=" << dose_margin << endl;
194     }
195
196     // build surface
197     vtkPolyData *data = buildPlane(reference,space_margin,dose_margin);
198     reference->Delete();
199
200     vtkAbstractCellLocator *locator = vtkCellLocator::New();
201     locator->SetDataSet(data);
202     data->Delete();
203     locator->CacheCellBoundsOn();
204     locator->AutomaticOn();
205     locator->BuildLocator();
206
207     // load target
208     vtkImageData *target = NULL;
209     {
210         vvImageReader::Pointer reader = vvImageReader::New();
211         reader->SetInputFilename(target_filename);
212         reader->Update();
213         vvImage::Pointer vvimage = reader->GetOutput();
214         target = vvimage->GetFirstVTKImageData();
215     }
216
217     // allocate output
218     vtkImageData *output = vtkImageData::New();
219     output->SetExtent(target->GetWholeExtent());
220     output->SetOrigin(target->GetOrigin());
221     output->SetSpacing(target->GetSpacing());
222     output->SetScalarTypeToFloat();
223     output->AllocateScalars();
224
225     // fill output
226     unsigned long total = 0;
227     unsigned long over_one = 0;
228     for (int kk=0; kk<target->GetNumberOfPoints(); kk++) {
229         double *point = target->GetPoint(kk);
230         double value = target->GetPointData()->GetScalars()->GetTuple1(kk);
231         assert(value>=0);
232         assert(point[2]==0);
233         point[2] = value;
234
235         point[0] /= space_margin;
236         point[1] /= space_margin;
237         point[2] /= dose_margin;
238
239         double closest_point[3] = {0,0,0};
240         vtkIdType cell_id = 0;
241         int foo = 0;
242         double squared_distance = 0;
243
244         locator->FindClosestPoint(point,closest_point,cell_id,foo,squared_distance);
245
246         double distance = sqrt(squared_distance);
247         output->GetPointData()->GetScalars()->SetTuple1(kk,distance);
248
249         if (value>1) over_one++;
250         total++;
251
252     }
253
254     if (verbose) {
255         cout << "total=" << total << endl;
256         cout << "over_one=" << over_one << endl;
257         cout << "ratio=" << static_cast<float>(over_one)/total << endl;
258     }
259
260     locator->Delete();
261     target->Delete();
262
263     saveImage(output,gamma_filename);
264     output->Delete();
265
266     return 0;
267 }
268