Convert an itk::Image to vtkImageData#

Synopsis#

Convert an itk::Image to vtkImageData in a 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:

vtkImageData (0x51828b0)
  Debug: Off
  Modified Time: 240
  Reference Count: 2
  Registered Events: (none)
  Information: 0x35d7810
  Data Released: False
  Global Release Data: Off
  UpdateTime: 241
  Field Data:
    Debug: Off
    Modified Time: 186
    Reference Count: 1
    Registered Events: (none)
    Number Of Arrays: 0
    Number Of Components: 0
    Number Of Tuples: 0
  Number Of Points: 39277
  Number Of Cells: 38880
  Cell Data:
    Debug: Off
    Modified Time: 194
    Reference Count: 1
    Registered Events: (none)
    Number Of Arrays: 0
    Number Of Components: 0
    Number Of Tuples: 0
    Copy Tuple Flags: ( 1 1 1 1 1 0 1 1 )
    Interpolate Flags: ( 1 1 1 1 1 0 0 1 )
    Pass Through Flags: ( 1 1 1 1 1 1 1 1 )
    Scalars: (none)
    Vectors: (none)
    Normals: (none)
    TCoords: (none)
    Tensors: (none)
    GlobalIds: (none)
    PedigreeIds: (none)
    EdgeFlag: (none)
  Point Data:
    Debug: Off
    Modified Time: 240
    Reference Count: 1
    Registered Events: (none)
    Number Of Arrays: 1
    Array 0 name = scalars
    Number Of Components: 1
    Number Of Tuples: 39277
    Copy Tuple Flags: ( 1 1 1 1 1 0 1 1 )
    Interpolate Flags: ( 1 1 1 1 1 0 0 1 )
    Pass Through Flags: ( 1 1 1 1 1 1 1 1 )
    Scalars:
      Debug: Off
      Modified Time: 240
      Reference Count: 1
      Registered Events: (none)
      Name: scalars
      Data type: unsigned char
      Size: 39277
      MaxId: 39276
      NumberOfComponents: 1
      Information: 0
      Name: scalars
      Number Of Components: 1
      Number Of Tuples: 39277
      Size: 39277
      MaxId: 39276
      LookupTable: (none)
      Array: 0x51a7150
    Vectors: (none)
    Normals: (none)
    TCoords: (none)
    Tensors: (none)
    GlobalIds: (none)
    PedigreeIds: (none)
    EdgeFlag: (none)
  Bounds:
    Xmin,Xmax: (0, 180)
    Ymin,Ymax: (0, 216)
    Zmin,Zmax: (0, 0)
  Compute Time: 254
  Spacing: (1, 1, 1)
  Origin: (0, 0, 0)
  Dimensions: (181, 217, 1)
  Increments: (0, 0, 0)
  Extent: (0, 180, 0, 216, 0, 0)

Code#

Python#

#!/usr/bin/env python

import sys
import itk
import argparse

from distutils.version import StrictVersion

if StrictVersion(itk.Version.GetITKVersion()) < StrictVersion("5.2.0"):
    print("ITK 5.2.0 is required.")
    sys.exit(1)

parser = argparse.ArgumentParser(description="Convert An itk Image to vtk Image Data.")
parser.add_argument("input_image")
args = parser.parse_args()

inputImage = itk.imread(args.input_image)

vtkImage = itk.vtk_image_from_image(inputImage)

print(vtkImage)

C++#

#include "itkImageFileReader.h"
#include "itkImageToVTKImageFilter.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>;

  const auto input = itk::ReadImage<ImageType>(inputFileName);

  using FilterType = itk::ImageToVTKImageFilter<ImageType>;
  auto filter = FilterType::New();
  filter->SetInput(input);

  try
  {
    filter->Update();
  }
  catch (const itk::ExceptionObject & error)
  {
    std::cerr << "Error: " << error << std::endl;
    return EXIT_FAILURE;
  }

  vtkImageData * myvtkImageData = filter->GetOutput();
  myvtkImageData->Print(std::cout);

  return EXIT_SUCCESS;
}

Classes demonstrated#

template<typename TInputImage>
class ImageToVTKImageFilter : public itk::ProcessObject

Converts an ITK image into a VTK image and plugs a itk data pipeline to a VTK datapipeline.

This class puts together an itkVTKImageExporter and a vtkImageImporter. 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 an itk::Image can be plugged as input and a vtkImage is produced as output.

ITK Sphinx Examples:

See itk::ImageToVTKImageFilter for additional documentation.