Convert vtkImageData to an itk::Image#
Synopsis#
Convert vtkImageData to an itk::Image in a processing pipeline. This transfers the image buffer data along with image size, pixel spacing, and origin. Since vtkImageData does not yet support an orientation matrix, the direction cosines are lost. This requires Module_ITKVtkGlue to be turned on in ITK’s CMake configuration.
Results#
Output:
Image (0x3b9d810)
RTTI typeinfo: itk::Image<unsigned char, 3u>
Reference Count: 2
Modified Time: 62
Debug: Off
Object Name:
Observers:
none
Source: (0x47e5460)
Source output name: Primary
Release Data: Off
Data Released: False
Global Release Data: Off
PipelineMTime: 57
UpdateMTime: 63
RealTimeStamp: 0 seconds
LargestPossibleRegion:
Dimension: 3
Index: [0, 0, 0]
Size: [181, 217, 1]
BufferedRegion:
Dimension: 3
Index: [0, 0, 0]
Size: [181, 217, 1]
RequestedRegion:
Dimension: 3
Index: [0, 0, 0]
Size: [181, 217, 1]
Spacing: [1, 1, 1]
Origin: [0, 0, 0]
Direction:
1 0 0
0 1 0
0 0 1
IndexToPointMatrix:
1 0 0
0 1 0
0 0 1
PointToIndexMatrix:
1 0 0
0 1 0
0 0 1
Inverse Direction:
1 0 0
0 1 0
0 0 1
PixelContainer:
ImportImageContainer (0x4875620)
RTTI typeinfo: itk::ImportImageContainer<unsigned long, unsigned char>
Reference Count: 1
Modified Time: 61
Debug: Off
Object Name:
Observers:
none
Pointer: 0x334e260
Container manages memory: false
Size: 39277
Capacity: 39277
Code#
Python#
#!/usr/bin/env python
import itk
import vtk
import argparse
parser = argparse.ArgumentParser(description="Convert vtk Image Data To An itk Image.")
parser.add_argument("input_image")
args = parser.parse_args()
Dimension = 2
PixelType = itk.UC
ImageType = itk.Image[PixelType, Dimension]
reader = vtk.vtkPNGReader()
reader.SetFileName(args.input_image)
reader.SetDataScalarTypeToUnsignedChar()
magnitude = vtk.vtkImageMagnitude()
magnitude.SetInputConnection(reader.GetOutputPort())
magnitude.Update()
vtkToItkFilter = itk.VTKImageToImageFilter[ImageType].New()
vtkToItkFilter.SetInput(magnitude.GetOutput())
vtkToItkFilter.Update()
myitkImage = vtkToItkFilter.GetOutput()
print(myitkImage)
C++#
#include "vtkSmartPointer.h"
#include "vtkPNGReader.h"
#include "vtkImageMagnitude.h"
#include "itkVTKImageToImageFilter.h"
int
main(int argc, char * argv[])
{
if (argc != 2)
{
std::cerr << "Usage: " << std::endl;
std::cerr << argv[0];
std::cerr << " <InputFileName>";
std::cerr << std::endl;
return EXIT_FAILURE;
}
const char * inputFileName = argv[1];
constexpr unsigned int Dimension = 2;
using PixelType = unsigned char;
using ImageType = itk::Image<PixelType, Dimension>;
vtkSmartPointer<vtkPNGReader> reader = vtkSmartPointer<vtkPNGReader>::New();
reader->SetFileName(inputFileName);
reader->SetDataScalarTypeToUnsignedChar();
vtkSmartPointer<vtkImageMagnitude> magnitude = vtkSmartPointer<vtkImageMagnitude>::New();
magnitude->SetInputConnection(reader->GetOutputPort());
magnitude->Update();
using FilterType = itk::VTKImageToImageFilter<ImageType>;
auto filter = FilterType::New();
filter->SetInput(magnitude->GetOutput());
try
{
filter->Update();
}
catch (const itk::ExceptionObject & error)
{
std::cerr << "Error: " << error << std::endl;
return EXIT_FAILURE;
}
ImageType::ConstPointer myitkImage = filter->GetOutput();
myitkImage->Print(std::cout);
return EXIT_SUCCESS;
}
Classes demonstrated#
-
template<typename TOutputImage>
class VTKImageToImageFilter : public itk::VTKImageImport<TOutputImage> Converts a VTK image into an ITK image and plugs a VTK data pipeline to an ITK datapipeline.
This class puts together an itk::VTKImageImport and a vtk::ImageExport. It takes care of the details related to the connection of ITK and VTK pipelines. The User will perceive this filter as an adaptor to which a vtkImageData can be plugged as input and an itk::Image is produced as output.
- ITK Sphinx Examples: