ITK  4.8.0
Insight Segmentation and Registration Toolkit
Examples/RegistrationITKv4/DeformableRegistration8.cxx
/*=========================================================================
*
* Copyright Insight Software Consortium
*
* 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.
*
*=========================================================================*/
// Software Guide : BeginLatex
//
// This example illustrates the use of the \doxygen{BSplineTransform}
// class for performing registration of two $3D$ images and for the case of
// multi-modality images. The image metric of choice in this case is the
// \doxygen{MattesMutualInformationImageToImageMetricv4}.
//
// \index{itk::BSplineTransform}
// \index{itk::BSplineTransform!DeformableRegistration}
// \index{itk::LBFGSBOptimizerv4}
//
// Software Guide : EndLatex
// Software Guide : BeginLatex
//
// The following are the most relevant headers to this example.
//
// \index{itk::BSplineTransform!header}
// \index{itk::LBFGSBOptimizerv4!header}
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
// Software Guide : EndCodeSnippet
// The following section of code implements a Command observer
// used to monitor the evolution of the registration process.
//
#include "itkCommand.h"
class CommandIterationUpdate : public itk::Command
{
public:
typedef CommandIterationUpdate Self;
itkNewMacro( Self );
protected:
CommandIterationUpdate() {};
public:
typedef itk::LBFGSBOptimizerv4 OptimizerType;
typedef const OptimizerType * OptimizerPointer;
void Execute(itk::Object *caller, const itk::EventObject & event) ITK_OVERRIDE
{
Execute( (const itk::Object *)caller, event);
}
void Execute(const itk::Object * object, const itk::EventObject & event) ITK_OVERRIDE
{
OptimizerPointer optimizer = static_cast< OptimizerPointer >( object );
if( !(itk::IterationEvent().CheckEvent( &event )) )
{
return;
}
std::cout << optimizer->GetCurrentIteration() << " ";
std::cout << optimizer->GetCurrentMetricValue() << " ";
std::cout << optimizer->GetInfinityNormOfProjectedGradient() << std::endl;
}
};
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] ";
std::cerr << " [filenameForFinalTransformParameters] [numberOfGridNodesInOneDimension]";
std::cerr << std::endl;
return EXIT_FAILURE;
}
const unsigned int ImageDimension = 3;
typedef float PixelType;
typedef itk::Image< PixelType, ImageDimension > MovingImageType;
// Software Guide : BeginLatex
//
// We instantiate now the type of the \code{BSplineTransform} using
// as template parameters the type for coordinates representation, the
// dimension of the space, and the order of the BSpline.
//
// \index{BSplineTransform!New}
// \index{BSplineTransform!Instantiation}
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
const unsigned int SpaceDimension = ImageDimension;
const unsigned int SplineOrder = 3;
typedef double CoordinateRepType;
CoordinateRepType,
SpaceDimension,
SplineOrder > TransformType;
// Software Guide : EndCodeSnippet
typedef itk::LBFGSBOptimizerv4 OptimizerType;
FixedImageType,
MovingImageType > MetricType;
FixedImageType,
MovingImageType > RegistrationType;
MetricType::Pointer metric = MetricType::New();
OptimizerType::Pointer optimizer = OptimizerType::New();
RegistrationType::Pointer registration = RegistrationType::New();
registration->SetMetric( metric );
registration->SetOptimizer( optimizer );
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();
// Software Guide : BeginLatex
//
// The transform object is constructed, initialized like previous example
// and passed to the registration method.
//
// \index{itk::ImageRegistrationMethodv4!SetInitialTransform()}
// \index{itk::ImageRegistrationMethodv4!InPlaceOn()}
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
TransformType::Pointer transform = TransformType::New();
// Software Guide : EndCodeSnippet
// Initialize the transform
unsigned int numberOfGridNodesInOneDimension = 5;
if( argc > 8 )
{
numberOfGridNodesInOneDimension = atoi( argv[8] );
}
typedef itk::BSplineTransformInitializer< TransformType,
FixedImageType> InitializerType;
InitializerType::Pointer transformInitializer = InitializerType::New();
TransformType::MeshSizeType meshSize;
meshSize.Fill( numberOfGridNodesInOneDimension - SplineOrder );
transformInitializer->SetTransform( transform );
transformInitializer->SetImage( fixedImage );
transformInitializer->SetTransformDomainMeshSize( meshSize );
transformInitializer->InitializeTransform();
// Set transform to identity
transform->SetIdentity();
// Software Guide : BeginCodeSnippet
registration->SetInitialTransform( transform );
registration->InPlaceOn();
// Software Guide : EndCodeSnippet
// Software Guide : BeginLatex
//
// Next we set the parameters of the LBFGSB Optimizer.
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
const unsigned int numParameters = transform->GetNumberOfParameters();
OptimizerType::BoundSelectionType boundSelect( numParameters );
OptimizerType::BoundValueType upperBound( numParameters );
OptimizerType::BoundValueType lowerBound( numParameters );
boundSelect.Fill( OptimizerType::UNBOUNDED );
upperBound.Fill( 0.0 );
lowerBound.Fill( 0.0 );
optimizer->SetBoundSelection( boundSelect );
optimizer->SetUpperBound( upperBound );
optimizer->SetLowerBound( lowerBound );
optimizer->SetCostFunctionConvergenceFactor( 1.e7 );
optimizer->SetGradientConvergenceTolerance( 1e-6 );
optimizer->SetNumberOfIterations( 200 );
optimizer->SetMaximumNumberOfFunctionEvaluations( 30 );
optimizer->SetMaximumNumberOfCorrections( 5 );
// Software Guide : EndCodeSnippet
// Create the Command observer and register it with the optimizer.
//
CommandIterationUpdate::Pointer observer = CommandIterationUpdate::New();
optimizer->AddObserver( itk::IterationEvent(), observer );
// A single level registration process is run using
// the shrink factor 1 and smoothing sigma 0.
//
const unsigned int numberOfLevels = 1;
RegistrationType::ShrinkFactorsArrayType shrinkFactorsPerLevel;
shrinkFactorsPerLevel.SetSize( numberOfLevels );
shrinkFactorsPerLevel[0] = 1;
RegistrationType::SmoothingSigmasArrayType smoothingSigmasPerLevel;
smoothingSigmasPerLevel.SetSize( numberOfLevels );
smoothingSigmasPerLevel[0] = 0;
registration->SetNumberOfLevels( numberOfLevels );
registration->SetSmoothingSigmasPerLevel( smoothingSigmasPerLevel );
registration->SetShrinkFactorsPerLevel( shrinkFactorsPerLevel );
// Software Guide : BeginLatex
//
// Next we set the parameters of the Mattes Mutual Information Metric.
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
metric->SetNumberOfHistogramBins( 50 );
// Software Guide : EndCodeSnippet
// Add time and memory probes
std::cout << std::endl << "Starting Registration" << std::endl;
try
{
memorymeter.Start( "Registration" );
chronometer.Start( "Registration" );
registration->Update();
chronometer.Stop( "Registration" );
memorymeter.Stop( "Registration" );
std::cout << "Optimizer stop condition = "
<< registration->GetOptimizer()->GetStopConditionDescription()
<< std::endl;
}
catch( itk::ExceptionObject & err )
{
std::cerr << "ExceptionObject caught !" << std::endl;
std::cerr << err << std::endl;
return EXIT_FAILURE;
}
// Report the time and memory taken by the registration
chronometer.Report( std::cout );
memorymeter.Report( std::cout );
OptimizerType::ParametersType finalParameters =
transform->GetParameters();
std::cout << "Last Transform Parameters" << std::endl;
std::cout << finalParameters << std::endl;
// Finally we use the last transform in order to resample the image.
//
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->SetOutputDirection( fixedImage->GetDirection() );
// This value is set to zero in order to make easier to perform
// regression testing in this example. However, for didactic
// exercise it will be better to set it to a medium gray value
// such as 100 or 128.
resample->SetDefaultPixelValue( 0 );
typedef signed short OutputPixelType;
FixedImageType,
OutputImageType > CastFilterType;
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 EXIT_FAILURE;
}
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 > 4 )
{
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 EXIT_FAILURE;
}
}
// Compute the difference image between the
// fixed and moving image before registration.
if( argc > 5 )
{
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 EXIT_FAILURE;
}
}
// Generate the explicit deformation field resulting from
// the registration.
if( argc > 6 )
{
typedef itk::Vector< float, ImageDimension > VectorPixelType;
typedef itk::Image< VectorPixelType, ImageDimension > DisplacementFieldImageType;
DisplacementFieldImageType,
CoordinateRepType > DisplacementFieldGeneratorType;
DisplacementFieldGeneratorType::Pointer dispfieldGenerator =
DisplacementFieldGeneratorType::New();
dispfieldGenerator->UseReferenceImageOn();
dispfieldGenerator->SetReferenceImage( fixedImage );
dispfieldGenerator->SetTransform( transform );
try
{
dispfieldGenerator->Update();
}
catch ( itk::ExceptionObject & err )
{
std::cerr << "Exception detected while generating deformation field";
std::cerr << " : " << err << std::endl;
return EXIT_FAILURE;
}
FieldWriterType::Pointer fieldWriter = FieldWriterType::New();
fieldWriter->SetInput( dispfieldGenerator->GetOutput() );
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;
}
}
// Optionally, save the transform parameters in a file
if( argc > 7 )
{
std::ofstream parametersFile;
parametersFile.open( argv[7] );
parametersFile << finalParameters << std::endl;
parametersFile.close();
}
return EXIT_SUCCESS;
}