Compute Histogram of Masked Region in Image#
Note
Wish List Still needs additional work to finish proper creation of example.
Synopsis#
Compute the histogram of a masked region of an image.
Results#
Note
Help Wanted Implementation of Results for sphinx examples containing this message. Reconfiguration of CMakeList.txt may be necessary. Write An Example <https://itk.org/ITKExamples/Documentation/Contribute/WriteANewExample.html>
Code#
C++#
#include "itkMaskedImageToHistogramFilter.h"
#include "itkImage.h"
#include "itkImageFileWriter.h"
#include "itkRGBPixel.h"
#include "itkImageRegionIteratorWithIndex.h"
#include "itkImageRegionIterator.h"
#include "itkRescaleIntensityImageFilter.h"
using RGBPixelType = itk::RGBPixel<unsigned char>;
using RGBImageType = itk::Image<RGBPixelType, 2>;
using UnsignedCharImageType = itk::Image<unsigned char, 2>;
static void
CreateImage(RGBImageType::Pointer image);
static void CreateHalfMask(itk::ImageRegion<2>, UnsignedCharImageType::Pointer mask);
int
main()
{
constexpr unsigned int MeasurementVectorSize = 3; // RGB
auto image = RGBImageType::New();
CreateImage(image);
auto mask = UnsignedCharImageType::New();
CreateHalfMask(image->GetLargestPossibleRegion(), mask);
using HistogramFilterType = itk::Statistics::MaskedImageToHistogramFilter<RGBImageType, UnsignedCharImageType>;
using HistogramMeasurementVectorType = HistogramFilterType::HistogramMeasurementVectorType;
using HistogramSizeType = HistogramFilterType::HistogramSizeType;
using HistogramType = HistogramFilterType::HistogramType;
auto histogramFilter = HistogramFilterType::New();
histogramFilter->SetInput(image);
histogramFilter->SetMaskImage(mask);
histogramFilter->SetAutoMinimumMaximum(true);
HistogramSizeType histogramSize(MeasurementVectorSize);
histogramSize[0] = 4; // number of bins for the Red channel
histogramSize[1] = 4; // number of bins for the Green channel
histogramSize[2] = 4; // number of bins for the Blue channel
histogramFilter->SetHistogramSize(histogramSize);
histogramFilter->SetMarginalScale(10); // Required (could this be set in the filter?)
histogramFilter->Update();
const HistogramType * histogram = histogramFilter->GetOutput();
HistogramType::ConstIterator histogramIterator = histogram->Begin();
// std::string filename = "/home/doriad/histogram.txt";
// std::ofstream fout(filename.c_str());
while (histogramIterator != histogram->End())
{
// std::cout << "Index = " << histogramIterator.GetMeasurementVector() << "Frequency = " <<
// histogramIterator.GetFrequency() << std::endl; std::cout << "Index = " << histogramIterator.GetIndex() <<
// "Frequency = " << histogramIterator.GetFrequency() << std::endl; fout << "Index = " <<
// histogram->GetIndex(histogramItr.GetMeasurementVector()) << "Frequency = " << histogramItr.GetFrequency() <<
// std::endl;
++histogramIterator;
}
// fout.close();
HistogramType::MeasurementVectorType mv(3);
mv[0] = 255;
mv[1] = 0;
mv[2] = 0;
std::cout << "Frequency = " << histogram->GetFrequency(histogram->GetIndex(mv)) << std::endl;
return EXIT_SUCCESS;
}
void
CreateImage(RGBImageType::Pointer image)
{
// Create a black image with a red square and a green square.
// This should produce a histogram with very strong spikes.
RGBImageType::SizeType size;
size[0] = 3;
size[1] = 3;
RGBImageType::IndexType start;
start[0] = 0;
start[1] = 0;
RGBImageType::RegionType region;
region.SetIndex(start);
region.SetSize(size);
image->SetRegions(region);
image->Allocate();
itk::ImageRegionIteratorWithIndex<RGBImageType> iterator(image, image->GetLargestPossibleRegion());
iterator.GoToBegin();
RGBPixelType redPixel;
redPixel.SetRed(255);
redPixel.SetGreen(0);
redPixel.SetBlue(0);
RGBPixelType blackPixel;
blackPixel.SetRed(0);
blackPixel.SetGreen(0);
blackPixel.SetBlue(0);
itk::ImageRegionIterator<RGBImageType> imageIterator(image, region);
while (!imageIterator.IsAtEnd())
{
imageIterator.Set(blackPixel);
++imageIterator;
}
RGBImageType::IndexType index;
index[0] = 0;
index[1] = 0;
image->SetPixel(index, redPixel);
index[0] = 1;
index[1] = 0;
image->SetPixel(index, redPixel);
}
void CreateHalfMask(itk::ImageRegion<2> region, UnsignedCharImageType::Pointer mask)
{
mask->SetRegions(region);
mask->Allocate();
mask->FillBuffer(0);
itk::Size<2> regionSize = region.GetSize();
itk::ImageRegionIterator<UnsignedCharImageType> imageIterator(mask, region);
// Make the left half of the mask white and the right half black
while (!imageIterator.IsAtEnd())
{
if (imageIterator.GetIndex()[0] > regionSize[0] / 2)
{
imageIterator.Set(0);
}
else
{
imageIterator.Set(1);
}
++imageIterator;
}
using RescaleFilterType = itk::RescaleIntensityImageFilter<UnsignedCharImageType, UnsignedCharImageType>;
auto rescaleFilter = RescaleFilterType::New();
rescaleFilter->SetInput(mask);
rescaleFilter->Update();
itk::WriteImage(rescaleFilter->GetOutput(), "mask.png");
}
Classes demonstrated#
-
template<typename TImage, typename TMaskImage>
class MaskedImageToHistogramFilter : public itk::Statistics::ImageToHistogramFilter<TImage> Generate a histogram from the masked pixels of an image.
This class expands the features of the ImageToHistogramFilter by adding a required MaskImage input image. Only the pixel in the input image where the MaskImage’s value is the MaskValue will be added to the computed histogram.
- ITK Sphinx Examples: