ITK
6.0.0
Insight Toolkit
Examples/RegistrationITKv4/ImageRegistration20.cxx
/*=========================================================================
*
* Copyright NumFOCUS
*
* 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
*
* https://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
#include "
itkImageRegistrationMethod.h
"
#include "
itkMeanSquaresImageToImageMetric.h
"
#include "
itkRegularStepGradientDescentOptimizer.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::RegularStepGradientDescentOptimizer
;
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() << 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
using
OptimizerType =
itk::RegularStepGradientDescentOptimizer
;
using
MetricType =
itk::MeanSquaresImageToImageMetric<FixedImageType, MovingImageType>
;
using
InterpolatorType =
itk::LinearInterpolateImageFunction<MovingImageType, double>
;
using
RegistrationType =
itk::ImageRegistrationMethod<FixedImageType, MovingImageType>
;
auto
metric =
MetricType::New
();
auto
optimizer =
OptimizerType::New
();
auto
interpolator =
InterpolatorType::New
();
auto
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
auto
transform =
TransformType::New
();
registration->SetTransform(transform);
// Software Guide : EndCodeSnippet
using
FixedImageReaderType =
itk::ImageFileReader<FixedImageType>
;
using
MovingImageReaderType =
itk::ImageFileReader<MovingImageType>
;
auto
fixedImageReader =
FixedImageReaderType::New
();
auto
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>;
auto
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.
//
auto
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
(
const
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
const
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>
;
auto
finalTransform =
TransformType::New
();
finalTransform->SetParameters(finalParameters);
finalTransform->SetFixedParameters(transform->GetFixedParameters());
auto
resampler =
ResampleFilterType::New
();
resampler->SetTransform(finalTransform);
resampler->SetInput(movingImageReader->GetOutput());
const
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>
;
auto
writer =
WriterType::New
();
auto
caster =
CastFilterType::New
();
writer->SetFileName(argv[3]);
caster->SetInput(resampler->GetOutput());
writer->SetInput(caster->GetOutput());
writer->Update();
using
DifferenceFilterType =
itk::SubtractImageFilter<FixedImageType, FixedImageType, FixedImageType>
;
auto
difference =
DifferenceFilterType::New
();
difference->SetInput1(fixedImageReader->GetOutput());
difference->SetInput2(resampler->GetOutput());
auto
writer2 =
WriterType::New
();
using
RescalerType =
itk::RescaleIntensityImageFilter<FixedImageType, OutputImageType>
;
auto
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>
;
auto
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;
}
Pointer
SmartPointer< Self > Pointer
Definition:
itkAddImageFilter.h:93
itk::CastImageFilter
Casts input pixels to output pixel type.
Definition:
itkCastImageFilter.h:100
itk::IdentityTransform
Implementation of an Identity Transform.
Definition:
itkIdentityTransform.h:50
itkRegularStepGradientDescentOptimizer.h
itkCenteredTransformInitializer.h
itkImageFileReader.h
itk::ImageRegistrationMethod
Base class for Image Registration Methods.
Definition:
itkImageRegistrationMethod.h:70
itk::SmartPointer< Self >
itkCastImageFilter.h
itkAffineTransform.h
itk::AffineTransform
Definition:
itkAffineTransform.h:101
itk::RegularStepGradientDescentOptimizer
Implement a gradient descent optimizer.
Definition:
itkRegularStepGradientDescentOptimizer.h:33
itk::ImageFileReader
Data source that reads image data from a single file.
Definition:
itkImageFileReader.h:75
itkMeanSquaresImageToImageMetric.h
itk::LinearInterpolateImageFunction
Linearly interpolate an image at specified positions.
Definition:
itkLinearInterpolateImageFunction.h:51
itk::Command
Superclass for callback/observer methods.
Definition:
itkCommand.h:45
itk::ImageFileWriter
Writes image data to a single file.
Definition:
itkImageFileWriter.h:90
itk::Command
class ITK_FORWARD_EXPORT Command
Definition:
itkObject.h:42
itkImageRegistrationMethod.h
itkSubtractImageFilter.h
itk::Command::Execute
virtual void Execute(Object *caller, const EventObject &event)=0
itkRescaleIntensityImageFilter.h
itk::SubtractImageFilter
Pixel-wise subtraction of two images.
Definition:
itkSubtractImageFilter.h:68
itk::MeanSquaresImageToImageMetric
TODO.
Definition:
itkMeanSquaresImageToImageMetric.h:41
itkImageFileWriter.h
itk::ExceptionObject
Standard exception handling object.
Definition:
itkExceptionObject.h:50
itk::ResampleImageFilter
Resample an image via a coordinate transform.
Definition:
itkResampleImageFilter.h:90
itk::Object
Base class for most ITK classes.
Definition:
itkObject.h:61
itk::RescaleIntensityImageFilter
Applies a linear transformation to the intensity levels of the input Image.
Definition:
itkRescaleIntensityImageFilter.h:133
itk::Image
Templated n-dimensional image class.
Definition:
itkImage.h:88
itk::EventObject
Abstraction of the Events used to communicating among filters and with GUIs.
Definition:
itkEventObject.h:58
New
static Pointer New()
AddImageFilter
Definition:
itkAddImageFilter.h:81
itkResampleImageFilter.h
itk::GTest::TypedefsAndConstructors::Dimension2::Dimension
constexpr unsigned int Dimension
Definition:
itkGTestTypedefsAndConstructors.h:44
itkCommand.h
Superclass
BinaryGeneratorImageFilter< TInputImage1, TInputImage2, TOutputImage > Superclass
Definition:
itkAddImageFilter.h:90
itk::CenteredTransformInitializer
CenteredTransformInitializer is a helper class intended to initialize the center of rotation and the ...
Definition:
itkCenteredTransformInitializer.h:61
Generated on
unknown
for ITK by
1.8.16