/*========================================================================= Program: Insight Segmentation & Registration Toolkit Module: $RCSfile: DeformableRegistration6.cxx,v $ Language: C++ Date: $Date: 2006/05/14 12:16:22 $ Version: $Revision: 1.12 $ Copyright (c) Insight Software Consortium. All rights reserved. See ITKCopyright.txt or http://www.itk.org/HTML/Copyright.htm for details. This software is distributed WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the above copyright notices for more information. =========================================================================*/ #if defined(_MSC_VER) #pragma warning ( disable : 4786 ) #endif // This example is based on the example program DeformableRegistration6.cxx // but using a modified Thin Plate Spline transform. It is intended as an // example to support the modified Kernel Transform classes presented in // Brooks, R., Arbel. T. "Improvements to the itk::KernelTransform and its // subclasses". #include "itkImageRegistrationMethod.h" #include "itkMeanSquaresImageToImageMetric.h" #include "itkLinearInterpolateImageFunction.h" #include "itkImage.h" #include "itkThinPlateSplineKernelTransform2.h" #include "itkLBFGSBOptimizer.h" // this is changed from the example.. this one is more reliable #include "itkImageFileReader.h" #include "itkImageFileWriter.h" #include "itkResampleImageFilter.h" #include "itkCastImageFilter.h" #include "itkSquaredDifferenceImageFilter.h" int main( int argc, char *argv[] ) { if( argc < 4 ) { std::cerr << "Missing Parameters " << std::endl; std::cerr << "Usage: " << argv[0]; std::cerr << " fixedPointsFile fixedImageFile movingImageFile outputImagefile "; std::cerr << " [differenceOutputfile] [differenceBeforeRegistration] "; std::cerr << " [deformationField] "; return 1; } const unsigned int ImageDimension = 2; typedef float PixelType; typedef itk::Image< PixelType, ImageDimension > FixedImageType; typedef itk::Image< PixelType, ImageDimension > MovingImageType; const unsigned int SpaceDimension = ImageDimension; typedef double CoordinateRepType; typedef itk::ThinPlateSplineKernelTransform2< CoordinateRepType, ImageDimension> TransformType; typedef itk::LBFGSBOptimizer OptimizerType; typedef itk::MeanSquaresImageToImageMetric< FixedImageType, MovingImageType > MetricType; typedef itk:: LinearInterpolateImageFunction< MovingImageType, double > InterpolatorType; typedef itk::ImageRegistrationMethod< FixedImageType, MovingImageType > RegistrationType; typedef itk::Point< CoordinateRepType, ImageDimension > PointType; typedef TransformType::PointSetType PointSetType; typedef PointSetType::Pointer PointSetPointer; typedef PointSetType::PointIdentifier PointIdType; MetricType::Pointer metric = MetricType::New(); OptimizerType::Pointer optimizer = OptimizerType::New(); InterpolatorType::Pointer interpolator = InterpolatorType::New(); RegistrationType::Pointer registration = RegistrationType::New(); registration->SetMetric( metric ); registration->SetOptimizer( optimizer ); registration->SetInterpolator( interpolator ); TransformType::Pointer transform = TransformType::New(); registration->SetTransform( transform ); typedef itk::ImageFileReader< FixedImageType > FixedImageReaderType; typedef itk::ImageFileReader< MovingImageType > MovingImageReaderType; FixedImageReaderType::Pointer fixedImageReader = FixedImageReaderType::New(); MovingImageReaderType::Pointer movingImageReader = MovingImageReaderType::New(); fixedImageReader->SetFileName( argv[2] ); movingImageReader->SetFileName( argv[3] ); FixedImageType::ConstPointer fixedImage = fixedImageReader->GetOutput(); registration->SetFixedImage( fixedImage ); registration->SetMovingImage( movingImageReader->GetOutput() ); fixedImageReader->Update(); FixedImageType::RegionType fixedRegion = fixedImage->GetBufferedRegion(); registration->SetFixedImageRegion( fixedRegion ); // Here we define the parameters (ie, the points) of the TPS PointSetType::Pointer landMarks = PointSetType::New(); PointType p1; PointSetType::PointsContainer::Pointer landMarkContainer = landMarks->GetPoints(); PointIdType id = itk::NumericTraits< PointIdType >::Zero; // Read in the list of landmarks std::ifstream infile; infile.open( argv[1] ); int pointCount=0; while (!infile.eof()) { infile >> p1[0] >> p1[1] ; landMarkContainer->InsertElement( id, p1 ); id++; pointCount++; } infile.close(); std::cout<SetSourceLandmarks(landMarks); transform->SetIdentity(); // copies the source landmarks to the target registration->SetInitialTransformParameters( transform->GetParameters() ); OptimizerType::BoundSelectionType boundSelect( pointCount*2); OptimizerType::BoundValueType upperBound( pointCount*2); OptimizerType::BoundValueType lowerBound( pointCount*2); boundSelect.Fill( 0 ); upperBound.Fill( 0.0 ); lowerBound.Fill( 0.0 ); optimizer->SetBoundSelection( boundSelect ); optimizer->SetUpperBound( upperBound ); optimizer->SetLowerBound( lowerBound ); optimizer->SetCostFunctionConvergenceFactor( 0.001 ); optimizer->SetProjectedGradientTolerance( 0.05 ); optimizer->SetMaximumNumberOfIterations( 500 ); optimizer->SetMaximumNumberOfEvaluations( 1000 ); optimizer->SetMaximumNumberOfCorrections( 12 ); std::cout << "Starting Registration " << std::endl; try { registration->StartRegistration(); } catch( itk::ExceptionObject & err ) { std::cerr << "ExceptionObject caught !" << std::endl; std::cerr << err << std::endl; return -1; } transform->SetParameters( registration->GetLastTransformParameters() ); typedef itk::ResampleImageFilter< MovingImageType, FixedImageType > ResampleFilterType; ResampleFilterType::Pointer resample = ResampleFilterType::New(); resample->SetTransform( transform ); resample->SetInput( movingImageReader->GetOutput() ); resample->SetSize( fixedImage->GetLargestPossibleRegion().GetSize() ); resample->SetOutputOrigin( fixedImage->GetOrigin() ); resample->SetOutputSpacing( fixedImage->GetSpacing() ); resample->SetDefaultPixelValue( 100 ); typedef unsigned char OutputPixelType; typedef itk::Image< OutputPixelType, ImageDimension > OutputImageType; typedef itk::CastImageFilter< FixedImageType, OutputImageType > CastFilterType; typedef itk::ImageFileWriter< OutputImageType > WriterType; WriterType::Pointer writer = WriterType::New(); CastFilterType::Pointer caster = CastFilterType::New(); writer->SetFileName( argv[4] ); caster->SetInput( resample->GetOutput() ); writer->SetInput( caster->GetOutput() ); try { writer->Update(); } catch( itk::ExceptionObject & err ) { std::cerr << "ExceptionObject caught !" << std::endl; std::cerr << err << std::endl; return -1; } typedef itk::SquaredDifferenceImageFilter< FixedImageType, FixedImageType, OutputImageType > DifferenceFilterType; DifferenceFilterType::Pointer difference = DifferenceFilterType::New(); WriterType::Pointer writer2 = WriterType::New(); writer2->SetInput( difference->GetOutput() ); // Compute the difference image between the // fixed and resampled moving image. if( argc >= 6 ) { difference->SetInput1( fixedImageReader->GetOutput() ); difference->SetInput2( resample->GetOutput() ); writer2->SetFileName( argv[5] ); try { writer2->Update(); } catch( itk::ExceptionObject & err ) { std::cerr << "ExceptionObject caught !" << std::endl; std::cerr << err << std::endl; return -1; } } // Compute the difference image between the // fixed and moving image before registration. if( argc >= 7 ) { writer2->SetFileName( argv[6] ); difference->SetInput1( fixedImageReader->GetOutput() ); difference->SetInput2( movingImageReader->GetOutput() ); try { writer2->Update(); } catch( itk::ExceptionObject & err ) { std::cerr << "ExceptionObject caught !" << std::endl; std::cerr << err << std::endl; return -1; } } // Generate the explicit deformation field resulting from // the registration. typedef itk::Vector< float, ImageDimension > VectorType; typedef itk::Image< VectorType, ImageDimension > DeformationFieldType; DeformationFieldType::Pointer field = DeformationFieldType::New(); field->SetRegions( fixedRegion ); field->SetOrigin( fixedImage->GetOrigin() ); field->SetSpacing( fixedImage->GetSpacing() ); field->Allocate(); typedef itk::ImageRegionIterator< DeformationFieldType > FieldIterator; FieldIterator fi( field, fixedRegion ); fi.GoToBegin(); TransformType::InputPointType fixedPoint; TransformType::OutputPointType movingPoint; DeformationFieldType::IndexType index; VectorType displacement; while( ! fi.IsAtEnd() ) { index = fi.GetIndex(); field->TransformIndexToPhysicalPoint( index, fixedPoint ); movingPoint = transform->TransformPoint( fixedPoint ); displacement[0] = movingPoint[0] - fixedPoint[0]; displacement[1] = movingPoint[1] - fixedPoint[1]; fi.Set( displacement ); ++fi; } typedef itk::ImageFileWriter< DeformationFieldType > FieldWriterType; FieldWriterType::Pointer fieldWriter = FieldWriterType::New(); fieldWriter->SetInput( field ); if( argc >= 8 ) { fieldWriter->SetFileName( argv[7] ); try { fieldWriter->Update(); } catch( itk::ExceptionObject & excp ) { std::cerr << "Exception thrown " << std::endl; std::cerr << excp << std::endl; return EXIT_FAILURE; } } return 0; }