ITK  5.4.0
Insight Toolkit
SphinxExamples/src/Segmentation/ConnectedComponents/LabelConnectComponentsInGrayscaleImage/Code.cxx
/*=========================================================================
*
* Copyright NumFOCUS
*
* Licensed under the Apache License, Version 2.0 (the "License");
* you may not use this file except in compliance with the License.
* You may obtain a copy of the License at
*
* https://www.apache.org/licenses/LICENSE-2.0.txt
*
* Unless required by applicable law or agreed to in writing, software
* distributed under the License is distributed on an "AS IS" BASIS,
* WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
* See the License for the specific language governing permissions and
* limitations under the License.
*
*=========================================================================*/
#include "itkImage.h"
#include "itksys/SystemTools.hxx"
#include <sstream>
#ifdef ENABLE_QUICKVIEW
# include "QuickView.h"
#endif
template <typename TImage>
static void
CreateImage(TImage * const image);
template <typename TImage, typename TLabelImage>
static void
SummarizeLabelStatistics(TImage * image, TLabelImage * labelImage);
int
main(int argc, char * argv[])
{
constexpr unsigned int Dimension = 2;
using PixelType = short;
using RGBPixelType = itk::RGBPixel<unsigned char>;
using LabelPixelType = unsigned int;
using LabelImageType = itk::Image<LabelPixelType, Dimension>;
PixelType distanceThreshold = 4;
if (argc < 2)
{
image = ImageType::New();
CreateImage(image.GetPointer());
}
else
{
if (argc > 2)
{
distanceThreshold = static_cast<PixelType>(atoi(argv[2]));
}
image = itk::ReadImage<ImageType>(argv[1]);
}
connected->SetInput(image);
connected->SetDistanceThreshold(distanceThreshold);
auto relabel = RelabelFilterType::New();
RelabelFilterType::ObjectSizeType minSize = 20;
if (argc > 3)
{
minSize = std::stoi(argv[3]);
}
relabel->SetInput(connected->GetOutput());
relabel->SetMinimumObjectSize(minSize);
relabel->Update();
SummarizeLabelStatistics(image.GetPointer(), relabel->GetOutput());
auto rgbFilter = RGBFilterType::New();
rgbFilter->SetInput(relabel->GetOutput());
#ifdef ENABLE_QUICKVIEW
QuickView viewer;
viewer.AddImage(
image.GetPointer(), true, argc > 1 ? itksys::SystemTools::GetFilenameName(argv[1]) : "Generated image");
std::stringstream desc;
desc << "Scalar Connected Components:\n# of Objects: " << relabel->GetNumberOfObjects() << " Threshold: "
connected->GetDistanceThreshold())
<< " Min Size: " << relabel->GetMinimumObjectSize();
viewer.AddRGBImage(rgbFilter->GetOutput(), true, desc.str());
viewer.Visualize();
#endif
return EXIT_SUCCESS;
}
template <typename TImage>
void
CreateImage(TImage * const image)
{
// Create an image with 2 connected components
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);
}
}
}
template <typename TImage, typename TLabelImage>
void
SummarizeLabelStatistics(TImage * image, TLabelImage * labelImage)
{
using LabelStatisticsImageFilterType = itk::LabelStatisticsImageFilter<TImage, TLabelImage>;
auto labelStatisticsImageFilter = LabelStatisticsImageFilterType::New();
labelStatisticsImageFilter->SetLabelInput(labelImage);
labelStatisticsImageFilter->SetInput(image);
labelStatisticsImageFilter->UseHistogramsOn(); // needed to compute median
labelStatisticsImageFilter->Update();
std::cout << "Number of labels: " << labelStatisticsImageFilter->GetNumberOfLabels() << std::endl;
std::cout << std::endl;
using LabelPixelType = typename LabelStatisticsImageFilterType::LabelPixelType;
for (auto vIt = labelStatisticsImageFilter->GetValidLabelValues().begin();
vIt != labelStatisticsImageFilter->GetValidLabelValues().end();
++vIt)
{
if (labelStatisticsImageFilter->HasLabel(*vIt))
{
LabelPixelType labelValue = *vIt;
std::cout << "Label: " << *vIt << std::endl;
std::cout << "\tmin: " << labelStatisticsImageFilter->GetMinimum(labelValue) << std::endl;
std::cout << "\tmax: " << labelStatisticsImageFilter->GetMaximum(labelValue) << std::endl;
std::cout << "\tmedian: " << labelStatisticsImageFilter->GetMedian(labelValue) << std::endl;
std::cout << "\tmean: " << labelStatisticsImageFilter->GetMean(labelValue) << std::endl;
std::cout << "\tsigma: " << labelStatisticsImageFilter->GetSigma(labelValue) << std::endl;
std::cout << "\tvariance: " << labelStatisticsImageFilter->GetVariance(labelValue) << std::endl;
std::cout << "\tsum: " << labelStatisticsImageFilter->GetSum(labelValue) << std::endl;
std::cout << "\tcount: " << labelStatisticsImageFilter->GetCount(labelValue) << std::endl;
std::cout << "\tregion: " << labelStatisticsImageFilter->GetRegion(labelValue) << std::endl;
}
}
}
Pointer
SmartPointer< Self > Pointer
Definition: itkAddImageFilter.h:93
itkScalarConnectedComponentImageFilter.h
itk::ScalarConnectedComponentImageFilter
A connected components filter that labels the objects in an arbitrary image. Two pixels are similar i...
Definition: itkScalarConnectedComponentImageFilter.h:103
itk::RGBPixel
Represent Red, Green and Blue components for color images.
Definition: itkRGBPixel.h:58
itkImageFileReader.h
itk::GTest::TypedefsAndConstructors::Dimension2::SizeType
ImageBaseType::SizeType SizeType
Definition: itkGTestTypedefsAndConstructors.h:49
itkImage.h
itkLabelStatisticsImageFilter.h
itkRelabelComponentImageFilter.h
itk::IndexValueType
long IndexValueType
Definition: itkIntTypes.h:90
itkLabelToRGBImageFilter.h
itk::GTest::TypedefsAndConstructors::Dimension2::IndexType
ImageBaseType::IndexType IndexType
Definition: itkGTestTypedefsAndConstructors.h:50
QuickView::AddRGBImage
void AddRGBImage(TImage *, bool FlipVertical=true, std::string Description="")
QuickView.h
QuickView::AddImage
void AddImage(TImage *, bool FlipVertical=true, std::string Description="")
itk::GTest::TypedefsAndConstructors::Dimension2::RegionType
ImageBaseType::RegionType RegionType
Definition: itkGTestTypedefsAndConstructors.h:54
itk::NumericTraits::PrintType
T PrintType
Definition: itkNumericTraits.h:69
itk::LabelStatisticsImageFilter
Given an intensity image and a label map, compute min, max, variance and mean of the pixels associate...
Definition: itkLabelStatisticsImageFilter.h:63
itk::RelabelComponentImageFilter
Relabel the components in an image such that consecutive labels are used.
Definition: itkRelabelComponentImageFilter.h:82
itk::LabelToRGBImageFilter
Apply a colormap to a label image.
Definition: itkLabelToRGBImageFilter.h:50
QuickView
A convenient class to render itk images with vtk.
Definition: QuickView.h:111
itk::Image
Templated n-dimensional image class.
Definition: itkImage.h:88
New
static Pointer New()
QuickView::Visualize
void Visualize(bool interact=true)
itk::GTest::TypedefsAndConstructors::Dimension2::Dimension
constexpr unsigned int Dimension
Definition: itkGTestTypedefsAndConstructors.h:44