ITK  5.0.0
Insight Segmentation and Registration Toolkit
Examples/RegistrationITKv3/ImageRegistration20.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 : BeginCommandLineArgs
// INPUTS: {brainweb1e1a10f20.mha}
// INPUTS: {brainweb1e1a10f20Rot10Tx15.mha}
// ARGUMENTS: ImageRegistration20Output.mhd
// Software Guide : EndCommandLineArgs
// Software Guide : BeginLatex
//
// This example illustrates the use of the \doxygen{AffineTransform}
// for performing registration in $3D$.
//
// \index{itk::AffineTransform}
//
// Software Guide : EndLatex
// Software Guide : BeginLatex
//
// Let's start by including the header file of the AffineTransform.
//
// \index{itk::AffineTransform!header}
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
// Software Guide : EndCodeSnippet
//
// The following piece of code implements an observer
// that will monitor the evolution of the registration process.
//
#include "itkCommand.h"
class CommandIterationUpdate : public itk::Command
{
public:
using Self = CommandIterationUpdate;
using Superclass = itk::Command;
using Pointer = itk::SmartPointer<Self>;
itkNewMacro( Self );
protected:
CommandIterationUpdate() {};
public:
using OptimizerPointer = const OptimizerType *;
void Execute(itk::Object *caller, const itk::EventObject & event) override
{
Execute( (const itk::Object *)caller, event);
}
void Execute(const itk::Object * object, const itk::EventObject & event) override
{
OptimizerPointer optimizer = static_cast< OptimizerPointer >( object );
if( ! itk::IterationEvent().CheckEvent( &event ) )
{
return;
}
std::cout << optimizer->GetCurrentIteration() << " ";
std::cout << optimizer->GetValue() << " ";
std::cout << optimizer->GetCurrentPosition() << 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 " << std::endl;
std::cerr << " outputImagefile [differenceBeforeRegistration] " << std::endl;
std::cerr << " [differenceAfterRegistration] " << std::endl;
std::cerr << " [stepLength] [maxNumberOfIterations] "<< std::endl;
return EXIT_FAILURE;
}
// Software Guide : BeginLatex
//
// We define then the types of the images to be registered.
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
constexpr unsigned int Dimension = 3;
using PixelType = float;
using FixedImageType = itk::Image< PixelType, Dimension >;
using MovingImageType = itk::Image< PixelType, Dimension >;
// Software Guide : EndCodeSnippet
// Software Guide : BeginLatex
//
// The transform type is instantiated using the code below. The template
// parameters of this class are the representation type of the space
// coordinates and the space dimension.
//
// \index{itk::AffineTransform!Instantiation}
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
using TransformType = itk::AffineTransform<
double,
Dimension >;
// Software Guide : EndCodeSnippet
FixedImageType,
MovingImageType >;
using InterpolatorType = itk:: LinearInterpolateImageFunction<
MovingImageType,
double >;
using RegistrationType = itk::ImageRegistrationMethod<
FixedImageType,
MovingImageType >;
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 );
// Software Guide : BeginLatex
//
// The transform object is constructed below and passed to the registration
// method.
//
// \index{itk::AffineTransform!New()}
// \index{itk::AffineTransform!Pointer}
// \index{itk::RegistrationMethod!SetTransform()}
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
TransformType::Pointer transform = TransformType::New();
registration->SetTransform( transform );
// Software Guide : EndCodeSnippet
using FixedImageReaderType = itk::ImageFileReader< FixedImageType >;
using MovingImageReaderType = itk::ImageFileReader< MovingImageType >;
FixedImageReaderType::Pointer fixedImageReader = FixedImageReaderType::New();
MovingImageReaderType::Pointer movingImageReader = MovingImageReaderType::New();
fixedImageReader->SetFileName( argv[1] );
movingImageReader->SetFileName( argv[2] );
registration->SetFixedImage( fixedImageReader->GetOutput() );
registration->SetMovingImage( movingImageReader->GetOutput() );
fixedImageReader->Update();
registration->SetFixedImageRegion(
fixedImageReader->GetOutput()->GetBufferedRegion() );
// Software Guide : BeginLatex
//
// In this example, we again use the
// \doxygen{CenteredTransformInitializer} helper class in order to compute
// a reasonable value for the initial center of rotation and the
// translation. The initializer is set to use the center of mass of each
// image as the initial correspondence correction.
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
using TransformInitializerType = itk::CenteredTransformInitializer<
TransformType, FixedImageType,
MovingImageType >;
TransformInitializerType::Pointer initializer
= TransformInitializerType::New();
initializer->SetTransform( transform );
initializer->SetFixedImage( fixedImageReader->GetOutput() );
initializer->SetMovingImage( movingImageReader->GetOutput() );
initializer->MomentsOn();
initializer->InitializeTransform();
// Software Guide : EndCodeSnippet
// Software Guide : BeginLatex
//
// Now we pass the parameters of the current transform as the initial
// parameters to be used when the registration process starts.
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
registration->SetInitialTransformParameters(
transform->GetParameters() );
// Software Guide : EndCodeSnippet
// Software Guide : BeginLatex
//
// Keeping in mind that the scale of units in scaling, rotation and
// translation are quite different, we take advantage of the scaling
// functionality provided by the optimizers. We know that the first $N
// \times N$ elements of the parameters array correspond to the rotation
// matrix factor, and the last $N$ are the components of the translation to
// be applied after multiplication with the matrix is performed.
//
// Software Guide : EndLatex
double translationScale = 1.0 / 1000.0;
if( argc > 8 )
{
translationScale = std::stod( argv[8] );
}
// Software Guide : BeginCodeSnippet
using OptimizerScalesType = OptimizerType::ScalesType;
OptimizerScalesType optimizerScales( transform->GetNumberOfParameters() );
optimizerScales[0] = 1.0;
optimizerScales[1] = 1.0;
optimizerScales[2] = 1.0;
optimizerScales[3] = 1.0;
optimizerScales[4] = 1.0;
optimizerScales[5] = 1.0;
optimizerScales[6] = 1.0;
optimizerScales[7] = 1.0;
optimizerScales[8] = 1.0;
optimizerScales[9] = translationScale;
optimizerScales[10] = translationScale;
optimizerScales[11] = translationScale;
optimizer->SetScales( optimizerScales );
// Software Guide : EndCodeSnippet
// Software Guide : BeginLatex
//
// We also set the usual parameters of the optimization method. In this
// case we are using an
// \doxygen{RegularStepGradientDescentOptimizer}. Below, we define the
// optimization parameters like initial step length, minimal step length
// and number of iterations. These last two act as stopping criteria for
// the optimization.
//
// Software Guide : EndLatex
double steplength = 0.1;
if( argc > 6 )
{
steplength = std::stod( argv[6] );
}
unsigned int maxNumberOfIterations = 300;
if( argc > 7 )
{
maxNumberOfIterations = std::stoi( argv[7] );
}
// Software Guide : BeginCodeSnippet
optimizer->SetMaximumStepLength( steplength );
optimizer->SetMinimumStepLength( 0.0001 );
optimizer->SetNumberOfIterations( maxNumberOfIterations );
// Software Guide : EndCodeSnippet
// Software Guide : BeginLatex
//
// We also set the optimizer to do minimization by calling the
// \code{MinimizeOn()} method.
//
// \index{itk::Regular\-Step\-Gradient\-Descent\-Optimizer!MinimizeOn()}
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
optimizer->MinimizeOn();
// Software Guide : EndCodeSnippet
// Create the Command observer and register it with the optimizer.
//
CommandIterationUpdate::Pointer observer = CommandIterationUpdate::New();
optimizer->AddObserver( itk::IterationEvent(), observer );
// Software Guide : BeginLatex
//
// Finally we trigger the execution of the registration method by calling
// the \code{Update()} method. The call is placed in a \code{try/catch}
// block in case any exceptions are thrown.
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
try
{
registration->Update();
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;
}
// Software Guide : EndCodeSnippet
// Software Guide : BeginLatex
//
// Once the optimization converges, we recover the parameters from the
// registration method. This is done with the
// \code{GetLastTransformParameters()} method. We can also recover the
// final value of the metric with the \code{GetValue()} method and the
// final number of iterations with the \code{GetCurrentIteration()}
// method.
//
// \index{itk::RegistrationMethod!GetValue()}
// \index{itk::RegistrationMethod!GetCurrentIteration()}
// \index{itk::RegistrationMethod!GetLastTransformParameters()}
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
OptimizerType::ParametersType finalParameters =
registration->GetLastTransformParameters();
const unsigned int numberOfIterations = optimizer->GetCurrentIteration();
const double bestValue = optimizer->GetValue();
// Software Guide : EndCodeSnippet
// Print out results
//
std::cout << "Result = " << std::endl;
std::cout << " Iterations = " << numberOfIterations << std::endl;
std::cout << " Metric value = " << bestValue << std::endl;
// The following code is used to dump output images to files.
// They illustrate the final results of the registration.
// We will resample the moving image and write out the difference image
// before and after registration. We will also rescale the intensities of the
// difference images, so that they look better!
using ResampleFilterType = itk::ResampleImageFilter<
MovingImageType,
FixedImageType >;
TransformType::Pointer finalTransform = TransformType::New();
finalTransform->SetParameters( finalParameters );
finalTransform->SetFixedParameters( transform->GetFixedParameters() );
ResampleFilterType::Pointer resampler = ResampleFilterType::New();
resampler->SetTransform( finalTransform );
resampler->SetInput( movingImageReader->GetOutput() );
FixedImageType::Pointer fixedImage = fixedImageReader->GetOutput();
resampler->SetSize( fixedImage->GetLargestPossibleRegion().GetSize() );
resampler->SetOutputOrigin( fixedImage->GetOrigin() );
resampler->SetOutputSpacing( fixedImage->GetSpacing() );
resampler->SetOutputDirection( fixedImage->GetDirection() );
resampler->SetDefaultPixelValue( 100 );
using OutputPixelType = unsigned char;
using CastFilterType = itk::CastImageFilter<
FixedImageType,
OutputImageType >;
WriterType::Pointer writer = WriterType::New();
CastFilterType::Pointer caster = CastFilterType::New();
writer->SetFileName( argv[3] );
caster->SetInput( resampler->GetOutput() );
writer->SetInput( caster->GetOutput() );
writer->Update();
using DifferenceFilterType = itk::SubtractImageFilter<
FixedImageType,
FixedImageType,
FixedImageType >;
DifferenceFilterType::Pointer difference = DifferenceFilterType::New();
difference->SetInput1( fixedImageReader->GetOutput() );
difference->SetInput2( resampler->GetOutput() );
WriterType::Pointer writer2 = WriterType::New();
using RescalerType = itk::RescaleIntensityImageFilter<
FixedImageType,
OutputImageType >;
RescalerType::Pointer intensityRescaler = RescalerType::New();
intensityRescaler->SetInput( difference->GetOutput() );
intensityRescaler->SetOutputMinimum( 0 );
intensityRescaler->SetOutputMaximum( 255 );
writer2->SetInput( intensityRescaler->GetOutput() );
resampler->SetDefaultPixelValue( 1 );
// Compute the difference image between the
// fixed and resampled moving image.
if( argc > 5 )
{
writer2->SetFileName( argv[5] );
writer2->Update();
}
using IdentityTransformType = itk::IdentityTransform< double, Dimension >;
IdentityTransformType::Pointer identity = IdentityTransformType::New();
// Compute the difference image between the
// fixed and moving image before registration.
if( argc > 4 )
{
resampler->SetTransform( identity );
writer2->SetFileName( argv[4] );
writer2->Update();
}
return EXIT_SUCCESS;
}