SimpleITK  
Segmentation/ConnectedThresholdImageFilter.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.
*
*=========================================================================*/
#if defined(_MSC_VER)
# pragma warning(disable : 4786)
#endif
#include "sitkImage.h"
#include <stdlib.h>
#include <iostream>
namespace sitk = itk::simple;
int
main(int argc, char * argv[])
{
//
// Check command line parameters
//
if (argc < 7)
{
std::cerr << "Missing Parameters " << std::endl;
std::cerr << "Usage: " << argv[0];
std::cerr << " inputImage outputImage lowerThreshold upperThreshold seedX seedY [seed2X seed2Y ... ]" << std::endl;
return 1;
}
//
// Read the image
//
reader.SetFileName(std::string(argv[1]));
sitk::Image image = reader.Execute();
//
// Blur using CurvatureFlowImageFilter
//
blurFilter.SetNumberOfIterations(5);
blurFilter.SetTimeStep(0.125);
image = blurFilter.Execute(image);
//
// Set up ConnectedThresholdImageFilter for segmentation
//
segmentationFilter.SetLower(atof(argv[3]));
segmentationFilter.SetUpper(atof(argv[4]));
segmentationFilter.SetReplaceValue(255);
for (int i = 5; i + 1 < argc; i += 2)
{
std::vector<unsigned int> seed = { (unsigned int)atoi(argv[i]), (unsigned int)atoi(argv[i + 1]) };
segmentationFilter.AddSeed(seed);
std::cout << "Adding a seed at: ";
for (unsigned int j = 0; j < seed.size(); ++j)
{
std::cout << seed[j] << " ";
}
std::cout << std::endl;
}
sitk::Image outImage = segmentationFilter.Execute(image);
//
// Write out the resulting file
//
writer.SetFileName(std::string(argv[2]));
writer.Execute(outImage);
return 0;
}
itk::simple::Image
The Image class for SimpleITK.
Definition: sitkImage.h:76
itk::simple::CurvatureFlowImageFilter::Execute
Image Execute(Image &&image1)
itk::simple::ConnectedThresholdImageFilter::SetLower
Self & SetLower(double Lower)
Definition: sitkConnectedThresholdImageFilter.h:74
itk::simple::ConnectedThresholdImageFilter::Execute
Image Execute(const Image &image1)
itk::simple::ConnectedThresholdImageFilter::SetReplaceValue
Self & SetReplaceValue(uint8_t ReplaceValue)
Definition: sitkConnectedThresholdImageFilter.h:94
sitkImageFileWriter.h
sitkImage.h
itk::simple::ImageFileWriter::Execute
Self & Execute(const Image &)
itk::simple::ImageFileReader
Read an image file and return a SimpleITK Image.
Definition: sitkImageFileReader.h:74
itk::simple::ConnectedThresholdImageFilter::AddSeed
Self & AddSeed(std::vector< unsigned int > point)
Add SeedList point.
Definition: sitkConnectedThresholdImageFilter.h:66
itk::simple::CurvatureFlowImageFilter
Denoise an image using curvature driven flow.
Definition: sitkCurvatureFlowImageFilter.h:71
itk::simple::CurvatureFlowImageFilter::SetNumberOfIterations
Self & SetNumberOfIterations(uint32_t NumberOfIterations)
Definition: sitkCurvatureFlowImageFilter.h:98
itk::simple::ImageFileReader::SetFileName
Self & SetFileName(const PathType &fn)
itk::simple::CurvatureFlowImageFilter::SetTimeStep
Self & SetTimeStep(double TimeStep)
Definition: sitkCurvatureFlowImageFilter.h:89
sitkImageFileReader.h
itk::simple::ConnectedThresholdImageFilter
Label pixels that are connected to a seed and lie within a range of values.
Definition: sitkConnectedThresholdImageFilter.h:41
sitkCurvatureFlowImageFilter.h
sitkConnectedThresholdImageFilter.h
itk::simple::ConnectedThresholdImageFilter::SetUpper
Self & SetUpper(double Upper)
Definition: sitkConnectedThresholdImageFilter.h:84
itk::simple::ImageFileReader::Execute
Image Execute() override
Set/Get The output PixelType of the image.
itk::simple::ImageFileWriter::SetFileName
Self & SetFileName(const PathType &fileName)
itk::simple
Definition: sitkAdditionalProcedures.h:28
itk::simple::ImageFileWriter
Write out a SimpleITK image to the specified file location.
Definition: sitkImageFileWriter.h:51