ITK  5.4.0
Insight Toolkit
Examples/Filtering/VesselnessMeasureImageFilter.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.
*
*=========================================================================*/
// The example takes an image (say MRA image), computes the vesselness measure
// of the image using the HessianRecursiveGaussianImageFilter and the
// Hessian3DToVesselnessMeasureImageFilter. The goal is to detect bright
// tubular structures in the image.
int
main(int argc, char * argv[])
{
if (argc < 3)
{
std::cerr << "Usage: inputImage outputImage [sigma] [alpha_1] [alpha_2]"
<< std::endl;
}
constexpr unsigned int Dimension = 3;
using InputPixelType = double;
using OutputPixelType = float;
using InputImageType = itk::Image<InputPixelType, Dimension>;
using OutputImageType = itk::Image<OutputPixelType, Dimension>;
using HessianFilterType =
using VesselnessMeasureFilterType =
auto hessianFilter = HessianFilterType::New();
auto vesselnessFilter = VesselnessMeasureFilterType::New();
auto reader = ReaderType::New();
auto writer = WriterType::New();
reader->SetFileName(argv[1]);
hessianFilter->SetInput(reader->GetOutput());
if (argc >= 4)
{
hessianFilter->SetSigma(static_cast<double>(std::stod(argv[3])));
}
vesselnessFilter->SetInput(hessianFilter->GetOutput());
writer->SetInput(vesselnessFilter->GetOutput());
writer->SetFileName(argv[2]);
if (argc >= 5)
{
vesselnessFilter->SetAlpha1(static_cast<double>(std::stod(argv[4])));
}
if (argc >= 6)
{
vesselnessFilter->SetAlpha2(static_cast<double>(std::stod(argv[5])));
}
writer->Update();
return EXIT_SUCCESS;
}
itkHessian3DToVesselnessMeasureImageFilter.h
itkImageFileReader.h
itkHessianRecursiveGaussianImageFilter.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
itkImageFileWriter.h
itk::HessianRecursiveGaussianImageFilter
Computes the Hessian matrix of an image by convolution with the Second and Cross derivatives of a Gau...
Definition: itkHessianRecursiveGaussianImageFilter.h:47
itk::Image
Templated n-dimensional image class.
Definition: itkImage.h:88
New
static Pointer New()
itk::Hessian3DToVesselnessMeasureImageFilter
Line filter to provide a vesselness measure for tubular objects from the hessian matrix.
Definition: itkHessian3DToVesselnessMeasureImageFilter.h:77
itk::GTest::TypedefsAndConstructors::Dimension2::Dimension
constexpr unsigned int Dimension
Definition: itkGTestTypedefsAndConstructors.h:44