ITK
5.2.0
Insight Toolkit
Examples/RegistrationITKv4/ImageRegistration15.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 One+One
// evolutionary optimizer.
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
#include "
itkImageRegistrationMethod.h
"
#include "
itkTranslationTransform.h
"
#include "
itkNormalizedMutualInformationHistogramImageToImageMetric.h
"
#include "
itkOnePlusOneEvolutionaryOptimizer.h
"
#include "
itkNormalVariateGenerator.h
"
#include "
itkImageFileReader.h
"
#include "
itkImageFileWriter.h
"
#include "
itkResampleImageFilter.h
"
#include "
itkCastImageFilter.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() { m_LastMetricValue = 0; }
public
:
using
OptimizerType =
itk::OnePlusOneEvolutionaryOptimizer
;
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
;
}
double
currentValue = optimizer->GetValue();
// Only print out when the Metric value changes
if
(std::fabs(m_LastMetricValue - currentValue) > 1
e
-7)
{
std::cout << optimizer->GetCurrentIteration() <<
" "
;
std::cout << currentValue <<
" "
;
std::cout << optimizer->GetFrobeniusNorm() <<
" "
;
std::cout << optimizer->GetCurrentPosition() << std::endl;
m_LastMetricValue = currentValue;
}
}
private
:
double
m_LastMetricValue;
};
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 [numberOfHistogramBins] "
;
std::cerr <<
"[initialRadius] [epsilon] [initialTx] [initialTy]"
<< 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
TransformType =
itk::TranslationTransform<double, Dimension>
;
using
OptimizerType =
itk::OnePlusOneEvolutionaryOptimizer
;
using
InterpolatorType =
itk::LinearInterpolateImageFunction<MovingImageType, double>
;
using
RegistrationType =
itk::ImageRegistrationMethod<FixedImageType, MovingImageType>
;
using
MetricType =
itk::NormalizedMutualInformationHistogramImageToImageMetric
<
FixedImageType,
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);
unsigned
int
numberOfHistogramBins = 32;
if
(argc > 4)
{
numberOfHistogramBins = std::stoi(argv[4]);
std::cout <<
"Using "
<< numberOfHistogramBins <<
" Histogram bins"
<< std::endl;
}
MetricType::HistogramType::SizeType
histogramSize;
histogramSize.
SetSize
(2);
histogramSize[0] = numberOfHistogramBins;
histogramSize[1] = numberOfHistogramBins;
metric->SetHistogramSize(histogramSize);
const
unsigned
int
numberOfParameters = transform->GetNumberOfParameters();
using
ScalesType = MetricType::ScalesType;
ScalesType scales(numberOfParameters);
scales.Fill(1.0);
metric->SetDerivativeStepLengthScales(scales);
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 > 8)
{
initialParameters[0] = std::stod(argv[7]);
initialParameters[1] = std::stod(argv[8]);
}
registration->SetInitialTransformParameters(initialParameters);
std::cout <<
"Initial transform parameters = "
;
std::cout << initialParameters << std::endl;
using
OptimizerScalesType = OptimizerType::ScalesType;
OptimizerScalesType optimizerScales(transform->GetNumberOfParameters());
FixedImageType::RegionType
region = fixedImage->GetLargestPossibleRegion();
FixedImageType::SizeType
size = region.
GetSize
();
FixedImageType::SpacingType spacing = fixedImage->GetSpacing();
optimizerScales[0] = 1.0 / (0.1 * size[0] * spacing[0]);
optimizerScales[1] = 1.0 / (0.1 * size[1] * spacing[1]);
optimizer->SetScales(optimizerScales);
using
GeneratorType =
itk::Statistics::NormalVariateGenerator
;
GeneratorType::Pointer generator = GeneratorType::New();
generator->Initialize(12345);
optimizer->MaximizeOn();
optimizer->SetNormalVariateGenerator(generator);
double
initialRadius = 0.01;
if
(argc > 5)
{
initialRadius = std::stod(argv[5]);
std::cout <<
"Using initial radius = "
<< initialRadius << std::endl;
}
optimizer->Initialize(initialRadius);
double
epsilon = 0.001;
if
(argc > 6)
{
epsilon = std::stod(argv[6]);
std::cout <<
"Using epsilon = "
<< epsilon << std::endl;
}
optimizer->SetEpsilon(epsilon);
optimizer->SetMaximumIteration(2000);
// 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];
unsigned
int
numberOfIterations = optimizer->GetCurrentIteration();
const
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 <<
" Iterations = "
<< numberOfIterations << std::endl;
std::cout <<
" Metric value = "
<< bestValue << std::endl;
using
ResampleFilterType =
itk::ResampleImageFilter<MovingImageType, FixedImageType>
;
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>
;
using
WriterType =
itk::ImageFileWriter<OutputImageType>
;
WriterType::Pointer writer = WriterType::New();
writer->SetFileName(argv[3]);
writer->SetInput(resample->GetOutput());
writer->Update();
// Software Guide : EndCodeSnippet
return
EXIT_SUCCESS;
}
itk::NormalizedMutualInformationHistogramImageToImageMetric
Computes normalized mutual information between two images to be registered using the histograms of th...
Definition:
itkNormalizedMutualInformationHistogramImageToImageMetric.h:52
itkImageFileReader.h
itk::GTest::TypedefsAndConstructors::Dimension2::SizeType
ImageBaseType::SizeType SizeType
Definition:
itkGTestTypedefsAndConstructors.h:49
itk::ImageRegistrationMethod
Base class for Image Registration Methods.
Definition:
itkImageRegistrationMethod.h:70
itk::SmartPointer< Self >
itkCastImageFilter.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::Statistics::NormalVariateGenerator
Normal random variate generator.
Definition:
itkNormalVariateGenerator.h:98
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::GTest::TypedefsAndConstructors::Dimension2::RegionType
ImageBaseType::RegionType RegionType
Definition:
itkGTestTypedefsAndConstructors.h:54
itk::TranslationTransform
Translation transformation of a vector space (e.g. space coordinates)
Definition:
itkTranslationTransform.h:43
itkImageRegistrationMethod.h
itkNormalizedMutualInformationHistogramImageToImageMetric.h
itkOnePlusOneEvolutionaryOptimizer.h
itk::Command::Execute
virtual void Execute(Object *caller, const EventObject &event)=0
itkImageFileWriter.h
itk::Size::SetSize
void SetSize(const SizeValueType val[VDimension])
Definition:
itkSize.h:179
itkNormalVariateGenerator.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:62
itk::Math::e
static constexpr double e
Definition:
itkMath.h:54
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
itk::OnePlusOneEvolutionaryOptimizer
1+1 evolutionary strategy optimizer
Definition:
itkOnePlusOneEvolutionaryOptimizer.h:71
itkResampleImageFilter.h
itk::GTest::TypedefsAndConstructors::Dimension2::Dimension
constexpr unsigned int Dimension
Definition:
itkGTestTypedefsAndConstructors.h:44
itkCommand.h
itk::Size::GetSize
const SizeValueType * GetSize() const
Definition:
itkSize.h:169
Generated on Thu Apr 1 2021 01:33:46 for ITK by
1.8.16