FAST  3.2.0
Framework for Heterogeneous Medical Image Computing and Visualization
ITKImageExporter.hpp
Go to the documentation of this file.
1 #ifndef ITKIMAGEEXPORTER_HPP_
2 #define ITKIMAGEEXPORTER_HPP_
3 
4 #include "FAST/ProcessObject.hpp"
5 #include "FAST/Data/Image.hpp"
6 
7 #include <itkImageSource.h>
8 #include <itkImageRegionIterator.h>
9 
10 namespace fast {
11 
12 template<class TImage>
13 class FAST_EXPORT ITKImageExporter: public itk::ImageSource<TImage>, public ProcessObject {
14  public:
17  typedef itk::ImageSource<TImage> Superclass;
18  typedef itk::SmartPointer<Self> Pointer;
19 
21  itkNewMacro(Self);
22 
24  itkTypeMacro(MyImageSource, ImageSource);
25 
26  std::string getNameOfClass() const { return "ITKImageExporter"; };
27  private:
29  void execute() {};
30 
31  // Is called by ITK
32  void GenerateData();
33 
34  template <class T>
35  void transferDataToITKImage(Image::pointer input);
36 
37 };
38 
39 } // end namespace fast
40 
41 template<class TImage>
43  createInputPort<Image>(0);
44 }
45 
46 template <class TImage>
47 template <class T>
49 
50  typename TImage::Pointer output = this->GetOutput();
51  typename TImage::RegionType region;
52  typename TImage::IndexType start;
53  start[0] = 0;
54  start[1] = 0;
55  if(input->getDimensions() == 3)
56  start[2] = 0;
57 
58  typename TImage::SizeType size;
59  size[0] = input->getWidth();
60  size[1] = input->getHeight();
61  if(input->getDimensions() == 3)
62  size[2] = input->getDepth();
63 
64  region.SetSize(size);
65  region.SetIndex(start);
66 
67  output->SetRegions(region);
68  output->Allocate();
69 
70  ImageAccess::pointer access = input->getImageAccess(ACCESS_READ);
71  T * fastPixelData = (T*)access->get();
72 
73  itk::ImageRegionIterator<TImage> imageIterator(output,
74  output->GetLargestPossibleRegion());
75  unsigned int width = input->getWidth();
76  if(input->getDimensions() == 2) {
77  while (!imageIterator.IsAtEnd()) {
78  unsigned int x = imageIterator.GetIndex()[0];
79  unsigned int y = imageIterator.GetIndex()[1];
80  imageIterator.Set(fastPixelData[x + y*width]);
81 
82  ++imageIterator;
83  }
84  } else {
85  unsigned int height = input->getHeight();
86 
87  while (!imageIterator.IsAtEnd()) {
88  unsigned int x = imageIterator.GetIndex()[0];
89  unsigned int y = imageIterator.GetIndex()[1];
90  unsigned int z = imageIterator.GetIndex()[2];
91  imageIterator.Set(fastPixelData[x + y*width + z*width*height]);
92 
93  ++imageIterator;
94  }
95  }
96 }
97 
98 template<class TImage>
100 
101  update(); // Update FAST pipeline
102 
103  Image::pointer input = getInputData<Image>();
104  if(input->getDimensions() != TImage::ImageDimension)
105  throw Exception("The dimension of the input and output images of the ITKImageExporter are unequal.");
106 
107  if(input->getNrOfChannels() != 1)
108  throw Exception("The ITKImageExporter currently doesn't support images with multiple components");
109 
110 
111  // Transfer data from mInput to vtk image
112  switch(input->getDataType()) {
113  fastSwitchTypeMacro(transferDataToITKImage<FAST_TYPE>(input))
114  }
115 
116 
117 }
118 
119 #endif // ITKIMAGEEXPORTER_HPP_
fast::ITKImageExporter::Self
ITKImageExporter Self
Definition: ITKImageExporter.hpp:16
Image.hpp
fast
Definition: AffineTransformation.hpp:7
fast::SpatialDataObject::pointer
std::shared_ptr< SpatialDataObject > pointer
Definition: SpatialDataObject.hpp:12
fastSwitchTypeMacro
#define fastSwitchTypeMacro(call)
Definition: DataTypes.hpp:55
fast::ITKImageExporter::Superclass
itk::ImageSource< TImage > Superclass
Definition: ITKImageExporter.hpp:17
ACCESS_READ
@ ACCESS_READ
Definition: Access.hpp:5
ProcessObject.hpp
fast::ITKImageExporter::Pointer
itk::SmartPointer< Self > Pointer
Definition: ITKImageExporter.hpp:18
fast::ITKImageExporter
Definition: ITKImageExporter.hpp:13
fast::ITKImageExporter::getNameOfClass
std::string getNameOfClass() const
Definition: ITKImageExporter.hpp:26
fast::ProcessObject
Definition: ProcessObject.hpp:22