/*========================================================================= Program: Insight Segmentation & Registration Toolkit Module: $RCSfile: DeformableRegistration6.cxx,v $ Language: C++ Date: $Date: 2005/02/08 03:55:57 $ Version: $Revision: 1.9 $ 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 // Software Guide : BeginLatex // // This example illustrates the use a modified // \doxygen{MeanSquaresImageToImageMetric} which includes efficiency optimization // when used in registration with \doxygen{BSplineDeformableTransform}. // This example is identical to DeformableRegistration6.cxx in all other apsects. // // Software Guide : EndLatex #include "itkImageRegistrationMethod.h" #include "itkMeanSquaresImageToImageMetric2.h" #include "itkLinearInterpolateImageFunction.h" #include "itkImage.h" #include "itkBSplineDeformableTransform.h" #include "itkLBFGSOptimizer.h" #include "itkImageFileReader.h" #include "itkImageFileWriter.h" #include "itkResampleImageFilter.h" #include "itkCastImageFilter.h" #include "itkSquaredDifferenceImageFilter.h" #include "itkBSplineResampleImageFunction.h" #include "itkIdentityTransform.h" #include "itkBSplineDecompositionImageFilter.h" #include "itkTimeProbesCollectorBase.h" // NOTE: the LBFGSOptimizer does not invoke events int main( int argc, char *argv[] ) { if( argc < 4 ) { std::cerr << "Missing Parameters " << std::endl; std::cerr << "Usage: " << argv[0]; std::cerr << " 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; // // We instantiate now the type of the \code{BSplineDeformableTransform} using // as template parameters the type for coordinates representation, the // dimension of the space, and the order of the BSpline. // // const unsigned int SpaceDimension = ImageDimension; const unsigned int SplineOrder = 3; typedef double CoordinateRepType; typedef itk::BSplineDeformableTransform< CoordinateRepType, SpaceDimension, SplineOrder > TransformType; typedef itk::LBFGSOptimizer OptimizerType; typedef itk::MeanSquaresImageToImageMetric2< FixedImageType, MovingImageType > MetricType; typedef itk:: LinearInterpolateImageFunction< MovingImageType, double > InterpolatorType; typedef itk::ImageRegistrationMethod< FixedImageType, MovingImageType > RegistrationType; 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 ); // We construct two transform objects, each one will be configured for a resolution level. // Notice than in this multi-resolution scheme we are not modifying the // resolution of the image, but rather the flexibility of the deformable // transform itself. // TransformType::Pointer transformLow = TransformType::New(); registration->SetTransform( transformLow ); typedef itk::ImageFileReader< FixedImageType > FixedImageReaderType; typedef itk::ImageFileReader< MovingImageType > MovingImageReaderType; FixedImageReaderType::Pointer fixedImageReader = FixedImageReaderType::New(); MovingImageReaderType::Pointer movingImageReader = MovingImageReaderType::New(); fixedImageReader->SetFileName( argv[1] ); movingImageReader->SetFileName( argv[2] ); 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 of the BSplineDeformableTransform grid. We // arbitrarily decide to use a grid with $5 \times 5$ nodes within the image. // The reader should note that the BSpline computation requires a // finite support region ( 1 grid node at the lower borders and 2 // grid nodes at upper borders). Therefore in this example, we set // the grid size to be $8 \times 8$ and place the grid origin such that // grid node (1,1) coinicides with the first pixel in the fixed image. // typedef TransformType::RegionType RegionType; RegionType bsplineRegionLow; RegionType::SizeType gridBorderSize; RegionType::SizeType totalGridSize; gridBorderSize.Fill( 3 ); // Border for spline order = 3 ( 1 lower, 2 upper ) // // Here we define the parameters of the BSpline transform at low resolution RegionType::SizeType gridLowSizeOnImage; gridLowSizeOnImage.Fill( 5 ); totalGridSize = gridLowSizeOnImage + gridBorderSize; RegionType bsplineRegion; bsplineRegion.SetSize( totalGridSize ); typedef TransformType::SpacingType SpacingType; SpacingType spacingLow = fixedImage->GetSpacing(); typedef TransformType::OriginType OriginType; OriginType originLow = fixedImage->GetOrigin();; FixedImageType::SizeType fixedImageSize = fixedRegion.GetSize(); for(unsigned int r=0; r(fixedImageSize[r] - 1) / static_cast(gridLowSizeOnImage[r] - 1) ); originLow[r] -= spacingLow[r]; } transformLow->SetGridSpacing( spacingLow ); transformLow->SetGridOrigin( originLow ); transformLow->SetGridRegion( bsplineRegion ); typedef TransformType::ParametersType ParametersType; const unsigned int numberOfParameters = transformLow->GetNumberOfParameters(); ParametersType parametersLow( numberOfParameters ); parametersLow.Fill( 0.0 ); transformLow->SetParameters( parametersLow ); // // We now pass the parameters of the current transform as the initial // parameters to be used when the registration process starts. // registration->SetInitialTransformParameters( transformLow->GetParameters() ); optimizer->SetGradientConvergenceTolerance( 0.05 ); optimizer->SetLineSearchAccuracy( 0.9 ); optimizer->SetDefaultStepLength( 1.5 ); optimizer->TraceOn(); optimizer->SetMaximumNumberOfFunctionEvaluations( 1000 ); // Add a time probe itk::TimeProbesCollectorBase collector; std::cout << "Starting Registration with low resolution transform" << std::endl; try { collector.Start( "Registration" ); registration->StartRegistration(); collector.Stop( "Registration" ); } catch( itk::ExceptionObject & err ) { std::cerr << "ExceptionObject caught !" << std::endl; std::cerr << err << std::endl; return -1; } // Once the registration has finished with the low resolution grid, we // proceed to instantiate a higher resolution // \code{BSplineDeformableTransform}. // TransformType::Pointer transformHigh = TransformType::New(); RegionType::SizeType gridHighSizeOnImage; gridHighSizeOnImage.Fill( 10 ); totalGridSize = gridHighSizeOnImage + gridBorderSize; bsplineRegion.SetSize( totalGridSize ); SpacingType spacingHigh = fixedImage->GetSpacing(); OriginType originHigh = fixedImage->GetOrigin();; for(unsigned int rh=0; rh(fixedImageSize[rh] - 1) / static_cast(gridHighSizeOnImage[rh] - 1) ); originHigh[rh] -= spacingHigh[rh]; } transformHigh->SetGridSpacing( spacingHigh ); transformHigh->SetGridOrigin( originHigh ); transformHigh->SetGridRegion( bsplineRegion ); ParametersType parametersHigh( transformHigh->GetNumberOfParameters() ); parametersHigh.Fill( 0.0 ); // // Now we need to initialize the BSpline coefficients of the higher resolution // transform. This is done by first computing the actual deformation field // at the higher resolution from the lower resolution BSpline coefficients. // Then a BSpline decomposition is done to obtain the BSpline coefficient of // the higher resolution transform. // unsigned int counter = 0; for ( unsigned int k = 0; k < SpaceDimension; k++ ) { typedef TransformType::ImageType ParametersImageType; typedef itk::ResampleImageFilter ResamplerType; ResamplerType::Pointer upsampler = ResamplerType::New(); typedef itk::BSplineResampleImageFunction FunctionType; FunctionType::Pointer function = FunctionType::New(); typedef itk::IdentityTransform IdentityTransformType; IdentityTransformType::Pointer identity = IdentityTransformType::New(); upsampler->SetInput( transformLow->GetCoefficientImage()[k] ); upsampler->SetInterpolator( function ); upsampler->SetTransform( identity ); upsampler->SetSize( transformHigh->GetGridRegion().GetSize() ); upsampler->SetOutputSpacing( transformHigh->GetGridSpacing() ); upsampler->SetOutputOrigin( transformHigh->GetGridOrigin() ); typedef itk::BSplineDecompositionImageFilter DecompositionType; DecompositionType::Pointer decomposition = DecompositionType::New(); decomposition->SetSplineOrder( SplineOrder ); decomposition->SetInput( upsampler->GetOutput() ); decomposition->Update(); ParametersImageType::Pointer newCoefficients = decomposition->GetOutput(); // copy the coefficients into the parameter array typedef itk::ImageRegionIterator Iterator; Iterator it( newCoefficients, transformHigh->GetGridRegion() ); while ( !it.IsAtEnd() ) { parametersHigh[ counter++ ] = it.Get(); ++it; } } transformHigh->SetParameters( parametersHigh ); // We now pass the parameters of the high resolution transform as the initial // parameters to be used in a second stage of the registration process. // std::cout << "Starting Registration with high resolution transform" << std::endl; registration->SetInitialTransformParameters( transformHigh->GetParameters() ); registration->SetTransform( transformHigh ); // // Typically, we will also want to tighten the optimizer parameters // when we move from lower to higher resolution grid. // optimizer->SetGradientConvergenceTolerance( 0.01 ); optimizer->SetDefaultStepLength( 0.25 ); try { collector.Start( "Registration" ); registration->StartRegistration(); collector.Stop( "Registration" ); } catch( itk::ExceptionObject & err ) { std::cerr << "ExceptionObject caught !" << std::endl; std::cerr << err << std::endl; return -1; } // Report the time taken by the registration collector.Report(); // Finally we use the last transform parameters in order to resample the image. // transformHigh->SetParameters( registration->GetLastTransformParameters() ); typedef itk::ResampleImageFilter< MovingImageType, FixedImageType > ResampleFilterType; ResampleFilterType::Pointer resample = ResampleFilterType::New(); resample->SetTransform( transformHigh ); 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[3] ); 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 >= 5 ) { difference->SetInput1( fixedImageReader->GetOutput() ); difference->SetInput2( resample->GetOutput() ); writer2->SetFileName( argv[4] ); 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 >= 6 ) { writer2->SetFileName( argv[5] ); 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 = transformHigh->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 >= 7 ) { fieldWriter->SetFileName( argv[6] ); try { fieldWriter->Update(); } catch( itk::ExceptionObject & excp ) { std::cerr << "Exception thrown " << std::endl; std::cerr << excp << std::endl; return EXIT_FAILURE; } } return 0; }