ITK  5.4.0
Insight Toolkit
SphinxExamples/src/Nonunit/Review/SegmentBloodVesselsWithMultiScaleHessianBasedMeasure/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.
*
*=========================================================================*/
int
main(int argc, char * argv[])
{
if (argc < 3)
{
std::cerr << "Usage: " << std::endl;
std::cerr << argv[0];
std::cerr << " <InputFileName> <OutputFileName>";
std::cerr << " [SigmaMinimum] [SigmaMaximum]";
std::cerr << " [NumberOfSigmaSteps]";
std::cerr << std::endl;
return EXIT_FAILURE;
}
const char * inputFileName = argv[1];
const char * outputFileName = argv[2];
double sigmaMinimum = 1.0;
if (argc > 3)
{
sigmaMinimum = std::stod(argv[3]);
}
double sigmaMaximum = 10.0;
if (argc > 4)
{
sigmaMaximum = std::stod(argv[4]);
}
unsigned int numberOfSigmaSteps = 10;
if (argc > 5)
{
numberOfSigmaSteps = std::stoi(argv[5]);
}
constexpr unsigned int Dimension = 2;
using PixelType = float;
const auto input = itk::ReadImage<ImageType>(inputFileName);
using HessianImageType = itk::Image<HessianPixelType, Dimension>;
auto objectnessFilter = ObjectnessFilterType::New();
objectnessFilter->SetBrightObject(false);
objectnessFilter->SetScaleObjectnessMeasure(false);
objectnessFilter->SetAlpha(0.5);
objectnessFilter->SetBeta(1.0);
objectnessFilter->SetGamma(5.0);
using MultiScaleEnhancementFilterType =
auto multiScaleEnhancementFilter = MultiScaleEnhancementFilterType::New();
multiScaleEnhancementFilter->SetInput(input);
multiScaleEnhancementFilter->SetHessianToMeasureFilter(objectnessFilter);
multiScaleEnhancementFilter->SetSigmaStepMethodToLogarithmic();
multiScaleEnhancementFilter->SetSigmaMinimum(sigmaMinimum);
multiScaleEnhancementFilter->SetSigmaMaximum(sigmaMaximum);
multiScaleEnhancementFilter->SetNumberOfSigmaSteps(numberOfSigmaSteps);
using OutputImageType = itk::Image<unsigned char, Dimension>;
auto rescaleFilter = RescaleFilterType::New();
rescaleFilter->SetInput(multiScaleEnhancementFilter->GetOutput());
try
{
itk::WriteImage(rescaleFilter->GetOutput(), outputFileName);
}
catch (const itk::ExceptionObject & error)
{
std::cerr << "Error: " << error << std::endl;
return EXIT_FAILURE;
}
return EXIT_SUCCESS;
}
itkImageFileReader.h
itkHessianToObjectnessMeasureImageFilter.h
itk::SymmetricSecondRankTensor
Represent a symmetric tensor of second rank.
Definition: itkSymmetricSecondRankTensor.h:75
itk::MultiScaleHessianBasedMeasureImageFilter
A filter to enhance structures using Hessian eigensystem-based measures in a multiscale framework.
Definition: itkMultiScaleHessianBasedMeasureImageFilter.h:87
itkRescaleIntensityImageFilter.h
itkImageFileWriter.h
itkMultiScaleHessianBasedMeasureImageFilter.h
itk::HessianToObjectnessMeasureImageFilter
A filter to enhance M-dimensional objects in N-dimensional images.
Definition: itkHessianToObjectnessMeasureImageFilter.h:61
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::GTest::TypedefsAndConstructors::Dimension2::Dimension
constexpr unsigned int Dimension
Definition: itkGTestTypedefsAndConstructors.h:44
itk::WriteImage
ITK_TEMPLATE_EXPORT void WriteImage(TImagePointer &&image, const std::string &filename, bool compress=false)
Definition: itkImageFileWriter.h:254