SimpleITK  2.0.0
ITKIntegration/ITKIntegration.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
// SimpleITK includes
#include "SimpleITK.h"
// ITK includes
#include "itkImage.h"
// create convenient namespace alias
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();
//
// Set up writer
//
writer.SetFileName( std::string( argv[2] ) );
// Blur using CurvatureFlowImageFilter
//
// Here we demonstrate the use of the ITK version of CurvatureFlowImageFilter
// instead of the SimpleITK version.
//
// First, define the type alias that correspond to the types of the input
// image. This requires foreknowlege of the data type of the input image.
//
const unsigned int Dimension = 2;
using InternalPixelType = float;
using InternalImageType = itk::Image< InternalPixelType, Dimension >;
//
// We must check the image dimension and the pixel type of the
// SimpleITK image match the ITK image we will cast to.s
//
if ( image.GetDimension() != Dimension )
{
std::cerr << "Input image is not a " << Dimension << " dimensional image as expected!" << std::endl;
return 1;
}
//
// The read sitk::Image could be any pixel type. Cast the image, to
// float so we know what type we have.
//
image = caster.Execute( image );
//
// Extract the itk image from the SimpleITK image
//
InternalImageType::Pointer itkImage =
dynamic_cast <InternalImageType*>( image.GetITKBase() );
//
// Always check the results of dynamic_casts
//
if ( itkImage.IsNull() )
{
std::cerr << "Unexpected error converting SimpleITK image to ITK image!" << std::endl;
return 1;
}
//
// Set up the blur filter and attach it to the pipeline.
//
BlurFilterType::Pointer blurFilter = BlurFilterType::New();
blurFilter->SetInput( itkImage );
blurFilter->SetNumberOfIterations( 5 );
blurFilter->SetTimeStep( 0.125 );
//
// Execute the Blur pipeline by calling Update() on the blur filter.
//
blurFilter->Update();
//
// Return to the simpleITK setting by making a SimpleITK image using the
// output of the blur filter.
//
sitk::Image blurredImage = sitk::Image( blurFilter->GetOutput() );
// Now that we have finished the ITK section, we return to the SimpleITK API
//
// 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(); ++i )
{
std::cout << seed[j] << " ";
}
std::cout << std::endl;
}
sitk::Image outImage = segmentationFilter.Execute(blurredImage);
//
// Write out the resulting file
//
writer.Execute(outImage);
return 0;
}
itk::CurvatureFlowImageFilter
itk::simple::Image
The Image class for SimpleITK.
Definition: sitkImage.h:75
itk::simple::ConnectedThresholdImageFilter::SetLower
Self & SetLower(double Lower)
Definition: sitkConnectedThresholdImageFilter.h:75
itk::simple::ImageFileWriter::SetFileName
Self & SetFileName(const std::string &fileName)
itk::simple::ConnectedThresholdImageFilter::Execute
Image Execute(const Image &image1)
itk::simple::ConnectedThresholdImageFilter::SetReplaceValue
Self & SetReplaceValue(uint8_t ReplaceValue)
Definition: sitkConnectedThresholdImageFilter.h:95
itk::simple::Image::GetITKBase
itk::DataObject * GetITKBase()
itk::simple::ImageFileWriter::Execute
Self & Execute(const Image &)
itk::simple::ImageFileReader
Read an image file and return a SimpleITK Image.
Definition: sitkImageFileReader.h:60
SimpleITK.h
itkImage.h
itk::simple::ConnectedThresholdImageFilter::AddSeed
Self & AddSeed(std::vector< unsigned int > point)
Add SeedList point.
Definition: sitkConnectedThresholdImageFilter.h:67
itk::simple::ImageFileReader::SetFileName
Self & SetFileName(const std::string &fn)
itk::simple::Image::GetDimension
unsigned int GetDimension() const
itk::simple::CastImageFilter
A hybrid cast image filter to convert images to other types of images.
Definition: sitkCastImageFilter.h:41
itk::simple::sitkFloat32
@ sitkFloat32
32 bit float
Definition: sitkPixelIDValues.h:100
itkCurvatureFlowImageFilter.h
itk::simple::ConnectedThresholdImageFilter
Label pixels that are connected to a seed and lie within a range of values.
Definition: sitkConnectedThresholdImageFilter.h:42
itk::simple::ConnectedThresholdImageFilter::SetUpper
Self & SetUpper(double Upper)
Definition: sitkConnectedThresholdImageFilter.h:85
itk::simple::ImageFileReader::Execute
Image Execute() override
Set/Get The output PixelType of the image.
itk::simple::CastImageFilter::SetOutputPixelType
Self & SetOutputPixelType(PixelIDValueEnum pixelID)
itk::Image
Definition: sitkPixelIDTypes.h:26
itk::simple
Definition: sitkAdditionalProcedures.h:29
Dimension
constexpr unsigned int Dimension
itk::simple::ImageFileWriter
Write out a SimpleITK image to the specified file location.
Definition: sitkImageFileWriter.h:48
itk::simple::CastImageFilter::Execute
Image Execute(const Image &)