Shape Attributes for Binary Image#

Synopsis#

Compute shape attributes for a binary image.

Results#

Output:

File "Generated image" has 2 labels.
Label: 1
BoundingBox: ImageRegion (0x7fab2f0003b8)
Dimension: 2
Index: [20, 30]
Size: [60, 70]

  NumberOfPixels: 4200
  PhysicalSize: 4200
  Centroid: [49.5, 64.5]
  NumberOfPixelsOnBorder: 0
  PerimeterOnBorder: 0
  FeretDiameter: 0
  PrincipalMoments: [299.917, 408.25]
  PrincipalAxes: 1 0
  0 1

  Elongation: 1.16671
  Perimeter: 245.385
  Roundness: 0.936229
  EquivalentSphericalRadius: 36.5637
  EquivalentSphericalPerimeter: 229.736
  EquivalentEllipsoidDiameter: [67.7015, 78.988]
  Flatness: 1.16671
  PerimeterOnBorderRatio: 0
Label: 2
BoundingBox: ImageRegion (0x7fab2ef00308)
Dimension: 2
Index: [100, 115]
Size: [30, 45]

  NumberOfPixels: 1350
  PhysicalSize: 1350
  Centroid: [114.5, 137]
  NumberOfPixelsOnBorder: 0
  PerimeterOnBorder: 0
  FeretDiameter: 0
  PrincipalMoments: [74.9167, 168.667]
  PrincipalAxes: 1 0
  0 1

  Elongation: 1.50046
  Perimeter: 141.098
  Roundness: 0.923103
  EquivalentSphericalRadius: 20.7296
  EquivalentSphericalPerimeter: 130.248
  EquivalentEllipsoidDiameter: [33.8461, 50.7849]
  Flatness: 1.50046
  PerimeterOnBorderRatio: 0

Code#

C++#

#include "itkImage.h"
#include "itkImageFileReader.h"
#include "itkConnectedComponentImageFilter.h"
#include "itkLabelImageToShapeLabelMapFilter.h"

template <typename TImage>
static void
CreateImage(TImage * image);

