ITK  5.2.0
Insight Toolkit
Examples/RegistrationITKv4/ImageRegistration16.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
*
* 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 how to do registration with a 2D Translation
// Transform, the Normalized Mutual Information metric and the Amoeba
// optimizer.
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
// 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:
using Self = CommandIterationUpdate;
using Superclass = itk::Command;
using Pointer = itk::SmartPointer<Self>;
itkNewMacro(Self);
protected:
CommandIterationUpdate() { m_IterationNumber = 0; }
public:
using OptimizerType = itk::AmoebaOptimizer;
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 << m_IterationNumber++ << " ";
std::cout << optimizer->GetCachedValue() << " ";
std::cout << optimizer->GetCachedCurrentPosition() << std::endl;
}
private:
unsigned long m_IterationNumber;
};
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::cerr << " outputImagefile ";
std::cerr << " [initialTx] [initialTy]";
std::cerr << "[useExplicitPDFderivatives ] " << std::endl;
return EXIT_FAILURE;
}
constexpr unsigned int Dimension = 2;
using PixelType = unsigned char;
using FixedImageType = itk::Image<PixelType, Dimension>;
using MovingImageType = itk::Image<PixelType, Dimension>;
using OptimizerType = itk::AmoebaOptimizer;
using InterpolatorType =
using RegistrationType =
using MetricType =
MovingImageType>;
TransformType::Pointer transform = TransformType::New();
OptimizerType::Pointer optimizer = OptimizerType::New();
InterpolatorType::Pointer interpolator = InterpolatorType::New();
RegistrationType::Pointer registration = RegistrationType::New();
registration->SetOptimizer(optimizer);
registration->SetTransform(transform);
registration->SetInterpolator(interpolator);
MetricType::Pointer metric = MetricType::New();
registration->SetMetric(metric);
// Software Guide : BeginLatex
//
// The metric requires two parameters to be selected: the number
// of bins used to compute the entropy and the number of spatial samples
// used to compute the density estimates. In typical application, 50
// histogram bins are sufficient and the metric is relatively insensitive
// to changes in the number of bins. The number of spatial samples
// to be used depends on the content of the image. If the images are
// smooth and do not contain much detail, then using approximately
// $1$ percent of the pixels will do. On the other hand, if the images
// are detailed, it may be necessary to use a much higher proportion,
// such as $20$ percent.
//
// \index{itk::Mattes\-Mutual\-Information\-Image\-To\-Image\-Metric!SetNumberOfHistogramBins()}
// \index{itk::Mattes\-Mutual\-Information\-Image\-To\-Image\-Metric!SetNumberOfSpatialSamples()}
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
metric->SetNumberOfHistogramBins(20);
metric->SetNumberOfSpatialSamples(10000);
// Software Guide : EndCodeSnippet
// For consistent results when regression testing.
metric->ReinitializeSeed(121212);
if (argc > 6)
{
// Define whether to calculate the metric derivative by explicitly
// computing the derivatives of the joint PDF with respect to the
// Transform parameters, or doing it by progressively accumulating
// contributions from each bin in the joint PDF.
metric->SetUseExplicitPDFDerivatives(std::stoi(argv[6]));
}
const unsigned int numberOfParameters = transform->GetNumberOfParameters();
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();
movingImageReader->Update();
FixedImageType::ConstPointer fixedImage = fixedImageReader->GetOutput();
registration->SetFixedImageRegion(fixedImage->GetBufferedRegion());
transform->SetIdentity();
using ParametersType = RegistrationType::ParametersType;
ParametersType initialParameters = transform->GetParameters();
initialParameters[0] = 0.0;
initialParameters[1] = 0.0;
if (argc > 5)
{
initialParameters[0] = std::stod(argv[4]);
initialParameters[1] = std::stod(argv[5]);
}
registration->SetInitialTransformParameters(initialParameters);
std::cout << "Initial transform parameters = ";
std::cout << initialParameters << std::endl;
// Software Guide : BeginLatex
//
// The AmoebaOptimizer moves a simplex around the cost surface. Here we
// set the initial size of the simplex (5 units in each of the parameters)
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
OptimizerType::ParametersType simplexDelta(numberOfParameters);
simplexDelta.Fill(5.0);
optimizer->AutomaticInitialSimplexOff();
optimizer->SetInitialSimplexDelta(simplexDelta);
// Software Guide : EndCodeSnippet
// Software Guide : BeginLatex
//
// We also adjust the tolerances on the optimizer to define convergence.
// Here, we used a tolerance on the parameters of 0.1 (which will be one
// tenth of image unit, in this case pixels). We also set the tolerance on
// the cost function value to define convergence. The metric we are using
// returns the value of Mutual Information. So we set the function
// convergence to be 0.001 bits (bits are the appropriate units for
// measuring Information).
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
optimizer->SetParametersConvergenceTolerance(0.1); // 1/10th pixel
optimizer->SetFunctionConvergenceTolerance(0.001); // 0.001 bits
// Software Guide : EndCodeSnippet
// Software Guide : BeginLatex
//
// In the case where the optimizer never succeeds in reaching the desired
// precision tolerance, it is prudent to establish a limit on the number of
// iterations to be performed. This maximum number is defined with the
// method \code{SetMaximumNumberOfIterations()}.
//
// \index{itk::Amoeba\-Optimizer!SetMaximumNumberOfIterations()}
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
optimizer->SetMaximumNumberOfIterations(200);
// Software Guide : EndCodeSnippet
// Create the Command observer and register it with the optimizer.
//
CommandIterationUpdate::Pointer observer = CommandIterationUpdate::New();
optimizer->AddObserver(itk::IterationEvent(), observer);
try
{
registration->Update();
std::cout << "Optimizer stop condition: "
<< registration->GetOptimizer()->GetStopConditionDescription()
<< std::endl;
}
catch (const itk::ExceptionObject & err)
{
std::cout << "ExceptionObject caught !" << std::endl;
std::cout << err << std::endl;
return EXIT_FAILURE;
}
ParametersType finalParameters = registration->GetLastTransformParameters();
const double finalTranslationX = finalParameters[0];
const double finalTranslationY = finalParameters[1];
double bestValue = optimizer->GetValue();
// Print out results
//
std::cout << "Result = " << std::endl;
std::cout << " Translation X = " << finalTranslationX << std::endl;
std::cout << " Translation Y = " << finalTranslationY << std::endl;
std::cout << " Metric value = " << bestValue << std::endl;
using ResampleFilterType =
TransformType::Pointer finalTransform = TransformType::New();
finalTransform->SetParameters(finalParameters);
finalTransform->SetFixedParameters(transform->GetFixedParameters());
ResampleFilterType::Pointer resample = ResampleFilterType::New();
resample->SetTransform(finalTransform);
resample->SetInput(movingImageReader->GetOutput());
resample->SetSize(fixedImage->GetLargestPossibleRegion().GetSize());
resample->SetOutputOrigin(fixedImage->GetOrigin());
resample->SetOutputSpacing(fixedImage->GetSpacing());
resample->SetOutputDirection(fixedImage->GetDirection());
resample->SetDefaultPixelValue(100);
using OutputImageType = itk::Image<PixelType, Dimension>;
WriterType::Pointer writer = WriterType::New();
writer->SetFileName(argv[3]);
writer->SetInput(resample->GetOutput());
writer->Update();
// Software Guide : EndCodeSnippet
return EXIT_SUCCESS;
}
itkImageFileReader.h
itk::ImageRegistrationMethod
Base class for Image Registration Methods.
Definition: itkImageRegistrationMethod.h:70
itk::SmartPointer< Self >
itkCastImageFilter.h
itkAmoebaOptimizer.h
itkTranslationTransform.h
itk::ImageFileReader
Data source that reads image data from a single file.
Definition: itkImageFileReader.h:75
itk::LinearInterpolateImageFunction
Linearly interpolate an image at specified positions.
Definition: itkLinearInterpolateImageFunction.h:50
itk::Command
Superclass for callback/observer methods.
Definition: itkCommand.h:45
itk::ImageFileWriter
Writes image data to a single file.
Definition: itkImageFileWriter.h:88
itk::Command
class ITK_FORWARD_EXPORT Command
Definition: itkObject.h:43
itk::TranslationTransform
Translation transformation of a vector space (e.g. space coordinates)
Definition: itkTranslationTransform.h:43
itkImageRegistrationMethod.h
itk::Command::Execute
virtual void Execute(Object *caller, const EventObject &event)=0
itkMersenneTwisterRandomVariateGenerator.h
itkImageFileWriter.h
itk::AmoebaOptimizer
Wrap of the vnl_amoeba algorithm.
Definition: itkAmoebaOptimizer.h:66
itk::ResampleImageFilter
Resample an image via a coordinate transform.
Definition: itkResampleImageFilter.h:90
itk::Object
Base class for most ITK classes.
Definition: itkObject.h:62
itk::Image
Templated n-dimensional image class.
Definition: itkImage.h:86
itk::EventObject
Abstraction of the Events used to communicating among filters and with GUIs.
Definition: itkEventObject.h:57
itkResampleImageFilter.h
itk::GTest::TypedefsAndConstructors::Dimension2::Dimension
constexpr unsigned int Dimension
Definition: itkGTestTypedefsAndConstructors.h:44
itkCommand.h
itkMattesMutualInformationImageToImageMetric.h
itk::MattesMutualInformationImageToImageMetric
Computes the mutual information between two images to be registered using the method of Mattes et al.
Definition: itkMattesMutualInformationImageToImageMetric.h:116