ITK  6.0.0
Insight Toolkit
Examples/RegistrationITKv4/DeformableRegistration12.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 for performing
// registration of two $2D$ images. The example code is for the most
// part identical to the code presented in
// Section~\ref{sec:DeformableRegistration8}. The major difference is
// that this example we set the image dimension to 2.
//
// \index{itk::BSplineTransform}
// \index{itk::BSplineTransform!DeformableRegistration}
// \index{itk::LBFGSBOptimizerv4}
//
//
// Software Guide : EndLatex
// Software Guide : BeginLatex
//
// The following are the most relevant headers to this example.
//
// \index{itk::BSplineTransform!header}
// \index{itk::LBFGSBOptimizerv4!header}
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
// Software Guide : EndCodeSnippet
// Software Guide : BeginLatex
//
// The parameter space of the \code{BSplineTransform} is composed by
// the set of all the deformations associated with the nodes of the BSpline
// grid. This large number of parameters makes possible to represent a wide
// variety of deformations, but it also has the price of requiring a
// significant amount of computation time.
//
// \index{itk::BSplineTransform!header}
//
// Software Guide : EndLatex
// 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;
itkNewMacro(Self);
protected:
CommandIterationUpdate() = default;
public:
using OptimizerType = itk::LBFGSBOptimizerv4;
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->GetCurrentMetricValue() << " ";
std::cout << optimizer->GetInfinityNormOfProjectedGradient() << 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] ";
std::cerr << " [filenameForFinalTransformParameters] ";
std::cerr << " [numberOfGridNodesInOneDimension] ";
std::cerr << std::endl;
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 now 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
using FixedImageReaderType = itk::ImageFileReader<FixedImageType>;
auto fixedImageReader = FixedImageReaderType::New();
fixedImageReader->SetFileName(argv[1]);
fixedImageReader->Update();
const FixedImageType::ConstPointer fixedImage =
fixedImageReader->GetOutput();
const FixedImageType::RegionType fixedRegion =
fixedImage->GetBufferedRegion();
using MovingImageReaderType = itk::ImageFileReader<MovingImageType>;
auto movingImageReader = MovingImageReaderType::New();
movingImageReader->SetFileName(argv[2]);
movingImageReader->Update();
const MovingImageType::ConstPointer movingImage =
movingImageReader->GetOutput();
// Software Guide : BeginCodeSnippet
constexpr unsigned int SpaceDimension = ImageDimension;
constexpr unsigned int SplineOrder = 3;
using CoordinateRepType = double;
using TransformType =
// Software Guide : EndCodeSnippet
using RegistrationType =
auto registration = RegistrationType::New();
// Software Guide : BeginLatex
//
// Final BSpline transform will be the output of the registration method.
//
// Software Guide : EndLatex
//
// Software Guide : BeginCodeSnippet
auto transform = TransformType::New();
// Software Guide : EndCodeSnippet
unsigned int numberOfGridNodesInOneDimension = 7;
if (argc > 8)
{
numberOfGridNodesInOneDimension = std::stoi(argv[8]);
}
// Software Guide : BeginCodeSnippet
TransformType::PhysicalDimensionsType fixedPhysicalDimensions;
TransformType::MeshSizeType meshSize;
TransformType::OriginType fixedOrigin;
for (unsigned int i = 0; i < SpaceDimension; ++i)
{
fixedOrigin[i] = fixedImage->GetOrigin()[i];
fixedPhysicalDimensions[i] =
fixedImage->GetSpacing()[i] *
static_cast<double>(
fixedImage->GetLargestPossibleRegion().GetSize()[i] - 1);
}
meshSize.Fill(numberOfGridNodesInOneDimension - SplineOrder);
transform->SetTransformDomainOrigin(fixedOrigin);
transform->SetTransformDomainPhysicalDimensions(fixedPhysicalDimensions);
transform->SetTransformDomainMeshSize(meshSize);
transform->SetTransformDomainDirection(fixedImage->GetDirection());
registration->SetInitialTransform(transform);
registration->InPlaceOn();
using ParametersType = TransformType::ParametersType;
const unsigned int numberOfParameters = transform->GetNumberOfParameters();
ParametersType parameters(numberOfParameters);
parameters.Fill(0.0);
transform->SetParameters(parameters);
// Software Guide : EndCodeSnippet
using MetricType =
MovingImageType>;
auto metric = MetricType::New();
metric->SetNumberOfHistogramBins(32);
metric->SetUseMovingImageGradientFilter(false);
metric->SetUseFixedImageGradientFilter(false);
metric->SetUseSampledPointSet(false);
using OptimizerType = itk::LBFGSBOptimizerv4;
auto optimizer = OptimizerType::New();
// Software Guide : BeginCodeSnippet
const unsigned int numParameters = transform->GetNumberOfParameters();
OptimizerType::BoundSelectionType boundSelect(numParameters);
OptimizerType::BoundValueType upperBound(numParameters);
OptimizerType::BoundValueType lowerBound(numParameters);
boundSelect.Fill(0);
upperBound.Fill(0.0);
lowerBound.Fill(0.0);
optimizer->SetBoundSelection(boundSelect);
optimizer->SetUpperBound(upperBound);
optimizer->SetLowerBound(lowerBound);
optimizer->SetCostFunctionConvergenceFactor(1.e7);
optimizer->SetGradientConvergenceTolerance(1e-35);
optimizer->SetNumberOfIterations(200);
optimizer->SetMaximumNumberOfFunctionEvaluations(200);
optimizer->SetMaximumNumberOfCorrections(7);
// Software Guide : EndCodeSnippet
// Create the Command observer and register it with the optimizer.
//
auto observer = CommandIterationUpdate::New();
optimizer->AddObserver(itk::IterationEvent(), observer);
// One level registration is performed using the shrink factor 1 and
// smoothing sigma 1
//
constexpr unsigned int numberOfLevels = 1;
RegistrationType::ShrinkFactorsArrayType shrinkFactorsPerLevel;
shrinkFactorsPerLevel.SetSize(1);
shrinkFactorsPerLevel[0] = 1;
RegistrationType::SmoothingSigmasArrayType smoothingSigmasPerLevel;
smoothingSigmasPerLevel.SetSize(1);
smoothingSigmasPerLevel[0] = 0;
registration->SetFixedImage(fixedImage);
registration->SetMovingImage(movingImage);
registration->SetMetric(metric);
registration->SetOptimizer(optimizer);
registration->SetNumberOfLevels(numberOfLevels);
registration->SetSmoothingSigmasPerLevel(smoothingSigmasPerLevel);
registration->SetShrinkFactorsPerLevel(shrinkFactorsPerLevel);
// Add time and memory probes
std::cout << std::endl << "Starting Registration" << std::endl;
try
{
memorymeter.Start("Registration");
chronometer.Start("Registration");
registration->Update();
chronometer.Stop("Registration");
memorymeter.Stop("Registration");
const OptimizerType::ConstPointer outputOptimizer =
dynamic_cast<const OptimizerType *>(registration->GetOptimizer());
if (outputOptimizer.IsNotNull())
{
std::cout << "Optimizer stop condition = "
<< outputOptimizer->GetStopConditionDescription()
<< std::endl;
}
else
{
std::cerr << "Output optimizer is null." << std::endl;
return EXIT_FAILURE;
}
}
catch (const itk::ExceptionObject & err)
{
std::cerr << "ExceptionObject caught !" << std::endl;
std::cerr << err << std::endl;
return EXIT_FAILURE;
}
// While the registration filter is run, it updates the output transform
// parameters with the final registration parameters
const OptimizerType::ParametersType finalParameters =
transform->GetParameters();
// Report the time and memory taken by the registration
chronometer.Report(std::cout);
memorymeter.Report(std::cout);
using ResampleFilterType =
auto resample = ResampleFilterType::New();
resample->SetTransform(transform);
resample->SetInput(movingImageReader->GetOutput());
resample->SetSize(fixedImage->GetLargestPossibleRegion().GetSize());
resample->SetOutputOrigin(fixedImage->GetOrigin());
resample->SetOutputSpacing(fixedImage->GetSpacing());
resample->SetOutputDirection(fixedImage->GetDirection());
// This value is set to zero in order to make easier to perform
// regression testing in this example. However, for didactic
// exercise it will be better to set it to a medium gray value
// such as 100 or 128.
resample->SetDefaultPixelValue(0);
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 > 4)
{
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 > 5)
{
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.
if (argc > 6)
{
using DisplacementFieldType = itk::Image<VectorType, ImageDimension>;
field->SetRegions(fixedRegion);
field->SetOrigin(fixedImage->GetOrigin());
field->SetSpacing(fixedImage->GetSpacing());
field->SetDirection(fixedImage->GetDirection());
field->Allocate();
FieldIterator fi(field, fixedRegion);
fi.GoToBegin();
TransformType::InputPointType fixedPoint;
TransformType::OutputPointType movingPoint;
VectorType displacement;
while (!fi.IsAtEnd())
{
index = fi.GetIndex();
field->TransformIndexToPhysicalPoint(index, fixedPoint);
movingPoint = transform->TransformPoint(fixedPoint);
displacement = movingPoint - fixedPoint;
fi.Set(displacement);
++fi;
}
auto fieldWriter = FieldWriterType::New();
fieldWriter->SetInput(field);
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;
}
}
// Optionally, save the transform parameters in a file
if (argc > 7)
{
std::ofstream parametersFile;
parametersFile.open(argv[7]);
parametersFile << finalParameters << std::endl;
parametersFile.close();
}
return EXIT_SUCCESS;
}
Pointer
SmartPointer< Self > Pointer
Definition: itkAddImageFilter.h:93
itk::CastImageFilter
Casts input pixels to output pixel type.
Definition: itkCastImageFilter.h:100
itkTimeProbesCollectorBase.h
itk::SquaredDifferenceImageFilter
Implements pixel-wise the computation of squared difference.
Definition: itkSquaredDifferenceImageFilter.h:82
ConstPointer
SmartPointer< const Self > ConstPointer
Definition: itkAddImageFilter.h:94
itk::LBFGSBOptimizerv4
Limited memory Broyden Fletcher Goldfarb Shannon minimization with simple bounds.
Definition: itkLBFGSBOptimizerv4.h:67
itk::ResourceProbesCollectorBase::Start
virtual void Start(const char *id)
itk::GTest::TypedefsAndConstructors::Dimension2::VectorType
ImageBaseType::SpacingType VectorType
Definition: itkGTestTypedefsAndConstructors.h:53
itk::MemoryProbesCollectorBase
Aggregates a set of memory probes.
Definition: itkMemoryProbesCollectorBase.h:37
itk::Vector
A templated class holding a n-Dimensional vector.
Definition: itkVector.h:62
itkImageFileReader.h
itk::ResourceProbesCollectorBase::Report
virtual void Report(std::ostream &os=std::cout, bool printSystemInfo=true, bool printReportHead=true, bool useTabs=false)
itk::SmartPointer< Self >
itkCastImageFilter.h
itkLBFGSBOptimizerv4.h
itkImageRegistrationMethodv4.h
itkMemoryProbesCollectorBase.h
itk::ImageFileReader
Data source that reads image data from a single file.
Definition: itkImageFileReader.h:75
itk::ImageRegionIterator
A multi-dimensional iterator templated over image type that walks a region of pixels.
Definition: itkImageRegionIterator.h:80
itk::GTest::TypedefsAndConstructors::Dimension2::IndexType
ImageBaseType::IndexType IndexType
Definition: itkGTestTypedefsAndConstructors.h:50
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::GTest::TypedefsAndConstructors::Dimension2::RegionType
ImageBaseType::RegionType RegionType
Definition: itkGTestTypedefsAndConstructors.h:54
itk::ResourceProbesCollectorBase::Stop
virtual void Stop(const char *id)
itk::Command::Execute
virtual void Execute(Object *caller, const EventObject &event)=0
itkImageFileWriter.h
itk::ExceptionObject
Standard exception handling object.
Definition: itkExceptionObject.h:50
itkSquaredDifferenceImageFilter.h
itk::ImageRegistrationMethodv4
Interface method for the current registration framework.
Definition: itkImageRegistrationMethodv4.h:117
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::Math::e
static constexpr double e
Definition: itkMath.h:56
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
itkCommand.h
Superclass
BinaryGeneratorImageFilter< TInputImage1, TInputImage2, TOutputImage > Superclass
Definition: itkAddImageFilter.h:90
itkMattesMutualInformationImageToImageMetricv4.h
itk::MattesMutualInformationImageToImageMetricv4
Computes the mutual information between two images to be registered using the method of Mattes et al.
Definition: itkMattesMutualInformationImageToImageMetricv4.h:103
itk::TimeProbesCollectorBase
Aggregates a set of time probes.
Definition: itkTimeProbesCollectorBase.h:38