Regional Minimal#
Synopsis#
Regional minimal image filter.
Results#
Code#
C++#
#include "itkImage.h"
#include "itkImageFileWriter.h"
#include "itkRegionalMinimaImageFilter.h"
using ImageType = itk::Image<unsigned char, 2>;
static void
CreateImage(ImageType::Pointer image);
int
main()
{
auto image = ImageType::New();
CreateImage(image);
using RegionalMinimaImageFilterType = itk::RegionalMinimaImageFilter<ImageType, ImageType>;
auto regionalMinimaImageFilter = RegionalMinimaImageFilterType::New();
regionalMinimaImageFilter->SetInput(image);
regionalMinimaImageFilter->Update();
itk::WriteImage(image, "input.png");
itk::WriteImage(regionalMinimaImageFilter->GetOutput(), "output.png");
return EXIT_SUCCESS;
}
void
CreateImage(ImageType::Pointer image)
{
ImageType::IndexType start;
start.Fill(0);
ImageType::SizeType size;
size[0] = 200;
size[1] = 300;
ImageType::RegionType region(start, size);
image->SetRegions(region);
image->Allocate();
// Make two intensity blobs
for (unsigned int r = 0; r < size[1]; ++r)
{
for (unsigned int c = 0; c < size[0]; ++c)
{
ImageType::IndexType pixelIndex;
pixelIndex[0] = c;
pixelIndex[1] = r;
double c1 = c - 100.0;
double c2 = c - 200.0;
double rr = r - 100.0;
// purposely use 270,257 since it is > 255
double v1 = 270.0 - std::sqrt(rr * rr + c1 * c1);
double v2 = 257.0 - std::sqrt(rr * rr + c2 * c2);
double maxv = v1;
if (maxv < v2)
maxv = v2;
double val = maxv;
// Clip
if (val < 0.0)
val = 0.0;
if (val > 255.0)
val = 255.0;
image->SetPixel(pixelIndex, static_cast<unsigned char>(val));
}
}
}
Classes demonstrated#
-
template<typename TInputImage, typename TOutputImage>
class RegionalMinimaImageFilter : public itk::ImageToImageFilter<TInputImage, TOutputImage> Produce a binary image where foreground is the regional minima of the input image.
Regional minima are flat zones surrounded by pixels of greater value.
If the input image is constant, the entire image can be considered as a minima or not. The SetFlatIsMinima() method let the user choose which behavior to use.
This class was contributed to the Insight Journal by
- Author
Gaetan Lehmann. Biologie du Developpement et de la Reproduction, INRA de Jouy-en-Josas, France. https://www.insight-journal.org/browse/publication/65
- See
RegionalMaximaImageFilter
- See
ValuedRegionalMinimaImageFilter
- See
HConcaveImageFilter
- ITK Sphinx Examples: