ITK  5.4.0
Insight Toolkit
SphinxExamples/src/Filtering/LabelMap/KeepRegionsThatMeetSpecific/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"
using ImageType = itk::Image<unsigned char, 2>;
using LabelImageType = itk::Image<unsigned char, 2>;
static void
CreateImage(ImageType::Pointer image);
int
main()
{
auto image = ImageType::New();
CreateImage(image);
// Create a ShapeLabelMap from the image
using BinaryImageToShapeLabelMapFilterType = itk::BinaryImageToShapeLabelMapFilter<ImageType>;
BinaryImageToShapeLabelMapFilterType::Pointer binaryImageToShapeLabelMapFilter =
binaryImageToShapeLabelMapFilter->SetInput(image);
binaryImageToShapeLabelMapFilter->Update();
// Remove label objects that have PERIMETER greater than 50
using ShapeOpeningLabelMapFilterType =
auto shapeOpeningLabelMapFilter = ShapeOpeningLabelMapFilterType::New();
shapeOpeningLabelMapFilter->SetInput(binaryImageToShapeLabelMapFilter->GetOutput());
shapeOpeningLabelMapFilter->SetLambda(50);
shapeOpeningLabelMapFilter->ReverseOrderingOn();
shapeOpeningLabelMapFilter->SetAttribute(ShapeOpeningLabelMapFilterType::LabelObjectType::PERIMETER);
shapeOpeningLabelMapFilter->Update();
// Create a label image
using LabelMapToLabelImageFilterType =
auto labelMapToLabelImageFilter = LabelMapToLabelImageFilterType::New();
labelMapToLabelImageFilter->SetInput(shapeOpeningLabelMapFilter->GetOutput());
labelMapToLabelImageFilter->Update();
using RGBPixelType = itk::RGBPixel<unsigned char>;
using RGBImageType = itk::Image<RGBPixelType, 2>;
// Color each label/object a different color
auto colormapImageFilter = RGBFilterType::New();
colormapImageFilter->SetInput(labelMapToLabelImageFilter->GetOutput());
colormapImageFilter->Update();
itk::WriteImage(colormapImageFilter->GetOutput(), "output.png");
return EXIT_SUCCESS;
}
void
CreateImage(ImageType::Pointer image)
{
// Create a black image with three white squares
start.Fill(0);
size.Fill(200);
region.SetSize(size);
region.SetIndex(start);
image->SetRegions(region);
image->Allocate();
itk::ImageRegionIterator<ImageType> imageIterator(image, image->GetLargestPossibleRegion());
// Make a square
while (!imageIterator.IsAtEnd())
{
if (((imageIterator.GetIndex()[0] > 5 && imageIterator.GetIndex()[0] < 10) &&
(imageIterator.GetIndex()[1] > 5 && imageIterator.GetIndex()[1] < 10)) ||
((imageIterator.GetIndex()[0] > 50 && imageIterator.GetIndex()[0] < 60) &&
(imageIterator.GetIndex()[1] > 50 && imageIterator.GetIndex()[1] < 60)) ||
((imageIterator.GetIndex()[0] > 100 && imageIterator.GetIndex()[0] < 130) &&
(imageIterator.GetIndex()[1] > 100 && imageIterator.GetIndex()[1] < 130)))
{
imageIterator.Set(255);
}
else
{
imageIterator.Set(0);
}
++imageIterator;
}
using WriterType = itk::ImageFileWriter<ImageType>;
auto writer = WriterType::New();
writer->SetFileName("input.png");
writer->SetInput(image);
writer->Update();
}
Pointer
SmartPointer< Self > Pointer
Definition: itkAddImageFilter.h:93
itk::RGBPixel
Represent Red, Green and Blue components for color images.
Definition: itkRGBPixel.h:58
itk::ScalarToRGBColormapImageFilter
Implements pixel-wise intensity->rgb mapping operation on one image.
Definition: itkScalarToRGBColormapImageFilter.h:130
itk::GTest::TypedefsAndConstructors::Dimension2::SizeType
ImageBaseType::SizeType SizeType
Definition: itkGTestTypedefsAndConstructors.h:49
itkImage.h
itkImageRegionIterator.h
itk::ScalarToRGBColormapImageFilterEnums::RGBColormapFilter::Jet
itkShapeOpeningLabelMapFilter.h
itk::Index::Fill
void Fill(IndexValueType value)
Definition: itkIndex.h:274
itk::Size::Fill
void Fill(SizeValueType value)
Definition: itkSize.h:213
itk::ImageRegionIterator
A multi-dimensional iterator templated over image type that walks a region of pixels.
Definition: itkImageRegionIterator.h:80
itk::GTest::TypedefsAndConstructors::Dimension2::IndexType
ImageBaseType::IndexType IndexType
Definition: itkGTestTypedefsAndConstructors.h:50
itk::BinaryImageToShapeLabelMapFilter
Converts a binary image to a label map and valuate the shape attributes.
Definition: itkBinaryImageToShapeLabelMapFilter.h:62
itk::ImageFileWriter
Writes image data to a single file.
Definition: itkImageFileWriter.h:88
itk::GTest::TypedefsAndConstructors::Dimension2::RegionType
ImageBaseType::RegionType RegionType
Definition: itkGTestTypedefsAndConstructors.h:54
itk::ShapeOpeningLabelMapFilter
Remove objects according to the value of their shape attribute.
Definition: itkShapeOpeningLabelMapFilter.h:49
itk::ImageRegionIterator::Set
void Set(const PixelType &value) const
Definition: itkImageRegionIterator.h:116
itkBinaryImageToShapeLabelMapFilter.h
itkImageFileWriter.h
itkScalarToRGBColormapImageFilter.h
itk::LabelMapToLabelImageFilter
Converts a LabelMap to a labeled image.
Definition: itkLabelMapToLabelImageFilter.h:46
itk::Image
Templated n-dimensional image class.
Definition: itkImage.h:88
New
static Pointer New()
itkLabelMapToLabelImageFilter.h
itk::ImageRegion::SetSize
void SetSize(const SizeType &size)
Definition: itkImageRegion.h:202
itk::WriteImage
ITK_TEMPLATE_EXPORT void WriteImage(TImagePointer &&image, const std::string &filename, bool compress=false)
Definition: itkImageFileWriter.h:254