{
public:
using Self = CommandIterationUpdate;
protected:
CommandIterationUpdate() = default;
public:
using OptimizerPointer = const OptimizerType *;
void
{
}
void
{
auto optimizer = static_cast<OptimizerPointer>(object);
if (!itk::IterationEvent().CheckEvent(&event))
{
return;
}
std::cout << optimizer->GetCurrentIteration() << " ";
std::cout << optimizer->GetValue() << " ";
std::cout << optimizer->GetCurrentPosition() << 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 ";
std::cerr << "outputImagefile [defaultPixelValue]" << std::endl;
std::cerr << "[checkerBoardAfter] [checkerBoardBefore]" << std::endl;
std::cerr << "[numberOfBins] [numberOfSamples]";
std::cerr << "[useExplicitPDFderivatives ] " << std::endl;
return EXIT_FAILURE;
}
using PixelType = float;
using RegistrationType = itk::
ImageRegistrationMethodv4<FixedImageType, MovingImageType, TransformType>;
using MetricType =
MovingImageType>;
auto optimizer = OptimizerType::New();
auto registration = RegistrationType::New();
registration->SetOptimizer(optimizer);
auto metric = MetricType::New();
registration->SetMetric(metric);
unsigned int numberOfBins = 24;
if (argc > 7)
{
numberOfBins = std::stoi(argv[7]);
}
metric->SetNumberOfHistogramBins(numberOfBins);
metric->SetUseMovingImageGradientFilter(false);
metric->SetUseFixedImageGradientFilter(false);
auto fixedImageReader = FixedImageReaderType::New();
auto movingImageReader = MovingImageReaderType::New();
fixedImageReader->SetFileName(argv[1]);
movingImageReader->SetFileName(argv[2]);
registration->SetFixedImage(fixedImageReader->GetOutput());
registration->SetMovingImage(movingImageReader->GetOutput());
optimizer->SetLearningRate(8.00);
optimizer->SetMinimumStepLength(0.001);
optimizer->SetNumberOfIterations(200);
optimizer->ReturnBestParametersAndValueOn();
optimizer->SetRelaxationFactor(0.8);
auto observer = CommandIterationUpdate::New();
optimizer->AddObserver(itk::IterationEvent(), observer);
constexpr unsigned int numberOfLevels = 1;
RegistrationType::ShrinkFactorsArrayType shrinkFactorsPerLevel;
shrinkFactorsPerLevel.SetSize(1);
shrinkFactorsPerLevel[0] = 1;
RegistrationType::SmoothingSigmasArrayType smoothingSigmasPerLevel;
smoothingSigmasPerLevel.SetSize(1);
smoothingSigmasPerLevel[0] = 0;
registration->SetNumberOfLevels(numberOfLevels);
registration->SetSmoothingSigmasPerLevel(smoothingSigmasPerLevel);
registration->SetShrinkFactorsPerLevel(shrinkFactorsPerLevel);
constexpr RegistrationType::MetricSamplingStrategyEnum samplingStrategy =
RegistrationType::MetricSamplingStrategyEnum::RANDOM;
constexpr double samplingPercentage = 0.20;
registration->SetMetricSamplingStrategy(samplingStrategy);
registration->SetMetricSamplingPercentage(samplingPercentage);
registration->MetricSamplingReinitializeSeed(121213);
try
{
registration->Update();
std::cout << "Optimizer stop condition: "
<< registration->GetOptimizer()->GetStopConditionDescription()
<< std::endl;
}
{
std::cerr << "ExceptionObject caught !" << std::endl;
std::cerr << err << std::endl;
return EXIT_FAILURE;
}
TransformType::ParametersType finalParameters =
registration->GetOutput()->Get()->GetParameters();
const double TranslationAlongX = finalParameters[0];
const double TranslationAlongY = finalParameters[1];
const unsigned int numberOfIterations = optimizer->GetCurrentIteration();
const double bestValue = optimizer->GetValue();
std::cout << std::endl;
std::cout << "Result = " << std::endl;
std::cout << " Translation X = " << TranslationAlongX << std::endl;
std::cout << " Translation Y = " << TranslationAlongY << std::endl;
std::cout << " Iterations = " << numberOfIterations << std::endl;
std::cout << " Metric value = " << bestValue << std::endl;
std::cout << " Stop Condition = "
<< optimizer->GetStopConditionDescription() << std::endl;
using ResampleFilterType =
auto resample = ResampleFilterType::New();
resample->SetTransform(registration->GetTransform());
resample->SetInput(movingImageReader->GetOutput());
const FixedImageType::Pointer fixedImage = fixedImageReader->GetOutput();
PixelType defaultPixelValue = 100;
if (argc > 4)
{
defaultPixelValue = std::stoi(argv[4]);
}
resample->SetSize(fixedImage->GetLargestPossibleRegion().GetSize());
resample->SetOutputOrigin(fixedImage->GetOrigin());
resample->SetOutputSpacing(fixedImage->GetSpacing());
resample->SetOutputDirection(fixedImage->GetDirection());
resample->SetDefaultPixelValue(defaultPixelValue);
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());
writer->Update();
auto checker = CheckerBoardFilterType::New();
checker->SetInput1(fixedImage);
checker->SetInput2(resample->GetOutput());
caster->SetInput(checker->GetOutput());
writer->SetInput(caster->GetOutput());
resample->SetDefaultPixelValue(0);
auto identityTransform = TransformType::New();
identityTransform->SetIdentity();
resample->SetTransform(identityTransform);
if (argc > 5)
{
writer->SetFileName(argv[5]);
writer->Update();
}
resample->SetTransform(registration->GetTransform());
if (argc > 6)
{
writer->SetFileName(argv[6]);
writer->Update();
}
return EXIT_SUCCESS;
}
Casts input pixels to output pixel type.
Combines two images in a checkerboard pattern.
Superclass for callback/observer methods.
virtual void Execute(Object *caller, const EventObject &event)=0
Abstraction of the Events used to communicating among filters and with GUIs.
Standard exception handling object.
Data source that reads image data from a single file.
Writes image data to a single file.
Templated n-dimensional image class.
Base class for most ITK classes.
Regular Step Gradient descent optimizer.
Resample an image via a coordinate transform.
Implements transparent reference counting.
BinaryGeneratorImageFilter< TInputImage1, TInputImage2, TOutputImage > Superclass
SmartPointer< Self > Pointer
constexpr unsigned int Dimension