View Component Vector Image as Scalar Image#

Synopsis#

View a component of a vector image as if it were a scalar image

Results#

Output:

1

Code#

C++#

#include "itkVectorImage.h"
#include "itkImage.h"
#include "itkVectorImageToImageAdaptor.h"
#include "itkImageRegionIterator.h"

using ScalarImageType = itk::Image<float, 2>;
using VectorImageType = itk::VectorImage<float, 2>;

void
CreateImage(VectorImageType::Pointer image);

int
main()
{
  auto image = VectorImageType::New();
  CreateImage(image);

  using ImageAdaptorType = itk::VectorImageToImageAdaptor<float, 2>;
  auto adaptor = ImageAdaptorType::New();
  adaptor->SetExtractComponentIndex(0);
  adaptor->SetImage(image);

  itk::Index<2> index;
  index[0] = 0;
  index[1] = 0;

  std::cout << adaptor->GetPixel(index) << std::endl;

  return EXIT_SUCCESS;
}

void
CreateImage(VectorImageType::Pointer image)
{
  VectorImageType::IndexType start;
  start.Fill(0);

  VectorImageType::SizeType size;
  size.Fill(2);

  VectorImageType::RegionType region;
  region.SetSize(size);
  region.SetIndex(start);

  image->SetRegions(region);
  image->SetNumberOfComponentsPerPixel(2);
  image->Allocate();

  using VariableVectorType = itk::VariableLengthVector<double>;
  VariableVectorType variableLengthVector;
  variableLengthVector.SetSize(2);

  itk::ImageRegionIterator<VectorImageType> imageIterator(image, image->GetLargestPossibleRegion());
  variableLengthVector[0] = 1;
  variableLengthVector[1] = 2;

  while (!imageIterator.IsAtEnd())
  {
    imageIterator.Set(variableLengthVector);

    ++imageIterator;
  }
}

Classes demonstrated#

template<typename TPixelType, unsigned int Dimension>
class VectorImageToImageAdaptor : public itk::ImageAdaptor<VectorImage<TPixelType, Dimension>, Accessor::VectorImageToImagePixelAccessor<TPixelType>>

Presents a VectorImage and extracts a component from it into an image.

The class is expected to be templated over a pixel type and dimension. These are the pixel types and dimension of the VectorImage.

The component to extract is set with SetExtractComponentIdx() method.

Note

This work is part of the National Alliance for Medical Image Computing (NAMIC), funded by the National Institutes of Health through the NIH Roadmap for Medical Research, Grant U54 EB005149.

ITK Sphinx Examples:

See itk::VectorImageToImageAdaptor for additional documentation.