ITK  6.0.0
Insight Toolkit
Examples/Filtering/GaussianBlurImageFunction.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"
int
main(int argc, char * argv[])
{
if (argc < 5)
{
std::cerr << "Usage: " << std::endl;
std::cerr << argv[0]
<< " inputImageFile outputImageFile sigma maxKernelWidth"
<< std::endl;
return EXIT_FAILURE;
}
using ImageType = itk::Image<float, 2>;
using ReaderType = itk::ImageFileReader<ImageType>;
using ConstIteratorType = itk::ImageRegionConstIterator<ImageType>;
auto reader = ReaderType::New();
reader->SetFileName(argv[1]);
reader->Update();
const ImageType * inputImage = reader->GetOutput();
ImageType::RegionType region = inputImage->GetBufferedRegion();
ConstIteratorType it(inputImage, region);
auto output = ImageType::New();
output->SetRegions(region);
output->SetOrigin(inputImage->GetOrigin());
output->SetSpacing(inputImage->GetSpacing());
output->Allocate();
IteratorType out(output, region);
auto gaussianFunction = GFunctionType::New();
gaussianFunction->SetInputImage(inputImage);
GFunctionType::ErrorArrayType setError;
setError.Fill(0.01);
gaussianFunction->SetMaximumError(setError);
gaussianFunction->SetSigma(std::stod(argv[3]));
gaussianFunction->SetMaximumKernelWidth(std::stoi(argv[4]));
it.GoToBegin();
out.GoToBegin();
while (!it.IsAtEnd())
{
out.Set(gaussianFunction->EvaluateAtIndex(it.GetIndex()));
++it;
++out;
}
using WriterType = itk::ImageFileWriter<ImageType>;
auto writer = WriterType::New();
writer->SetFileName(argv[2]);
writer->SetInput(output);
writer->Update();
return EXIT_SUCCESS;
}
itkImageFileReader.h
itkGaussianBlurImageFunction.h
itkImage.h
itkImageRegionIterator.h
itk::ImageFileReader
Data source that reads image data from a single file.
Definition: itkImageFileReader.h:75
itk::ImageRegionIterator
A multi-dimensional iterator templated over image type that walks a region of pixels.
Definition: itkImageRegionIterator.h:80
itk::GaussianBlurImageFunction
Compute the convolution of a neighborhood operator with the image at a specific location in space,...
Definition: itkGaussianBlurImageFunction.h:39
itk::ImageFileWriter
Writes image data to a single file.
Definition: itkImageFileWriter.h:90
itk::GTest::TypedefsAndConstructors::Dimension2::RegionType
ImageBaseType::RegionType RegionType
Definition: itkGTestTypedefsAndConstructors.h:54
itkImageFileWriter.h
itk::ImageRegionConstIterator< ImageType >
itk::Image
Templated n-dimensional image class.
Definition: itkImage.h:88
New
static Pointer New()