ITK  5.4.0
Insight Toolkit
SphinxExamples/src/Filtering/ImageFilterBase/ApplyKernelToEveryPixel/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 FloatImageType = itk::Image<float, 2>;
void
CreateImage(FloatImageType::Pointer image);
void
CastRescaleAndWrite(FloatImageType::Pointer image, const std::string & filename);
int
main()
{
auto image = FloatImageType::New();
CreateImage(image);
CastRescaleAndWrite(image, "input.png");
using SobelOperatorType = itk::SobelOperator<float, 2>;
SobelOperatorType sobelOperator;
itk::Size<2> radius;
radius.Fill(1); // a radius of 1x1 creates a 3x3 operator
sobelOperator.SetDirection(0); // Create the operator for the X axis derivative
sobelOperator.CreateToRadius(radius);
using NeighborhoodOperatorImageFilterType = itk::NeighborhoodOperatorImageFilter<FloatImageType, FloatImageType>;
filter->SetOperator(sobelOperator);
filter->SetInput(image);
filter->Update();
CastRescaleAndWrite(filter->GetOutput(), "output.png");
return EXIT_SUCCESS;
}
void
CreateImage(FloatImageType::Pointer image)
{
start.Fill(0);
size.Fill(100);
FloatImageType::RegionType region(start, size);
image->SetRegions(region);
image->Allocate();
image->FillBuffer(0);
// Make a square
for (unsigned int r = 20; r < 80; ++r)
{
for (unsigned int c = 20; c < 80; ++c)
{
pixelIndex[0] = r;
pixelIndex[1] = c;
image->SetPixel(pixelIndex, 15);
}
}
}
void
CastRescaleAndWrite(FloatImageType::Pointer image, const std::string & filename)
{
using UnsignedCharImageType = itk::Image<unsigned char, 2>;
auto rescaleFilter = RescaleFilterType::New();
rescaleFilter->SetInput(image);
rescaleFilter->SetOutputMinimum(0);
rescaleFilter->SetOutputMaximum(255);
rescaleFilter->Update();
itk::WriteImage(rescaleFilter->GetOutput(), filename);
}
Pointer
SmartPointer< Self > Pointer
Definition: itkAddImageFilter.h:93
itk::SobelOperator
A NeighborhoodOperator for performing a directional Sobel edge-detection operation at a pixel locatio...
Definition: itkSobelOperator.h:98
itk::Size< 2 >
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
itkSobelOperator.h
itkRescaleIntensityImageFilter.h
itkImageFileWriter.h
itkNeighborhoodOperatorImageFilter.h
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::NeighborhoodOperatorImageFilter
Applies a single NeighborhoodOperator to an image region.
Definition: itkNeighborhoodOperatorImageFilter.h:51
itk::WriteImage
ITK_TEMPLATE_EXPORT void WriteImage(TImagePointer &&image, const std::string &filename, bool compress=false)
Definition: itkImageFileWriter.h:254