ITK  6.0.0
Insight Toolkit
Examples/RegistrationITKv4/DeformableRegistration6.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 : BeginLatex
//
// This example illustrates the use of the \doxygen{BSplineTransform}
// class in a multi-resolution scheme. Here we run 3 levels of resolutions.
// The first level of registration is performed with the spline grid of
// low resolution. Then, a common practice is to increase the resolution
// of the B-spline mesh (or, analogously, the control point grid size)
// at each level.
//
// For this purpose, we introduce the concept of transform adaptors.
// Each level of each stage is defined by a transform adaptor
// which describes how to adapt the transform to the current level by
// increasing the resolution from the previous level.
// Here, we used \doxygen{BSplineTransformParametersAdaptor} class
// to adapt the BSpline transform parameters at each resolution level.
// Note that for many transforms, such as affine, the
// concept of an adaptor may be nonsensical since the number of transform
// parameters does not change between resolution levels.
//
// This examples use the \doxygen{LBFGS2Optimizerv4}, which is the new
// implementation of the quasi-Newton unbounded limited-memory
// Broyden Fletcher Goldfarb Shannon (LBFGS) optimizer. The unbounded
// version does not require specification of the bounds of the
// parameters space, since the number of parameters change at each
// B-Spline resolution this implementation is preferred.
//
// Since this example is quite similar to the previous example on the use
// of the \code{BSplineTransform} we omit most of the details already
// discussed and will focus on the aspects related to the multi-resolution
// approach.
//
// \index{itk::BSplineTransform}
// \index{itk::BSplineTransform!DeformableRegistration}
// \index{itk::LBFGS2Optimizerv4}
// \index{itk::BSplineTransformParametersAdaptor}
//
//
// Software Guide : EndLatex
// Software Guide : BeginLatex
//
// We include the header files for the transform, optimizer and adaptor.
//
// \index{itk::BSplineTransform!header}
// \index{itk::LBFGS2Optimizer!header}
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
// Software Guide : EndCodeSnippet
class CommandIterationUpdate : public itk::Command
{
public:
using Self = CommandIterationUpdate;
itkNewMacro(Self);
protected:
CommandIterationUpdate() = default;
public:
using OptimizerType = itk::LBFGS2Optimizerv4;
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::endl;
std::cout << "\t" << optimizer->GetCurrentGradientNorm() << " "
<< optimizer->GetCurrentParameterNorm() << " "
<< optimizer->GetCurrentStepSize() << 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] ";
return EXIT_FAILURE;
}
constexpr unsigned int ImageDimension = 2;
using PixelType = float;
using FixedImageType = itk::Image<PixelType, ImageDimension>;
using MovingImageType = itk::Image<PixelType, ImageDimension>;
// Software Guide : BeginLatex
//
// We instantiate 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;
constexpr unsigned int SplineOrder = 3;
using CoordinateRepType = double;
using TransformType =
// Software Guide : EndCodeSnippet
using OptimizerType = itk::LBFGS2Optimizerv4;
using MetricType =
using RegistrationType =
auto metric = MetricType::New();
auto optimizer = OptimizerType::New();
auto registration = RegistrationType::New();
registration->SetMetric(metric);
registration->SetOptimizer(optimizer);
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]);
FixedImageType::ConstPointer fixedImage = fixedImageReader->GetOutput();
registration->SetFixedImage(fixedImage);
registration->SetMovingImage(movingImageReader->GetOutput());
fixedImageReader->Update();
// Software Guide : BeginLatex
//
// We construct the transform object, initialize its parameters and
// connect that to the registration object.
//
// \index{itk::RegistrationMethod!SetTransform()}
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
auto outputBSplineTransform = TransformType::New();
// Initialize the fixed parameters of transform (grid size, etc).
//
using InitializerType =
auto transformInitializer = InitializerType::New();
unsigned int numberOfGridNodesInOneDimension = 8;
TransformType::MeshSizeType meshSize;
meshSize.Fill(numberOfGridNodesInOneDimension - SplineOrder);
transformInitializer->SetTransform(outputBSplineTransform);
transformInitializer->SetImage(fixedImage);
transformInitializer->SetTransformDomainMeshSize(meshSize);
transformInitializer->InitializeTransform();
// Set transform to identity
//
using ParametersType = TransformType::ParametersType;
const unsigned int numberOfParameters =
outputBSplineTransform->GetNumberOfParameters();
ParametersType parameters(numberOfParameters);
parameters.Fill(0.0);
outputBSplineTransform->SetParameters(parameters);
registration->SetInitialTransform(outputBSplineTransform);
registration->InPlaceOn();
// Software Guide : EndCodeSnippet
// Software Guide : BeginLatex
//
// The registration process is run in three levels. The shrink factors
// and smoothing sigmas are set for each level.
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
constexpr unsigned int numberOfLevels = 3;
RegistrationType::ShrinkFactorsArrayType shrinkFactorsPerLevel;
shrinkFactorsPerLevel.SetSize(numberOfLevels);
shrinkFactorsPerLevel[0] = 3;
shrinkFactorsPerLevel[1] = 2;
shrinkFactorsPerLevel[2] = 1;
RegistrationType::SmoothingSigmasArrayType smoothingSigmasPerLevel;
smoothingSigmasPerLevel.SetSize(numberOfLevels);
smoothingSigmasPerLevel[0] = 2;
smoothingSigmasPerLevel[1] = 1;
smoothingSigmasPerLevel[2] = 0;
registration->SetNumberOfLevels(numberOfLevels);
registration->SetSmoothingSigmasPerLevel(smoothingSigmasPerLevel);
registration->SetShrinkFactorsPerLevel(shrinkFactorsPerLevel);
// Software Guide : EndCodeSnippet
// Software Guide : BeginLatex
//
// Create the transform adaptors to modify the flexibility
// of the deformable transform for each level of this
// multi-resolution scheme.
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
RegistrationType::TransformParametersAdaptorsContainerType adaptors;
// First, get fixed image physical dimensions
TransformType::PhysicalDimensionsType fixedPhysicalDimensions;
for (unsigned int i = 0; i < SpaceDimension; ++i)
{
fixedPhysicalDimensions[i] =
fixedImage->GetSpacing()[i] *
static_cast<double>(
fixedImage->GetLargestPossibleRegion().GetSize()[i] - 1);
}
// Create the transform adaptors specific to B-splines
for (unsigned int level = 0; level < numberOfLevels; ++level)
{
using ShrinkFilterType =
auto shrinkFilter = ShrinkFilterType::New();
shrinkFilter->SetShrinkFactors(shrinkFactorsPerLevel[level]);
shrinkFilter->SetInput(fixedImage);
shrinkFilter->Update();
// A good heuristic is to double the b-spline mesh resolution at each
// level
//
TransformType::MeshSizeType requiredMeshSize;
for (unsigned int d = 0; d < ImageDimension; ++d)
{
requiredMeshSize[d] = meshSize[d] << level;
}
using BSplineAdaptorType =
auto bsplineAdaptor = BSplineAdaptorType::New();
bsplineAdaptor->SetTransform(outputBSplineTransform);
bsplineAdaptor->SetRequiredTransformDomainMeshSize(requiredMeshSize);
bsplineAdaptor->SetRequiredTransformDomainOrigin(
shrinkFilter->GetOutput()->GetOrigin());
bsplineAdaptor->SetRequiredTransformDomainDirection(
shrinkFilter->GetOutput()->GetDirection());
bsplineAdaptor->SetRequiredTransformDomainPhysicalDimensions(
fixedPhysicalDimensions);
adaptors.push_back(bsplineAdaptor);
}
registration->SetTransformParametersAdaptorsPerLevel(adaptors);
// Software Guide : EndCodeSnippet
// Scale estimator
using ScalesEstimatorType =
auto scalesEstimator = ScalesEstimatorType::New();
scalesEstimator->SetMetric(metric);
scalesEstimator->SetTransformForward(true);
scalesEstimator->SetSmallParameterVariation(1.0);
// Set Optimizer
optimizer->SetScalesEstimator(scalesEstimator);
optimizer->SetSolutionAccuracy(1e-4);
optimizer->SetHessianApproximationAccuracy(5);
optimizer->SetMaximumIterations(100);
optimizer->SetMaximumLineSearchEvaluations(10);
// Connect an observer
auto observer = CommandIterationUpdate::New();
optimizer->AddObserver(itk::IterationEvent(), observer);
std::cout << "Starting Registration " << std::endl;
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;
}
// Finally we use the last transform in order to resample the image.
//
using ResampleFilterType =
auto resample = ResampleFilterType::New();
resample->SetTransform(outputBSplineTransform);
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 OutputPixelType = unsigned char;
using CastFilterType =
auto writer = WriterType::New();
auto caster = CastFilterType::New();
writer->SetFileName(argv[3]);
caster->SetInput(resample->GetOutput());
writer->SetInput(caster->GetOutput());
try
{
writer->Update();
}
catch (const itk::ExceptionObject & err)
{
std::cerr << "ExceptionObject caught !" << std::endl;
std::cerr << err << std::endl;
return EXIT_FAILURE;
}
using DifferenceFilterType =
FixedImageType,
OutputImageType>;
auto difference = DifferenceFilterType::New();
auto writer2 = WriterType::New();
writer2->SetInput(difference->GetOutput());
// Compute the difference image between the
// fixed and resampled moving image.
if (argc >= 5)
{
difference->SetInput1(fixedImageReader->GetOutput());
difference->SetInput2(resample->GetOutput());
writer2->SetFileName(argv[4]);
try
{
writer2->Update();
}
catch (const 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 >= 6)
{
writer2->SetFileName(argv[5]);
difference->SetInput1(fixedImageReader->GetOutput());
difference->SetInput2(movingImageReader->GetOutput());
try
{
writer2->Update();
}
catch (const 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.
using VectorPixelType = itk::Vector<float, ImageDimension>;
using DisplacementFieldImageType =
using DisplacementFieldGeneratorType =
itk::TransformToDisplacementFieldFilter<DisplacementFieldImageType,
CoordinateRepType>;
auto dispfieldGenerator = DisplacementFieldGeneratorType::New();
dispfieldGenerator->UseReferenceImageOn();
dispfieldGenerator->SetReferenceImage(fixedImage);
dispfieldGenerator->SetTransform(outputBSplineTransform);
try
{
dispfieldGenerator->Update();
}
catch (const itk::ExceptionObject & err)
{
std::cerr << "Exception detected while generating deformation field";
std::cerr << " : " << err << std::endl;
return EXIT_FAILURE;
}
auto fieldWriter = FieldWriterType::New();
fieldWriter->SetInput(dispfieldGenerator->GetOutput());
if (argc >= 7)
{
fieldWriter->SetFileName(argv[6]);
try
{
fieldWriter->Update();
}
catch (const itk::ExceptionObject & excp)
{
std::cerr << "Exception thrown " << std::endl;
std::cerr << excp << std::endl;
return EXIT_FAILURE;
}
}
return EXIT_SUCCESS;
}
Pointer
SmartPointer< Self > Pointer
Definition: itkAddImageFilter.h:93
itkBSplineTransformParametersAdaptor.h
itk::CastImageFilter
Casts input pixels to output pixel type.
Definition: itkCastImageFilter.h:100
itk::SquaredDifferenceImageFilter
Implements pixel-wise the computation of squared difference.
Definition: itkSquaredDifferenceImageFilter.h:82
itk::BSplineTransformParametersAdaptor
BSplineTransformParametersAdaptor adapts a BSplineTransform to the new specified fixed parameters.
Definition: itkBSplineTransformParametersAdaptor.h:64
ConstPointer
SmartPointer< const Self > ConstPointer
Definition: itkAddImageFilter.h:94
itkTransformToDisplacementFieldFilter.h
itk::Vector
A templated class holding a n-Dimensional vector.
Definition: itkVector.h:62
itkImageFileReader.h
itkLBFGS2Optimizerv4.h
itk::SmartPointer< Self >
itk::LBFGS2Optimizerv4
LBFGS2Optimizerv4Template< double > LBFGS2Optimizerv4
Definition: itkLBFGS2Optimizerv4.h:575
itkCastImageFilter.h
itk::TransformToDisplacementFieldFilter
Generate a displacement field from a coordinate transform.
Definition: itkTransformToDisplacementFieldFilter.h:55
itk::BSplineTransformInitializer
BSplineTransformInitializer is a helper class intended to initialize the control point grid such that...
Definition: itkBSplineTransformInitializer.h:41
itkImageRegistrationMethodv4.h
itkMeanSquaresImageToImageMetricv4.h
itk::ImageFileReader
Data source that reads image data from a single file.
Definition: itkImageFileReader.h:75
itk::Command
Superclass for callback/observer methods.
Definition: itkCommand.h:45
itk::BSplineTransform
Deformable transform using a BSpline representation.
Definition: itkBSplineTransform.h:103
itk::ImageFileWriter
Writes image data to a single file.
Definition: itkImageFileWriter.h:90
itkBSplineTransform.h
itk::Command
class ITK_FORWARD_EXPORT Command
Definition: itkObject.h:42
itk::ShrinkImageFilter
Reduce the size of an image by an integer factor in each dimension.
Definition: itkShrinkImageFilter.h:68
itk::Command::Execute
virtual void Execute(Object *caller, const EventObject &event)=0
itkIdentityTransform.h
itkImageFileWriter.h
itkSquaredDifferenceImageFilter.h
itk::ImageRegistrationMethodv4
Interface method for the current registration framework.
Definition: itkImageRegistrationMethodv4.h:117
itk::MeanSquaresImageToImageMetricv4
Class implementing a mean squares metric.
Definition: itkMeanSquaresImageToImageMetricv4.h:46
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
itkBSplineTransformInitializer.h
itk::Math::e
static constexpr double e
Definition: itkMath.h:56
itk::RegistrationParameterScalesFromPhysicalShift
Registration helper class for estimating scales of transform parameters a step sizes,...
Definition: itkRegistrationParameterScalesFromPhysicalShift.h:35
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
Superclass
BinaryGeneratorImageFilter< TInputImage1, TInputImage2, TOutputImage > Superclass
Definition: itkAddImageFilter.h:90