VMTK
vtkvmtkITKFilterUtilities.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3 Program: VMTK
4 Module: vtkvmtkITKFilterUtilities.h
5 Language: C++
6 
7  Copyright (c) Luca Antiga, David Steinman. All rights reserved.
8  See LICENSE file for details.
9 
10  Portions of this code are covered under the VTK copyright.
11  See VTKCopyright.txt or http://www.kitware.com/VTKCopyright.htm
12  for details.
13 
14  Portions of this code are covered under the ITK copyright.
15  See ITKCopyright.txt or http://www.itk.org/HTML/Copyright.htm
16  for details.
17 
18  This software is distributed WITHOUT ANY WARRANTY; without even
19  the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
20  PURPOSE. See the above copyright notices for more information.
21 
22 =========================================================================*/
23 
33 #ifndef __vtkvmtkITKFilterUtilities_h
34 #define __vtkvmtkITKFilterUtilities_h
35 
37 #include "vtkvmtkWin32Header.h"
38 
39 #include "vtkImageData.h"
40 #include "itkImage.h"
41 #include "itkCommand.h"
42 
43 class VTK_VMTK_SEGMENTATION_EXPORT vtkvmtkITKFilterUtilities
44 {
45 public:
46 
47  template<typename TImage>
48  static void
49  VTKToITKImage(vtkImageData* input, typename TImage::Pointer output) {
50 
51  typedef TImage ImageType;
52  typedef typename ImageType::Pointer ImagePointer;
53  typedef typename ImageType::PixelType PixelType;
54 
55  int dims[3];
56  input->GetDimensions(dims);
57  double spacing[3];
58  input->GetSpacing(spacing);
59  double origin[3];
60  input->GetOrigin(origin);
61  int extent[6];
62  input->GetExtent(extent);
63 
64  output->GetPixelContainer()->SetImportPointer(static_cast<PixelType*>(input->GetScalarPointer()),dims[0]*dims[1]*dims[2],false);
65  typename ImageType::RegionType region;
66  typename ImageType::IndexType index;
67  typename ImageType::SizeType size;
68  //index[0] = index[1] = index[2] = 0;
69  //size[0] = dims[0];
70  //size[1] = dims[1];
71  //size[2] = dims[2];
72  index[0] = extent[0];
73  index[1] = extent[2];
74  index[2] = extent[4];
75  size[0] = extent[1] - extent[0] + 1;
76  size[1] = extent[3] - extent[2] + 1;
77  size[2] = extent[5] - extent[4] + 1;
78  region.SetIndex(index);
79  region.SetSize(size);
80  output->SetLargestPossibleRegion(region);
81  output->SetBufferedRegion(region);
82  output->SetSpacing(spacing);
83  output->SetOrigin(origin);
84  }
85 
86  template<typename TImage>
87  static void
88  VTKToITKVectorImage(vtkImageData* input, typename TImage::Pointer output) {
89 
90  typedef TImage ImageType;
91  typedef typename ImageType::Pointer ImagePointer;
92  typedef typename ImageType::PixelType PixelType;
93  typedef typename ImageType::InternalPixelType InternalPixelType;
94 
95  int dims[3];
96  input->GetDimensions(dims);
97  double spacing[3];
98  input->GetSpacing(spacing);
99  int components = input->GetNumberOfScalarComponents();
100  double origin[3];
101  input->GetOrigin(origin);
102  int extent[6];
103  input->GetExtent(extent);
104 
105  output->GetPixelContainer()->SetImportPointer(static_cast<InternalPixelType*>(input->GetScalarPointer()),dims[0]*dims[1]*dims[2]*components,false);
106  typename ImageType::RegionType region;
107  typename ImageType::IndexType index;
108  typename ImageType::SizeType size;
109  //index[0] = index[1] = index[2] = 0;
110  //size[0] = dims[0];
111  //size[1] = dims[1];
112  //size[2] = dims[2];
113  index[0] = extent[0];
114  index[1] = extent[2];
115  index[2] = extent[4];
116  size[0] = extent[1] - extent[0] + 1;
117  size[1] = extent[3] - extent[2] + 1;
118  size[2] = extent[5] - extent[4] + 1;
119  region.SetIndex(index);
120  region.SetSize(size);
121  output->SetLargestPossibleRegion(region);
122  output->SetBufferedRegion(region);
123  output->SetSpacing(spacing);
124  output->SetOrigin(origin);
125  output->SetVectorLength(components);
126  }
127 
128  template<typename TImage>
129  static void
130  ITKToVTKImage(typename TImage::Pointer input, vtkImageData* output) {
131 
132  typedef TImage ImageType;
133  typedef typename ImageType::Pointer ImagePointer;
134  typedef typename ImageType::PixelType PixelType;
135  typedef typename ImageType::PointType PointType;
136  typedef typename ImageType::SpacingType SpacingType;
137  typedef typename ImageType::RegionType RegionType;
138  typedef typename ImageType::IndexType IndexType;
139  typedef typename ImageType::SizeType SizeType;
140 
141  PointType origin = input->GetOrigin();
142  SpacingType spacing = input->GetSpacing();
143 
144  double outputOrigin[3];
145  double outputSpacing[3];
146 
147  outputOrigin[0] = origin[0];
148  outputOrigin[1] = origin[1];
149  outputOrigin[2] = origin[2];
150 
151  outputSpacing[0] = spacing[0];
152  outputSpacing[1] = spacing[1];
153  outputSpacing[2] = spacing[2];
154 
155  output->SetOrigin(outputOrigin);
156  output->SetSpacing(outputSpacing);
157 
158  RegionType region = input->GetBufferedRegion();
159  IndexType index = region.GetIndex();
160  SizeType size = region.GetSize();
161 
162  //int dimensions[3];
163  //dimensions[0] = size[0];
164  //dimensions[1] = size[1];
165  //dimensions[2] = size[2];
166  int extent[6];
167  extent[0] = index[0];
168  extent[1] = index[0] + size[0] - 1;
169  extent[2] = index[1];
170  extent[3] = index[1] + size[1] - 1;
171  extent[4] = index[2];
172  extent[5] = index[2] + size[2] - 1;
173 
174  int components = input->GetNumberOfComponentsPerPixel();
175  int dataType = output->GetScalarType(); // WARNING: we delegate setting type to caller
176 
177  //output->SetDimensions(dimensions);
178  output->SetExtent(extent);
179  output->AllocateScalars(dataType,components);
180 
181  memcpy(static_cast<PixelType*>(output->GetScalarPointer()),input->GetBufferPointer(),input->GetBufferedRegion().GetNumberOfPixels()*sizeof(PixelType));
182  }
183 
184  static void
185  ProgressCallback(itk::Object *o, const itk::EventObject &, void *data)
186  {
187  ((vtkAlgorithm*)data)->UpdateProgress(dynamic_cast<const itk::ProcessObject*>(o)->GetProgress());
188  }
189 
190  static void
191  ConnectProgress(itk::Object* obj, vtkAlgorithm* alg)
192  {
193  itk::CStyleCommand::Pointer progressCommand = itk::CStyleCommand::New();
194  progressCommand->SetCallback(vtkvmtkITKFilterUtilities::ProgressCallback);
195  progressCommand->SetClientData(alg);
196  obj->AddObserver(itk::ProgressEvent(),progressCommand);
197  }
198 
199 protected:
202 
203 private:
204  vtkvmtkITKFilterUtilities(const vtkvmtkITKFilterUtilities&); // Not implemented.
205  void operator=(const vtkvmtkITKFilterUtilities&); // Not implemented.
206 };
207 
208 #if 0
209 template< class TInputPixel, class TOutputPixel>
210 void
211 vtkvmtkSimpleImageToImageITKFilter<TInputPixel,TOutputPixel>::
212 SimpleExecute(vtkImageData *input, vtkImageData *output)
213 {
214  int inputDims[3];
215  input->GetDimensions(inputDims);
216  double inputSpacing[3];
217  input->GetSpacing(inputSpacing);
218 
219  InputImagePointer inImage = InputImageType::New();
220  inImage->GetPixelContainer()->SetImportPointer(static_cast<InputPixelType*>(input->GetScalarPointer()),inputDims[0]*inputDims[1]*inputDims[2],false);
221  typename InputImageType::RegionType inputRegion;
222  typename InputImageType::IndexType inputIndex;
223  typename InputImageType::SizeType inputSize;
224  inputIndex[0] = inputIndex[1] = inputIndex[2] = 0;
225  inputSize[0] = inputDims[0];
226  inputSize[1] = inputDims[1];
227  inputSize[2] = inputDims[2];
228  inputRegion.SetIndex(inputIndex);
229  inputRegion.SetSize(inputSize);
230  inImage->SetLargestPossibleRegion(inputRegion);
231  inImage->SetBufferedRegion(inputRegion);
232  inImage->SetSpacing(inputSpacing);
233 
234  int outputDims[3];
235  output->GetDimensions(outputDims);
236  double outputSpacing[3];
237  output->GetSpacing(outputSpacing);
238 
239  OutputImagePointer outImage = OutputImageType::New();
240  outImage->GetPixelContainer()->SetImportPointer(static_cast<OutputPixelType*>(output->GetScalarPointer()),outputDims[0]*outputDims[1]*outputDims[2],false);
241  typename OutputImageType::RegionType outputRegion;
242  typename OutputImageType::IndexType outputIndex;
243  typename OutputImageType::SizeType outputSize;
244  outputIndex[0] = outputIndex[1] = outputIndex[2] = 0;
245  outputSize[0] = outputDims[0];
246  outputSize[1] = outputDims[1];
247  outputSize[2] = outputDims[2];
248  outputRegion.SetIndex(outputIndex);
249  outputRegion.SetSize(outputSize);
250  outImage->SetLargestPossibleRegion(outputRegion);
251  outImage->SetBufferedRegion(outputRegion);
252  outImage->SetSpacing(outputSpacing);
253 
254  this->SimpleExecuteITK(inImage,outImage);
255 
256  memcpy(static_cast<OutputPixelType*>(output->GetScalarPointer()),outImage->GetBufferPointer(),outImage->GetBufferedRegion().GetNumberOfPixels()*sizeof(OutputPixelType));
257 }
258 
259 typedef vtkvmtkSimpleImageToImageITKFilter<float,float> vtkvmtkSimpleImageToImageITKFilterFF;
260 #endif
261 
262 #endif
263 
static void ProgressCallback(itk::Object *o, const itk::EventObject &, void *data)
static void ITKToVTKImage(typename TImage::Pointer input, vtkImageData *output)
static void VTKToITKVectorImage(vtkImageData *input, typename TImage::Pointer output)
static void VTKToITKImage(vtkImageData *input, typename TImage::Pointer output)
Abstract class for wrapping ITK filters.
static void ConnectProgress(itk::Object *obj, vtkAlgorithm *alg)