ITK  5.4.0
Insight Toolkit
SphinxExamples/src/Filtering/BinaryMathematicalMorphology/ThinImage/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>;
static void
CreateImage(ImageType::Pointer image);
int
main(int argc, char * argv[])
{
auto image = ImageType::New();
if (argc == 2)
{
image = itk::ReadImage<ImageType>(argv[1]);
}
else
{
CreateImage(image);
itk::WriteImage(image, "Input.png");
}
using BinaryThinningImageFilterType = itk::BinaryThinningImageFilter<ImageType, ImageType>;
auto binaryThinningImageFilter = BinaryThinningImageFilterType::New();
binaryThinningImageFilter->SetInput(image);
binaryThinningImageFilter->Update();
// Rescale the image so that it can be seen (the output is 0 and 1, we want 0 and 255)
auto rescaler = RescaleType::New();
rescaler->SetInput(binaryThinningImageFilter->GetOutput());
rescaler->SetOutputMinimum(0);
rescaler->SetOutputMaximum(255);
rescaler->Update();
itk::WriteImage(rescaler->GetOutput(), "Output.png");
return EXIT_SUCCESS;
}
void
CreateImage(ImageType::Pointer image)
{
// Create an image
start.Fill(0);
size.Fill(100);
ImageType::RegionType region(start, size);
image->SetRegions(region);
image->Allocate();
image->FillBuffer(0);
// Draw a 5 pixel wide line
for (unsigned int i = 20; i < 80; ++i)
{
for (unsigned int j = 50; j < 55; ++j)
{
index[0] = i;
index[1] = j;
image->SetPixel(index, 255);
}
}
}
itkBinaryThinningImageFilter.h
Pointer
SmartPointer< Self > Pointer
Definition: itkAddImageFilter.h:93
itk::Index
Represent a n-dimensional index in a n-dimensional image.
Definition: itkIndex.h:70
itkImageFileReader.h
itk::GTest::TypedefsAndConstructors::Dimension2::SizeType
ImageBaseType::SizeType SizeType
Definition: itkGTestTypedefsAndConstructors.h:49
itkImage.h
itk::Index::Fill
void Fill(IndexValueType value)
Definition: itkIndex.h:274
itk::Size::Fill
void Fill(SizeValueType value)
Definition: itkSize.h:213
itk::GTest::TypedefsAndConstructors::Dimension2::IndexType
ImageBaseType::IndexType IndexType
Definition: itkGTestTypedefsAndConstructors.h:50
itk::GTest::TypedefsAndConstructors::Dimension2::RegionType
ImageBaseType::RegionType RegionType
Definition: itkGTestTypedefsAndConstructors.h:54
itkRescaleIntensityImageFilter.h
itkImageFileWriter.h
itk::BinaryThinningImageFilter
This filter computes one-pixel-wide edges of the input image.
Definition: itkBinaryThinningImageFilter.h:62
itk::RescaleIntensityImageFilter
Applies a linear transformation to the intensity levels of the input Image.
Definition: itkRescaleIntensityImageFilter.h:133
itk::Image
Templated n-dimensional image class.
Definition: itkImage.h:88
New
static Pointer New()
itk::WriteImage
ITK_TEMPLATE_EXPORT void WriteImage(TImagePointer &&image, const std::string &filename, bool compress=false)
Definition: itkImageFileWriter.h:254