]> Creatis software - clitk.git/blob - common/clitkDicomRTStruct2ImageFilter.cxx
GateAsciiImageIO is now cross-platform using itksys::RegularExpression
[clitk.git] / common / clitkDicomRTStruct2ImageFilter.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://www.centreleonberard.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 <iterator>
21 #include <algorithm>
22
23 // clitk
24 #include "clitkDicomRTStruct2ImageFilter.h"
25 #include "clitkImageCommon.h"
26
27 // vtk
28 #include <vtkPolyDataToImageStencil.h>
29 #include <vtkSmartPointer.h>
30 #include <vtkImageStencil.h>
31 #include <vtkLinearExtrusionFilter.h>
32 #include <vtkMetaImageWriter.h>
33
34
35 //--------------------------------------------------------------------
36 clitk::DicomRTStruct2ImageFilter::DicomRTStruct2ImageFilter()
37 {
38   mROI = NULL;
39   mWriteOutput = false;
40   mCropMask = true;
41 }
42 //--------------------------------------------------------------------
43
44
45 //--------------------------------------------------------------------
46 clitk::DicomRTStruct2ImageFilter::~DicomRTStruct2ImageFilter()
47 {
48
49 }
50 //--------------------------------------------------------------------
51
52 bool clitk::DicomRTStruct2ImageFilter::ImageInfoIsSet() const
53 {
54   return mSize.size() && mSpacing.size() && mOrigin.size();
55 }
56
57 //--------------------------------------------------------------------
58 void clitk::DicomRTStruct2ImageFilter::SetROI(clitk::DicomRT_ROI * roi)
59 {
60   mROI = roi;
61 }
62 //--------------------------------------------------------------------
63
64
65 //--------------------------------------------------------------------
66 void clitk::DicomRTStruct2ImageFilter::SetCropMaskEnabled(bool b)
67 {
68   mCropMask = b;
69 }
70 //--------------------------------------------------------------------
71
72
73 //--------------------------------------------------------------------
74 void clitk::DicomRTStruct2ImageFilter::SetOutputImageFilename(std::string s)
75 {
76   mOutputFilename = s;
77   mWriteOutput = true;
78 }
79 //--------------------------------------------------------------------
80
81
82 //--------------------------------------------------------------------
83 void clitk::DicomRTStruct2ImageFilter::SetImageFilename(std::string f)
84 {
85   itk::ImageIOBase::Pointer header = clitk::readImageHeader(f);
86   if (header->GetNumberOfDimensions() < 3) {
87     std::cerr << "Error. Please provide a 3D image instead of " << f << std::endl;
88     exit(0);
89   }
90   if (header->GetNumberOfDimensions() > 3) {
91     std::cerr << "Warning dimension > 3 are ignored" << std::endl;
92   }
93   mSpacing.resize(3);
94   mOrigin.resize(3);
95   mSize.resize(3);
96   for(unsigned int i=0; i<3; i++) {
97     mSpacing[i] = header->GetSpacing(i);
98     mOrigin[i] = header->GetOrigin(i);
99     mSize[i] = header->GetDimensions(i);
100   }
101 }
102 //--------------------------------------------------------------------
103
104 void clitk::DicomRTStruct2ImageFilter::SetOutputOrigin(const double* origin)
105 {
106   std::copy(origin,origin+3,std::back_inserter(mOrigin));
107 }
108 //--------------------------------------------------------------------
109 void clitk::DicomRTStruct2ImageFilter::SetOutputSpacing(const double* spacing)
110 {
111   std::copy(spacing,spacing+3,std::back_inserter(mSpacing));
112 }
113 //--------------------------------------------------------------------
114 void clitk::DicomRTStruct2ImageFilter::SetOutputSize(const unsigned long* size)
115 {
116   std::copy(size,size+3,std::back_inserter(mSize));
117 }
118
119 //--------------------------------------------------------------------
120 void clitk::DicomRTStruct2ImageFilter::Update()
121 {
122   if (!mROI) {
123     std::cerr << "Error. No ROI set, please use SetROI." << std::endl;
124     exit(0);
125   }
126   if (!ImageInfoIsSet()) {
127     std::cerr << "Error. Please provide image info (spacing/origin) with SetImageFilename" << std::endl;
128     exit(0);
129   }
130
131   // Get Mesh
132   vtkPolyData * mesh = mROI->GetMesh();
133
134   // Get bounds
135   double *bounds=mesh->GetBounds();
136   // for(int i=0; i<6; i++){
137 //     DD(bounds[i]);
138 //   }
139
140   // Compute origin
141   std::vector<double> origin;
142   origin.resize(3);
143   origin[0] = floor((bounds[0]-mOrigin[0])/mSpacing[0]-2)*mSpacing[0]+mOrigin[0];
144   origin[1] = floor((bounds[2]-mOrigin[1])/mSpacing[1]-2)*mSpacing[1]+mOrigin[1];
145   origin[2] = floor((bounds[4]-mOrigin[2])/mSpacing[2]-2)*mSpacing[2]+mOrigin[2];
146
147   // Compute extend
148   std::vector<double> extend;
149   extend.resize(3);
150   extend[0] = ceil((bounds[1]-origin[0])/mSpacing[0]+4);
151   extend[1] = ceil((bounds[3]-origin[1])/mSpacing[1]+4);
152   extend[2] = ceil((bounds[5]-origin[2])/mSpacing[2]+4);
153
154   // If no crop, set initial image size/origin
155   if (!mCropMask) {
156     for(int i=0; i<3; i++) {
157       origin[i] = mOrigin[i];
158       extend[i] = mSize[i]-1;
159     }
160   }
161
162   // Create new output image
163   mBinaryImage = vtkSmartPointer<vtkImageData>::New();
164   mBinaryImage->SetScalarTypeToUnsignedChar();
165   mBinaryImage->SetOrigin(&origin[0]);
166   mBinaryImage->SetSpacing(&mSpacing[0]);
167   mBinaryImage->SetExtent(0, extend[0],
168                           0, extend[1],
169                           0, extend[2]);
170   mBinaryImage->AllocateScalars();
171
172   // for(int i=0; i<3; i++){
173   //     DD(origin[i]);
174   //     DD(extend[i]);
175   //     DD(mBinaryImage->GetDimensions()[i]);
176   //   }
177   memset(mBinaryImage->GetScalarPointer(), 0,
178          mBinaryImage->GetDimensions()[0]*mBinaryImage->GetDimensions()[1]*mBinaryImage->GetDimensions()[2]*sizeof(unsigned char));
179
180   // Extrude
181   vtkSmartPointer<vtkLinearExtrusionFilter> extrude=vtkSmartPointer<vtkLinearExtrusionFilter>::New();
182   extrude->SetInput(mesh);
183   ///We extrude in the -slice_spacing direction to respect the FOCAL convention (NEEDED !)
184   extrude->SetVector(0, 0, -mSpacing[2]);
185
186   // Binarization
187   vtkSmartPointer<vtkPolyDataToImageStencil> sts=vtkSmartPointer<vtkPolyDataToImageStencil>::New();
188   //The following line is extremely important
189   //http://www.nabble.com/Bug-in-vtkPolyDataToImageStencil--td23368312.html#a23370933
190   sts->SetTolerance(0);
191   sts->SetInformationInput(mBinaryImage);
192   sts->SetInput(extrude->GetOutput());
193   //sts->SetInput(mesh);
194
195   vtkSmartPointer<vtkImageStencil> stencil=vtkSmartPointer<vtkImageStencil>::New();
196   stencil->SetStencil(sts->GetOutput());
197   stencil->SetInput(mBinaryImage);
198   stencil->ReverseStencilOn();
199   stencil->Update();
200
201   /*
202   vtkSmartPointer<vtkMetaImageWriter> w = vtkSmartPointer<vtkMetaImageWriter>::New();
203   w->SetInput(stencil->GetOutput());
204   w->SetFileName("binary2.mhd");
205   w->Write();
206   */
207
208   mBinaryImage->ShallowCopy(stencil->GetOutput());
209
210   if (mWriteOutput) {
211     typedef itk::Image<unsigned char, 3> ImageType;
212     typedef itk::VTKImageToImageFilter<ImageType> ConnectorType;
213     ConnectorType::Pointer connector = ConnectorType::New();
214     connector->SetInput(GetOutput());
215     connector->Update();
216     clitk::writeImage<ImageType>(connector->GetOutput(), mOutputFilename);
217   }
218 }
219 //--------------------------------------------------------------------
220
221
222
223 //--------------------------------------------------------------------
224 vtkImageData * clitk::DicomRTStruct2ImageFilter::GetOutput()
225 {
226   assert(mBinaryImage);
227   return mBinaryImage;
228 }
229 //--------------------------------------------------------------------
230