]> Creatis software - clitk.git/blob - common/clitkDicomRT_ROI_ConvertToImageFilter.cxx
- correct crop 2D bug
[clitk.git] / common / clitkDicomRT_ROI_ConvertToImageFilter.cxx
1 /*=========================================================================
2   Program:         vv http://www.creatis.insa-lyon.fr/rio/vv
3   Main authors :   XX XX XX
4
5   Authors belongs to: 
6   - University of LYON           http://www.universite-lyon.fr/
7   - Léon Bérard cancer center    http://oncora1.lyon.fnclcc.fr
8   - CREATIS CNRS laboratory      http://www.creatis.insa-lyon.fr
9
10   This software is distributed WITHOUT ANY WARRANTY; without even
11   the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
12   PURPOSE.  See the copyright notices for more information.
13
14   It is distributed under dual licence
15   - BSD       http://www.opensource.org/licenses/bsd-license.php
16   - CeCILL-B  http://www.cecill.info/licences/Licence_CeCILL-B_V1-en.html
17
18   =========================================================================*/
19
20 #include "clitkDicomRT_ROI_ConvertToImageFilter.h" 
21 #include <vtkPolyDataToImageStencil.h>
22 #include <vtkSmartPointer.h>
23 #include <vtkImageStencil.h>
24 #include <vtkLinearExtrusionFilter.h>
25 #include <itkVTKImageToImageFilter.h>
26 #include "clitkImageCommon.h"
27
28 //--------------------------------------------------------------------
29 clitk::DicomRT_ROI_ConvertToImageFilter::DicomRT_ROI_ConvertToImageFilter() {
30   mROI = NULL;
31   mImageInfoIsSet = false;
32   mWriteOutput = false;
33   mCropMask = true;
34 }
35 //--------------------------------------------------------------------
36
37
38 //--------------------------------------------------------------------
39 clitk::DicomRT_ROI_ConvertToImageFilter::~DicomRT_ROI_ConvertToImageFilter() {
40   
41 }
42 //--------------------------------------------------------------------
43
44
45 //--------------------------------------------------------------------
46 void clitk::DicomRT_ROI_ConvertToImageFilter::SetROI(clitk::DicomRT_ROI * roi) {
47   mROI = roi;
48 }
49 //--------------------------------------------------------------------
50
51
52 //--------------------------------------------------------------------
53 void clitk::DicomRT_ROI_ConvertToImageFilter::SetCropMaskEnabled(bool b) {
54   mCropMask = b;
55 }
56 //--------------------------------------------------------------------
57
58
59 //--------------------------------------------------------------------
60 void clitk::DicomRT_ROI_ConvertToImageFilter::SetOutputImageFilename(std::string s) {
61   mOutputFilename = s;
62   mWriteOutput = true;
63 }
64 //--------------------------------------------------------------------
65
66
67 //--------------------------------------------------------------------
68 void clitk::DicomRT_ROI_ConvertToImageFilter::SetImageFilename(std::string f) {
69   itk::ImageIOBase::Pointer header = clitk::readImageHeader(f);
70   if (header->GetNumberOfDimensions() < 3) {
71     std::cerr << "Error. Please provide a 3D image instead of " << f << std::endl;
72     exit(0);
73   }
74   if (header->GetNumberOfDimensions() > 3) {
75     std::cerr << "Warning dimension > 3 are ignored" << std::endl;
76   }
77   mSpacing.resize(3);
78   mOrigin.resize(3);
79   mSize.resize(3);
80   for(unsigned int i=0; i<3; i++) {
81     mSpacing[i] = header->GetSpacing(i);
82     mOrigin[i] = header->GetOrigin(i);
83     mSize[i] = header->GetDimensions(i);
84   }
85   mImageInfoIsSet = true;
86 }
87 //--------------------------------------------------------------------
88
89
90 //--------------------------------------------------------------------
91 void clitk::DicomRT_ROI_ConvertToImageFilter::Update() {
92   if (!mROI) {
93     std::cerr << "Error. No ROI set, please use SetROI." << std::endl;
94     exit(0);
95   }
96   if (!mImageInfoIsSet) {
97     std::cerr << "Error. Please provide image info (spacing/origin) with SetImageFilename" << std::endl;
98     exit(0);
99   }
100   // DD("Update");
101   
102   // Get Mesh
103   vtkPolyData * mesh = mROI->GetMesh();
104   DD(mesh->GetNumberOfCells());
105   
106   // Get bounds
107   double *bounds=mesh->GetBounds(); 
108   // for(int i=0; i<6; i++){
109 //     DD(bounds[i]);
110 //   }
111
112   // Compute origin
113   std::vector<double> origin; 
114   origin.resize(3);
115   origin[0] = floor((bounds[0]-mOrigin[0])/mSpacing[0]-2)*mSpacing[0]+mOrigin[0];
116   origin[1] = floor((bounds[2]-mOrigin[1])/mSpacing[1]-2)*mSpacing[1]+mOrigin[1];
117   origin[2] = floor((bounds[4]-mOrigin[2])/mSpacing[2]-2)*mSpacing[2]+mOrigin[2];
118   
119   // Compute extend
120   std::vector<double> extend; 
121   extend.resize(3);
122   extend[0] = ceil((bounds[1]-origin[0])/mSpacing[0]+4);
123   extend[1] = ceil((bounds[3]-origin[1])/mSpacing[1]+4);
124   extend[2] = ceil((bounds[5]-origin[2])/mSpacing[2]+4);
125   
126   // If no crop, set initial image size/origin
127   if (!mCropMask) {
128     for(int i=0; i<3; i++) {
129       origin[i] = mOrigin[i];
130       extend[i] = mSize[i]-1;
131     }
132   }
133   
134   // Create new output image
135   mBinaryImage = vtkImageData::New();
136   mBinaryImage->SetScalarTypeToUnsignedChar();
137   mBinaryImage->SetOrigin(&origin[0]);
138   mBinaryImage->SetSpacing(&mSpacing[0]);
139   mBinaryImage->SetExtent(0, extend[0], 
140                           0, extend[1], 
141                           0, extend[2]);
142   mBinaryImage->AllocateScalars();
143
144   // for(int i=0; i<3; i++){
145   //     DD(origin[i]);
146   //     DD(extend[i]);
147   //     DD(mBinaryImage->GetDimensions()[i]);
148   //   }
149   memset(mBinaryImage->GetScalarPointer(), 0,
150          mBinaryImage->GetDimensions()[0]*mBinaryImage->GetDimensions()[1]*mBinaryImage->GetDimensions()[2]*sizeof(unsigned char));
151   
152   // Extrude
153   vtkSmartPointer<vtkLinearExtrusionFilter> extrude=vtkSmartPointer<vtkLinearExtrusionFilter>::New();
154   extrude->SetInput(mesh);
155   ///We extrude in the -slice_spacing direction to respect the FOCAL convention // ?????????????
156   extrude->SetVector(0, 0, -mSpacing[2]);
157
158   // Binarization  
159   vtkSmartPointer<vtkPolyDataToImageStencil> sts=vtkSmartPointer<vtkPolyDataToImageStencil>::New();
160   //The following line is extremely important
161   //http://www.nabble.com/Bug-in-vtkPolyDataToImageStencil--td23368312.html#a23370933
162   sts->SetTolerance(0);
163   sts->SetInformationInput(mBinaryImage);
164   sts->SetInput(extrude->GetOutput());
165   //sts->SetInput(mesh);
166
167   vtkSmartPointer<vtkImageStencil> stencil=vtkSmartPointer<vtkImageStencil>::New();
168   stencil->SetStencil(sts->GetOutput());  
169   stencil->SetInput(mBinaryImage);
170   stencil->ReverseStencilOn();
171   stencil->Update();
172   mBinaryImage->ShallowCopy(stencil->GetOutput());
173   
174   if (mWriteOutput) {
175      typedef itk::Image<unsigned char, 3> ImageType;
176      typedef itk::VTKImageToImageFilter<ImageType> ConnectorType;
177      ConnectorType::Pointer connector = ConnectorType::New();
178      connector->SetInput(GetOutput());
179      connector->Update();
180      clitk::writeImage<ImageType>(connector->GetOutput(), mOutputFilename);  
181   }
182 }
183 //--------------------------------------------------------------------
184
185
186     
187 //--------------------------------------------------------------------
188 vtkImageData * clitk::DicomRT_ROI_ConvertToImageFilter::GetOutput() {
189   return mBinaryImage;
190 }
191 //--------------------------------------------------------------------