ITK  5.2.0
Insight Toolkit
Examples/Filtering/LaplacianImageFilter.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
*
* http://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.
*
*=========================================================================*/
int
main(int argc, char * argv[])
{
if (argc < 3)
{
std::cerr << "Usage: " << std::endl;
std::cerr << argv[0] << " inputImage outputImage " << std::endl;
return EXIT_FAILURE;
}
const char * inputFilename = argv[1];
const char * outputFilename = argv[2];
using CharPixelType = unsigned char; // IO
using RealPixelType = double; // Operations
constexpr unsigned int Dimension = 2;
using CharImageType = itk::Image<CharPixelType, Dimension>;
using RealImageType = itk::Image<RealPixelType, Dimension>;
using CastToRealFilterType =
using CastToCharFilterType =
using RescaleFilter =
using LaplacianFilter =
using ZeroCrossingFilter =
// Setting the IO
ReaderType::Pointer reader = ReaderType::New();
WriterType::Pointer writer = WriterType::New();
CastToRealFilterType::Pointer toReal = CastToRealFilterType::New();
CastToCharFilterType::Pointer toChar = CastToCharFilterType::New();
RescaleFilter::Pointer rescale = RescaleFilter::New();
// Setting the ITK pipeline filter
LaplacianFilter::Pointer lapFilter = LaplacianFilter::New();
ZeroCrossingFilter::Pointer zeroFilter = ZeroCrossingFilter::New();
reader->SetFileName(inputFilename);
writer->SetFileName(outputFilename);
// The output of an edge filter is 0 or 1
rescale->SetOutputMinimum(0);
rescale->SetOutputMaximum(255);
toReal->SetInput(reader->GetOutput());
toChar->SetInput(rescale->GetOutput());
writer->SetInput(toChar->GetOutput());
// Edge Detection by Laplacian Image Filter:
lapFilter->SetInput(toReal->GetOutput());
zeroFilter->SetInput(lapFilter->GetOutput());
rescale->SetInput(zeroFilter->GetOutput());
try
{
writer->Update();
}
catch (const itk::ExceptionObject & err)
{
std::cout << "ExceptionObject caught !" << std::endl;
std::cout << err << std::endl;
return EXIT_FAILURE;
}
return EXIT_SUCCESS;
}
itk::CastImageFilter
Casts input pixels to output pixel type.
Definition: itkCastImageFilter.h:104
itk::ZeroCrossingImageFilter
This filter finds the closest pixel to the zero-crossings (sign changes) in a signed itk::Image.
Definition: itkZeroCrossingImageFilter.h:64
itk::LaplacianImageFilter
This filter computes the Laplacian of a scalar-valued image.
Definition: itkLaplacianImageFilter.h:63
itkImageFileReader.h
itkCastImageFilter.h
itkLaplacianImageFilter.h
itk::ImageFileReader
Data source that reads image data from a single file.
Definition: itkImageFileReader.h:75
itk::ImageFileWriter
Writes image data to a single file.
Definition: itkImageFileWriter.h:88
itkRescaleIntensityImageFilter.h
itkImageFileWriter.h
itkZeroCrossingBasedEdgeDetectionImageFilter.h
itk::RescaleIntensityImageFilter
Applies a linear transformation to the intensity levels of the input Image.
Definition: itkRescaleIntensityImageFilter.h:154
itk::Image
Templated n-dimensional image class.
Definition: itkImage.h:86
itk::GTest::TypedefsAndConstructors::Dimension2::Dimension
constexpr unsigned int Dimension
Definition: itkGTestTypedefsAndConstructors.h:44