ITK
5.0.0
Insight Segmentation and Registration Toolkit
Main Page
Related Pages
Modules
Namespaces
Classes
Files
Examples
Examples/RegistrationITKv4/ImageRegistration9.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: {BrainProtonDensitySliceBorder20.png}
// INPUTS: {BrainProtonDensitySliceR10X13Y17.png}
// OUTPUTS: {ImageRegistration9Output.png}
// OUTPUTS: {ImageRegistration9DifferenceBefore.png}
// OUTPUTS: {ImageRegistration9DifferenceAfter.png}
// ARGUMENTS: 1.0 300
// Software Guide : EndCommandLineArgs
// Software Guide : BeginLatex
//
// This example illustrates the use of the \doxygen{AffineTransform}
// for performing registration in $2D$. The example code is, for the most part,
// identical to that in \ref{sec:InitializingRegistrationWithMoments}.
// The main difference is the use of the AffineTransform here instead of the
// \doxygen{Euler2DTransform}. We will focus on the most
// relevant changes in the current code and skip the basic elements already
// explained in previous examples.
//
// \index{itk::AffineTransform}
//
// Software Guide : EndLatex
#include "
itkImageRegistrationMethodv4.h
"
#include "
itkMeanSquaresImageToImageMetricv4.h
"
#include "
itkRegularStepGradientDescentOptimizerv4.h
"
#include "
itkCenteredTransformInitializer.h
"
// Software Guide : BeginLatex
//
// Let's start by including the header file of the AffineTransform.
//
// \index{itk::AffineTransform!header}
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
#include "
itkAffineTransform.h
"
// Software Guide : EndCodeSnippet
#include "
itkImageFileReader.h
"
#include "
itkImageFileWriter.h
"
#include "
itkResampleImageFilter.h
"
#include "
itkCastImageFilter.h
"
#include "
itkSubtractImageFilter.h
"
#include "
itkRescaleIntensityImageFilter.h
"
//
// 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() =
default
;
public
:
using
OptimizerType =
itk::RegularStepGradientDescentOptimizerv4<double>
;
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
{
auto
optimizer =
static_cast<
OptimizerPointer
>
( object );
if
( ! itk::IterationEvent().CheckEvent( &event ) )
{
return
;
}
std::cout << optimizer->GetCurrentIteration() <<
" "
;
std::cout << optimizer->GetValue() <<
" "
;
std::cout << optimizer->GetCurrentPosition();
// Print the angle for the trace plot
vnl_matrix<double> p(2, 2);
p[0][0] = (double) optimizer->GetCurrentPosition()[0];
p[0][1] = (double) optimizer->GetCurrentPosition()[1];
p[1][0] = (double) optimizer->GetCurrentPosition()[2];
p[1][1] = (double) optimizer->GetCurrentPosition()[3];
vnl_svd<double> svd(p);
vnl_matrix<double> r(2, 2);
r = svd.U() * vnl_transpose(svd.V());
double
angle = std::asin(r[1][0]);
std::cout <<
" AffineAngle: "
<< angle * 180.0 /
itk::Math::pi
<< 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 then define the types of the images to be registered.
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
constexpr
unsigned
int
Dimension
= 2;
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
using
OptimizerType =
itk::RegularStepGradientDescentOptimizerv4<double>
;
using
MetricType =
itk::MeanSquaresImageToImageMetricv4
<
FixedImageType,
MovingImageType >;
using
RegistrationType =
itk::ImageRegistrationMethodv4
<
FixedImageType,
MovingImageType,
TransformType >;
MetricType::Pointer metric = MetricType::New();
OptimizerType::Pointer optimizer = OptimizerType::New();
RegistrationType::Pointer registration = RegistrationType::New();
registration->SetMetric( metric );
registration->SetOptimizer( optimizer );
// Software Guide : BeginLatex
//
// The transform object is constructed below and is initialized before the registration
// process starts.
//
// \index{itk::AffineTransform!New()}
// \index{itk::AffineTransform!Pointer}
// \index{itk::RegistrationMethodv4!SetTransform()}
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
TransformType::Pointer transform = TransformType::New();
// 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() );
// Software Guide : BeginLatex
//
// In this example, we again use the
// \doxygen{CenteredTransformInitializer} helper class in order to compute
// reasonable values for the initial center of rotation and the
// translations. 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 transform object to the registration filter, and it will be grafted to
// the output transform of the registration filter by updating its parameters during the
// the registration process.
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
registration->SetInitialTransform( transform );
registration->InPlaceOn();
// 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] = translationScale;
optimizerScales[5] = 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{RegularStepGradientDescentOptimizerv4} as before. Below, we define the
// optimization parameters like learning rate (initial step length), minimum
// step length and number of iterations. These last two act as stopping
// criteria for the optimization.
//
// Software Guide : EndLatex
double
steplength = 1.0;
if
( argc > 6 )
{
steplength = std::stod( argv[6] );
}
unsigned
int
maxNumberOfIterations = 300;
if
( argc > 7 )
{
maxNumberOfIterations = std::stoi( argv[7] );
}
// Software Guide : BeginCodeSnippet
optimizer->SetLearningRate( steplength );
optimizer->SetMinimumStepLength( 0.0001 );
optimizer->SetNumberOfIterations( maxNumberOfIterations );
// Software Guide : EndCodeSnippet
// Create the Command observer and register it with the optimizer.
//
CommandIterationUpdate::Pointer observer = CommandIterationUpdate::New();
optimizer->AddObserver( itk::IterationEvent(), observer );
// One level registration process without shrinking and smoothing.
//
constexpr
unsigned
int
numberOfLevels = 1;
RegistrationType::ShrinkFactorsArrayType shrinkFactorsPerLevel;
shrinkFactorsPerLevel.SetSize( 1 );
shrinkFactorsPerLevel[0] = 1;
RegistrationType::SmoothingSigmasArrayType smoothingSigmasPerLevel;
smoothingSigmasPerLevel.SetSize( 1 );
smoothingSigmasPerLevel[0] = 0;
registration->SetNumberOfLevels ( numberOfLevels );
registration->SetSmoothingSigmasPerLevel( smoothingSigmasPerLevel );
registration->SetShrinkFactorsPerLevel( shrinkFactorsPerLevel );
// 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 the 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. 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::RegistrationMethodv4!GetValue()}
// \index{itk::RegistrationMethodv4!GetCurrentIteration()}
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
const
TransformType::ParametersType finalParameters =
registration->GetOutput()->Get()->GetParameters();
const
double
finalRotationCenterX = transform->GetCenter()[0];
const
double
finalRotationCenterY = transform->GetCenter()[1];
const
double
finalTranslationX = finalParameters[4];
const
double
finalTranslationY = finalParameters[5];
const
unsigned
int
numberOfIterations = optimizer->GetCurrentIteration();
const
double
bestValue = optimizer->GetValue();
// Software Guide : EndCodeSnippet
// Print out results
//
std::cout <<
"Result = "
<< std::endl;
std::cout <<
" Center X = "
<< finalRotationCenterX << std::endl;
std::cout <<
" Center Y = "
<< finalRotationCenterY << std::endl;
std::cout <<
" Translation X = "
<< finalTranslationX << std::endl;
std::cout <<
" Translation Y = "
<< finalTranslationY << std::endl;
std::cout <<
" Iterations = "
<< numberOfIterations << std::endl;
std::cout <<
" Metric value = "
<< bestValue << std::endl;
//Compute the rotation angle and scaling from SVD of the matrix
// \todo Find a way to figure out if the scales are along X or along Y.
// VNL returns the eigenvalues ordered from largest to smallest.
vnl_matrix<double> p(2, 2);
p[0][0] = (double) finalParameters[0];
p[0][1] = (double) finalParameters[1];
p[1][0] = (double) finalParameters[2];
p[1][1] = (double) finalParameters[3];
vnl_svd<double> svd(p);
vnl_matrix<double> r(2, 2);
r = svd.U() * vnl_transpose(svd.V());
double
angle = std::asin(r[1][0]);
const
double
angleInDegrees = angle * 180.0 /
itk::Math::pi
;
std::cout <<
" Scale 1 = "
<< svd.W(0) << std::endl;
std::cout <<
" Scale 2 = "
<< svd.W(1) << std::endl;
std::cout <<
" Angle (degrees) = "
<< angleInDegrees << std::endl;
// Software Guide : BeginLatex
//
// Let's execute this example over two of the images provided in
// \code{Examples/Data}:
//
// \begin{itemize}
// \item \code{BrainProtonDensitySliceBorder20.png}
// \item \code{BrainProtonDensitySliceR10X13Y17.png}
// \end{itemize}
//
// The second image is the result of intentionally rotating the first
// image by $10$ degrees and then translating by $(-13,-17)$. Both images
// have unit-spacing and are shown in Figure
// \ref{fig:FixedMovingImageRegistration9}. We execute the code using the
// following parameters: step length=1.0, translation scale= 0.0001 and
// maximum number of iterations = 300. With these images and parameters
// the registration takes $92$ iterations and produces
//
// \begin{center}
// \begin{verbatim}
// 90 44.0851 [0.9849, -0.1729, 0.1725, 0.9848, 12.4541, 16.0759] AffineAngle: 9.9494
// \end{verbatim}
// \end{center}
//
// These results are interpreted as
//
// \begin{itemize}
// \item Iterations = 92
// \item Final Metric = 44.0386
// \item Center = $( 111.204, 131.591 )$ millimeters
// \item Translation = $( 12.4542, 16.076 )$ millimeters
// \item Affine scales = $(1.00014, .999732)$
// \end{itemize}
//
// The second component of the matrix values is usually associated with
// $\sin{\theta}$. We obtain the rotation through SVD of the affine
// matrix. The value is $9.9494$ degrees, which is approximately the
// intentional misalignment of $10.0$ degrees.
//
// \begin{figure}
// \center
// \includegraphics[width=0.44\textwidth]{BrainProtonDensitySliceBorder20}
// \includegraphics[width=0.44\textwidth]{BrainProtonDensitySliceR10X13Y17}
// \itkcaption[AffineTransform registration]{Fixed and moving images
// provided as input to the registration method using the AffineTransform.}
// \label{fig:FixedMovingImageRegistration9}
// \end{figure}
//
//
// \begin{figure}
// \center
// \includegraphics[width=0.32\textwidth]{ImageRegistration9Output}
// \includegraphics[width=0.32\textwidth]{ImageRegistration9DifferenceBefore}
// \includegraphics[width=0.32\textwidth]{ImageRegistration9DifferenceAfter}
// \itkcaption[AffineTransform output images]{The resampled moving image
// (left), and the difference between the fixed and moving images before (center)
// and after (right) registration with the
// AffineTransform transform.}
// \label{fig:ImageRegistration9Outputs}
// \end{figure}
//
// Figure \ref{fig:ImageRegistration9Outputs} shows the output of the
// registration. The right most image of this figure shows the squared
// magnitude difference between the fixed image and the resampled
// moving image.
//
// \begin{figure}
// \center
// \includegraphics[height=0.32\textwidth]{ImageRegistration9TraceMetric}
// \includegraphics[height=0.32\textwidth]{ImageRegistration9TraceAngle}
// \includegraphics[height=0.32\textwidth]{ImageRegistration9TraceTranslations}
// \itkcaption[AffineTransform output plots]{Metric values,
// rotation angle and translations during the registration using the
// AffineTransform transform.}
// \label{fig:ImageRegistration9Plots}
// \end{figure}
//
// Figure \ref{fig:ImageRegistration9Plots} shows the plots of the main
// output parameters of the registration process. The metric values at every
// iteration are shown on the left plot. The angle values are shown on the middle plot,
// while the translation components of the registration are presented
// on the right plot. Note that the final total offset of the transform
// is to be computed as a combination of the shift due to rotation plus the
// explicit translation set on the transform.
//
// Software Guide : EndLatex
// 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 >;
ResampleFilterType::Pointer resampler = ResampleFilterType::New();
resampler->SetTransform( transform );
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
OutputImageType =
itk::Image< OutputPixelType, Dimension >
;
using
CastFilterType =
itk::CastImageFilter
<
FixedImageType,
OutputImageType >;
using
WriterType =
itk::ImageFileWriter< 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;
}
Generated on Sun Mar 24 2019 02:34:45 for ITK by
1.8.5