ITK  5.2.0
Insight Toolkit
SphinxExamples/src/Registration/Common/Perform2DTranslationRegistrationWithMeanSquares/Code.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.
*
*=========================================================================*/
int
main(int argc, char * argv[])
{
if (argc != 6)
{
std::cerr << "Missing Parameters " << std::endl;
std::cerr << "Usage: " << argv[0];
std::cerr << " fixedImageFile movingImageFile ";
std::cerr << "outputImagefile [differenceImageAfter]";
std::cerr << "[differenceImageBefore]" << std::endl;
return EXIT_FAILURE;
}
const char * fixedImageFile = argv[1];
const char * movingImageFile = argv[2];
const char * outputImageFile = argv[3];
const char * differenceImageAfterFile = argv[4];
const char * differenceImageBeforeFile = argv[5];
constexpr unsigned int Dimension = 2;
using PixelType = float;
using FixedImageType = itk::Image<PixelType, Dimension>;
using MovingImageType = itk::Image<PixelType, Dimension>;
using FixedImageReaderType = itk::ImageFileReader<FixedImageType>;
auto fixedImageReader = FixedImageReaderType::New();
fixedImageReader->SetFileName(fixedImageFile);
using MovingImageReaderType = itk::ImageFileReader<MovingImageType>;
auto movingImageReader = MovingImageReaderType::New();
movingImageReader->SetFileName(movingImageFile);
auto initialTransform = TransformType::New();
auto optimizer = OptimizerType::New();
optimizer->SetLearningRate(4);
optimizer->SetMinimumStepLength(0.001);
optimizer->SetRelaxationFactor(0.5);
optimizer->SetNumberOfIterations(200);
auto metric = MetricType::New();
auto registration = RegistrationType::New();
registration->SetInitialTransform(initialTransform);
registration->SetMetric(metric);
registration->SetOptimizer(optimizer);
registration->SetFixedImage(fixedImageReader->GetOutput());
registration->SetMovingImage(movingImageReader->GetOutput());
auto movingInitialTransform = TransformType::New();
TransformType::ParametersType initialParameters(movingInitialTransform->GetNumberOfParameters());
initialParameters[0] = 0.0;
initialParameters[1] = 0.0;
movingInitialTransform->SetParameters(initialParameters);
registration->SetMovingInitialTransform(movingInitialTransform);
auto identityTransform = TransformType::New();
identityTransform->SetIdentity();
registration->SetFixedInitialTransform(identityTransform);
constexpr unsigned int numberOfLevels = 1;
registration->SetNumberOfLevels(numberOfLevels);
RegistrationType::ShrinkFactorsArrayType shrinkFactorsPerLevel;
shrinkFactorsPerLevel.SetSize(1);
shrinkFactorsPerLevel[0] = 1;
registration->SetShrinkFactorsPerLevel(shrinkFactorsPerLevel);
RegistrationType::SmoothingSigmasArrayType smoothingSigmasPerLevel;
smoothingSigmasPerLevel.SetSize(1);
smoothingSigmasPerLevel[0] = 0;
registration->SetSmoothingSigmasPerLevel(smoothingSigmasPerLevel);
try
{
registration->Update();
std::cout << "Optimizer stop condition: " << registration->GetOptimizer()->GetStopConditionDescription()
<< std::endl;
}
catch (itk::ExceptionObject & err)
{
std::cerr << "ExceptionObject caught !" << std::endl;
std::cerr << err << std::endl;
return EXIT_FAILURE;
}
auto transform = registration->GetTransform();
auto finalParameters = transform->GetParameters();
auto translationAlongX = finalParameters[0];
auto translationAlongY = finalParameters[1];
auto numberOfIterations = optimizer->GetCurrentIteration();
auto bestValue = optimizer->GetValue();
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;
using CompositeTransformType = itk::CompositeTransform<double, Dimension>;
auto outputCompositeTransform = CompositeTransformType::New();
outputCompositeTransform->AddTransform(movingInitialTransform);
outputCompositeTransform->AddTransform(registration->GetModifiableTransform());
auto resampler = ResampleFilterType::New();
resampler->SetInput(movingImageReader->GetOutput());
resampler->SetTransform(outputCompositeTransform);
auto fixedImage = fixedImageReader->GetOutput();
resampler->SetUseReferenceImage(true);
resampler->SetReferenceImage(fixedImage);
resampler->SetDefaultPixelValue(100);
using OutputPixelType = unsigned char;
using OutputImageType = itk::Image<OutputPixelType, Dimension>;
auto caster = CastFilterType::New();
caster->SetInput(resampler->GetOutput());
auto writer = WriterType::New();
writer->SetFileName(outputImageFile);
writer->SetInput(caster->GetOutput());
writer->Update();
auto difference = DifferenceFilterType::New();
difference->SetInput1(fixedImageReader->GetOutput());
difference->SetInput2(resampler->GetOutput());
auto intensityRescaler = RescalerType::New();
intensityRescaler->SetInput(difference->GetOutput());
intensityRescaler->SetOutputMinimum(itk::NumericTraits<OutputPixelType>::min());
intensityRescaler->SetOutputMaximum(itk::NumericTraits<OutputPixelType>::max());
resampler->SetDefaultPixelValue(1);
writer->SetInput(intensityRescaler->GetOutput());
writer->SetFileName(differenceImageAfterFile);
writer->Update();
resampler->SetTransform(identityTransform);
writer->SetFileName(differenceImageBeforeFile);
writer->Update();
return EXIT_SUCCESS;
}
itk::CastImageFilter
Casts input pixels to output pixel type.
Definition: itkCastImageFilter.h:104
itk::CompositeTransform
This class contains a list of transforms and concatenates them by composition.
Definition: itkCompositeTransform.h:87
itkRegularStepGradientDescentOptimizerv4.h
itkImageFileReader.h
itkCastImageFilter.h
itkImageRegistrationMethodv4.h
itkMeanSquaresImageToImageMetricv4.h
itkTranslationTransform.h
itk::ImageFileReader
Data source that reads image data from a single file.
Definition: itkImageFileReader.h:75
itk::RegularStepGradientDescentOptimizerv4
Regular Step Gradient descent optimizer.
Definition: itkRegularStepGradientDescentOptimizerv4.h:47
itk::ImageFileWriter
Writes image data to a single file.
Definition: itkImageFileWriter.h:88
itk::TranslationTransform
Translation transformation of a vector space (e.g. space coordinates)
Definition: itkTranslationTransform.h:43
itkSubtractImageFilter.h
itkRescaleIntensityImageFilter.h
itk::SubtractImageFilter
Pixel-wise subtraction of two images.
Definition: itkSubtractImageFilter.h:68
itkImageFileWriter.h
itk::NumericTraits
Define additional traits for native types such as int or float.
Definition: itkNumericTraits.h:58
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::RescaleIntensityImageFilter
Applies a linear transformation to the intensity levels of the input Image.
Definition: itkRescaleIntensityImageFilter.h:154
itk::Image
Templated n-dimensional image class.
Definition: itkImage.h:86
itkResampleImageFilter.h
itk::GTest::TypedefsAndConstructors::Dimension2::Dimension
constexpr unsigned int Dimension
Definition: itkGTestTypedefsAndConstructors.h:44