int
main(int argc, char * argv[])
{
  constexpr unsigned int Dimension = 2;
  using PixelType = unsigned char;
  using LabelType = unsigned short;
  using InputImageType = itk::Image<PixelType, Dimension>;
  using OutputImageType = itk::Image<LabelType, Dimension>;
  using ShapeLabelObjectType = itk::ShapeLabelObject<LabelType, Dimension>;
  using LabelMapType = itk::LabelMap<ShapeLabelObjectType>;

  std::string             fileName;
  InputImageType::Pointer image;
  if (argc < 2)
  {
    image = InputImageType::New();
    CreateImage(image.GetPointer());
    fileName = "Generated image";
  }
  else
  {
    fileName = argv[1];
    image = itk::ReadImage<InputImageType>(fileName);
  }


  using ConnectedComponentImageFilterType = itk::ConnectedComponentImageFilter<InputImageType, OutputImageType>;
  using I2LType = itk::LabelImageToShapeLabelMapFilter<OutputImageType, LabelMapType>;

  auto connected = ConnectedComponentImageFilterType::New();
  connected->SetInput(image);
  connected->Update();

  using I2LType = itk::LabelImageToShapeLabelMapFilter<OutputImageType, LabelMapType>;
  auto i2l = I2LType::New();
  i2l->SetInput(connected->GetOutput());
  i2l->SetComputePerimeter(true);
  i2l->Update();

  LabelMapType * labelMap = i2l->GetOutput();
  std::cout << "File "
            << "\"" << fileName << "\""
            << " has " << labelMap->GetNumberOfLabelObjects() << " labels." << std::endl;

  // Retrieve all attributes
  for (unsigned int n = 0; n < labelMap->GetNumberOfLabelObjects(); ++n)
  {
    ShapeLabelObjectType * labelObject = labelMap->GetNthLabelObject(n);
    std::cout << "Label: " << itk::NumericTraits<LabelMapType::LabelType>::PrintType(labelObject->GetLabel())
              << std::endl;
    std::cout << "    BoundingBox: " << labelObject->GetBoundingBox() << std::endl;
    std::cout << "    NumberOfPixels: " << labelObject->GetNumberOfPixels() << std::endl;
    std::cout << "    PhysicalSize: " << labelObject->GetPhysicalSize() << std::endl;
    std::cout << "    Centroid: " << labelObject->GetCentroid() << std::endl;
    std::cout << "    NumberOfPixelsOnBorder: " << labelObject->GetNumberOfPixelsOnBorder() << std::endl;
    std::cout << "    PerimeterOnBorder: " << labelObject->GetPerimeterOnBorder() << std::endl;
    std::cout << "    FeretDiameter: " << labelObject->GetFeretDiameter() << std::endl;
    std::cout << "    PrincipalMoments: " << labelObject->GetPrincipalMoments() << std::endl;
    std::cout << "    PrincipalAxes: " << labelObject->GetPrincipalAxes() << std::endl;
    std::cout << "    Elongation: " << labelObject->GetElongation() << std::endl;
    std::cout << "    Perimeter: " << labelObject->GetPerimeter() << std::endl;
    std::cout << "    Roundness: " << labelObject->GetRoundness() << std::endl;
    std::cout << "    EquivalentSphericalRadius: " << labelObject->GetEquivalentSphericalRadius() << std::endl;
    std::cout << "    EquivalentSphericalPerimeter: " << labelObject->GetEquivalentSphericalPerimeter() << std::endl;
    std::cout << "    EquivalentEllipsoidDiameter: " << labelObject->GetEquivalentEllipsoidDiameter() << std::endl;
    std::cout << "    Flatness: " << labelObject->GetFlatness() << std::endl;
    std::cout << "    PerimeterOnBorderRatio: " << labelObject->GetPerimeterOnBorderRatio() << std::endl;
  }

  return EXIT_SUCCESS;
}

template <typename TImage>
void
CreateImage(TImage * const image)
{
  // Create an image with 2 objects
  typename TImage::IndexType start = { { 0, 0 } };
  start[0] = 0;
  start[1] = 0;

  typename TImage::SizeType size;
  unsigned int              NumRows = 200;
  unsigned int              NumCols = 300;
  size[0] = NumRows;
  size[1] = NumCols;

  typename TImage::RegionType region(start, size);

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

  // Make a square
  for (typename TImage::IndexValueType r = 20; r < 80; ++r)
  {
    for (typename TImage::IndexValueType c = 30; c < 100; ++c)
    {
      typename TImage::IndexType pixelIndex = { { r, c } };

      image->SetPixel(pixelIndex, 255);
    }
  }

  // Make another square
  for (typename TImage::IndexValueType r = 100; r < 130; ++r)
  {
    for (typename TImage::IndexValueType c = 115; c < 160; ++c)
    {
      typename TImage::IndexType pixelIndex = { { r, c } };

      image->SetPixel(pixelIndex, 255);
    }
  }
}

Classes demonstrated#

template<typename TInputImage, typename TOutputImage = LabelMap<ShapeLabelObject<typename TInputImage::PixelType, TInputImage::ImageDimension>>>
class LabelImageToShapeLabelMapFilter : public itk::ImageToImageFilter<TInputImage, TOutputImage>

Converts a label image to a label map and valuates the shape attributes.

A convenient class that converts a label image to a label map and valuates the shape attribute at once.

This implementation was taken from the Insight Journal paper: https://www.insight-journal.org/browse/publication/176

Author

Gaetan Lehmann. Biologie du Developpement et de la Reproduction, INRA de Jouy-en-Josas, France.

See

ShapeLabelObject, LabelShapeOpeningImageFilter, LabelStatisticsOpeningImageFilter

ITK Sphinx Examples:

See itk::LabelImageToShapeLabelMapFilter for additional documentation.