ITK
6.0.0
Insight Toolkit
Examples/RegistrationITKv4/DeformableRegistration13.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 is almost identical to
// Section~\ref{sec:DeformableRegistration12}, with the difference that it
// illustrates who to use the RegularStepGradientDescentOptimizer for a
// deformable registration task.
//
// \index{itk::BSplineTransform}
// \index{itk::BSplineTransform!DeformableRegistration}
// \index{itk::RegularStepGradientDescentOptimizer}
//
//
// Software Guide : EndLatex
#include "
itkImageRegistrationMethod.h
"
#include "
itkMattesMutualInformationImageToImageMetric.h
"
#include "
itkTimeProbesCollectorBase.h
"
#include "
itkMemoryProbesCollectorBase.h
"
// Software Guide : BeginLatex
//
// The following are the most relevant headers to this example.
//
// \index{itk::BSplineTransform!header}
// \index{itk::RegularStepGradientDescentOptimizer!header}
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
#include "
itkBSplineTransform.h
"
#include "
itkRegularStepGradientDescentOptimizer.h
"
// Software Guide : EndCodeSnippet
#include "
itkImageFileReader.h
"
#include "
itkImageFileWriter.h
"
#include "
itkResampleImageFilter.h
"
#include "
itkCastImageFilter.h
"
#include "
itkSquaredDifferenceImageFilter.h
"
#include "
itkMersenneTwisterRandomVariateGenerator.h
"
// 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() =
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 <<
"Iteration : "
;
std::cout << optimizer->GetCurrentIteration() <<
" "
;
std::cout << optimizer->GetValue() <<
" "
;
std::cout << 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 <<
" [useExplicitPDFderivatives ] [useCachingBSplineWeights ] "
;
std::cerr <<
" [filenameForFinalTransformParameters] "
;
std::cerr << std::endl;
return
EXIT_FAILURE;
}
// For consistent results when regression testing.
itk::Statistics::MersenneTwisterRandomVariateGenerator::GetInstance
()
->SetSeed(121212);
constexpr
unsigned
int
ImageDimension = 2;
using
PixelType =
unsigned
char;
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. We also
// instantiate the type of the optimizer.
//
// \index{BSplineTransform!New}
// \index{BSplineTransform!Instantiation}
// \index{RegularStepGradientDescentOptimizer!Instantiation}
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
const
unsigned
int
SpaceDimension = ImageDimension;
constexpr
unsigned
int
SplineOrder = 3;
using
CoordinateRepType = double;
using
TransformType =
itk::BSplineTransform<CoordinateRepType, SpaceDimension, SplineOrder>
;
using
OptimizerType =
itk::RegularStepGradientDescentOptimizer
;
// Software Guide : EndCodeSnippet
using
MetricType =
itk::MattesMutualInformationImageToImageMetric
<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);
auto
transform =
TransformType::New
();
registration->SetTransform(transform);
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();
FixedImageType::RegionType
fixedRegion = fixedImage->GetBufferedRegion();
registration->SetFixedImageRegion(fixedRegion);
unsigned
int
numberOfGridNodesInOneDimension = 7;
// 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());
using
ParametersType = TransformType::ParametersType;
const
unsigned
int
numberOfParameters = transform->GetNumberOfParameters();
ParametersType parameters(numberOfParameters);
parameters.Fill(0.0);
transform->SetParameters(parameters);
registration->SetInitialTransformParameters(transform->GetParameters());
// Software Guide : EndCodeSnippet
// Software Guide : BeginLatex
//
// Next we set the parameters of the RegularStepGradientDescentOptimizer.
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
optimizer->SetMaximumStepLength(10.0);
optimizer->SetMinimumStepLength(0.01);
optimizer->SetRelaxationFactor(0.7);
optimizer->SetNumberOfIterations(200);
// Software Guide : EndCodeSnippet
// Create the Command observer and register it with the optimizer.
//
auto
observer =
CommandIterationUpdate::New
();
optimizer->AddObserver(itk::IterationEvent(), observer);
metric->SetNumberOfHistogramBins(50);
const
auto
numberOfSamples =
static_cast<unsigned int>(fixedRegion.GetNumberOfPixels() * 60.0 / 100.0);
metric->SetNumberOfSpatialSamples(numberOfSamples);
if
(argc > 7)
{
// 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[7]));
}
if
(argc > 8)
{
// Define whether to cache the BSpline weights and indexes corresponding
// to each one of the samples used to compute the metric. Enabling caching
// will make the algorithm run faster but it will have a cost on the
// amount of memory that needs to be allocated. This option is only
// relevant when using the BSplineTransform.
metric->SetUseCachingOfBSplineWeights(std::stoi(argv[8]));
}
// Add a time probe
itk::TimeProbesCollectorBase
chronometer;
itk::MemoryProbesCollectorBase
memorymeter;
std::cout << std::endl <<
"Starting Registration"
<< std::endl;
try
{
memorymeter.
Start
(
"Registration"
);
chronometer.
Start
(
"Registration"
);
registration->Update();
chronometer.
Stop
(
"Registration"
);
memorymeter.
Stop
(
"Registration"
);
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;
}
OptimizerType::ParametersType finalParameters =
registration->GetLastTransformParameters();
// Report the time and memory taken by the registration
chronometer.
Report
(std::cout);
memorymeter.
Report
(std::cout);
transform->SetParameters(finalParameters);
using
ResampleFilterType =
itk::ResampleImageFilter<MovingImageType, FixedImageType>
;
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
OutputImageType =
itk::Image<OutputPixelType, ImageDimension>
;
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(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 =
itk::SquaredDifferenceImageFilter
<FixedImageType,
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
VectorType
=
itk::Vector<float, ImageDimension>
;
using
DisplacementFieldType =
itk::Image<VectorType, ImageDimension>
;
auto
field =
DisplacementFieldType::New
();
field->SetRegions(fixedRegion);
field->SetOrigin(fixedImage->GetOrigin());
field->SetSpacing(fixedImage->GetSpacing());
field->SetDirection(fixedImage->GetDirection());
field->Allocate();
using
FieldIterator =
itk::ImageRegionIterator<DisplacementFieldType>
;
FieldIterator fi(field, fixedRegion);
fi.GoToBegin();
TransformType::InputPointType fixedPoint;
TransformType::OutputPointType movingPoint;
DisplacementFieldType::IndexType
index;
VectorType
displacement;
while
(!fi.IsAtEnd())
{
index = fi.GetIndex();
field->TransformIndexToPhysicalPoint(index, fixedPoint);
movingPoint = transform->TransformPoint(fixedPoint);
displacement = movingPoint - fixedPoint;
fi.Set(displacement);
++fi;
}
using
FieldWriterType =
itk::ImageFileWriter<DisplacementFieldType>
;
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 > 9)
{
std::ofstream parametersFile;
parametersFile.open(argv[9]);
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
itkRegularStepGradientDescentOptimizer.h
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::ImageRegistrationMethod
Base class for Image Registration Methods.
Definition:
itkImageRegistrationMethod.h:70
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
itk::RegularStepGradientDescentOptimizer
Implement a gradient descent optimizer.
Definition:
itkRegularStepGradientDescentOptimizer.h:33
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::LinearInterpolateImageFunction
Linearly interpolate an image at specified positions.
Definition:
itkLinearInterpolateImageFunction.h:51
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
itkImageRegistrationMethod.h
itk::ResourceProbesCollectorBase::Stop
virtual void Stop(const char *id)
itk::Command::Execute
virtual void Execute(Object *caller, const EventObject &event)=0
itkMersenneTwisterRandomVariateGenerator.h
itkImageFileWriter.h
itkSquaredDifferenceImageFilter.h
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::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::Statistics::MersenneTwisterRandomVariateGenerator::GetInstance
static Pointer GetInstance()
itkCommand.h
Superclass
BinaryGeneratorImageFilter< TInputImage1, TInputImage2, TOutputImage > Superclass
Definition:
itkAddImageFilter.h:90
itkMattesMutualInformationImageToImageMetric.h
itk::TimeProbesCollectorBase
Aggregates a set of time probes.
Definition:
itkTimeProbesCollectorBase.h:38
itk::MattesMutualInformationImageToImageMetric
Computes the mutual information between two images to be registered using the method of Mattes et al.
Definition:
itkMattesMutualInformationImageToImageMetric.h:117
Generated on
unknown
for ITK by
1.8.16