SimpleITK  
ImageRegistrationMethod4/ImageRegistrationMethod4.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.
*
*=========================================================================*/
// This one header will include all SimpleITK filters and external
// objects.
#include <SimpleITK.h>
#include <iostream>
#include <stdlib.h>
#include <iomanip>
namespace sitk = itk::simple;
class IterationUpdate : public sitk::Command
{
public:
IterationUpdate(const sitk::ImageRegistrationMethod & m)
: m_Method(m)
{}
void
Execute() override
{
// use sitk's output operator for std::vector etc..
using sitk::operator<<;
// stash the stream state
std::ios state(NULL);
state.copyfmt(std::cout);
std::cout << std::fixed << std::setfill(' ') << std::setprecision(5);
std::cout << std::setw(3) << m_Method.GetOptimizerIteration();
std::cout << " = " << std::setw(10) << m_Method.GetMetricValue();
std::cout << " : " << m_Method.GetOptimizerPosition() << std::endl;
std::cout.copyfmt(state);
}
private:
};
int
main(int argc, char * argv[])
{
if (argc < 4)
{
std::cout << "Usage: " << argv[0]
<< " <fixedImageFile> <movingImageFile> <outputTransformFile> [numberOfBins] [samplingPercentage]"
<< std::endl;
return 1;
}
sitk::Image fixed = sitk::ReadImage(argv[1], sitk::sitkFloat32);
sitk::Image moving = sitk::ReadImage(argv[2], sitk::sitkFloat32);
unsigned int numberOfBins = 24;
double samplingPercentage = 0.10;
if (argc > 4)
{
numberOfBins = atoi(argv[4]);
}
if (argc > 5)
{
samplingPercentage = atof(argv[5]);
}
R.SetMetricSamplingPercentage(samplingPercentage);
IterationUpdate cmd(R);
sitk::Transform outTx = R.Execute(fixed, moving);
std::cout << "-------" << std::endl;
std::cout << outTx.ToString() << std::endl;
std::cout << "Optimizer stop condition: " << R.GetOptimizerStopConditionDescription() << std::endl;
std::cout << " Iteration: " << R.GetOptimizerIteration() << std::endl;
std::cout << " Metric value: " << R.GetMetricValue() << std::endl;
sitk::WriteTransform(outTx, argv[3]);
resampler.SetReferenceImage(fixed);
resampler.SetDefaultPixelValue(100);
resampler.SetTransform(outTx);
sitk::Image out = resampler.Execute(moving);
sitk::Image simg1 = sitk::Cast(sitk::RescaleIntensity(fixed), sitk::sitkUInt8);
sitk::Image simg2 = sitk::Cast(sitk::RescaleIntensity(out), sitk::sitkUInt8);
sitk::Image cimg = sitk::Compose(simg1, simg2, sitk::Divide(sitk::Add(simg1, simg2), 2));
// Show the composition image if SITK_NOSHOW environment variable is not set
if (getenv("SITK_NOSHOW") == nullptr)
{
sitk::Show(cimg, "ImageRegistration4 Composition");
}
return 0;
}
An implementation of the Command design pattern for callback.
Definition sitkCommand.h:45
virtual void Execute()
An interface method to the modular ITKv4 registration framework.
void SetMetricSamplingStrategy(MetricSamplingStrategyType strategy)
Set sampling strategy for sample generation.
void SetInterpolator(InterpolatorEnum Interpolator)
Set and get the interpolator to use.
Transform Execute(const Image &fixed, const Image &moving)
Optimize the configured registration problem.
void SetMetricAsMattesMutualInformation(unsigned int numberOfHistogramBins=50)
Use the mutual information between two images to be registered using the method of Mattes et al.
void SetMetricSamplingPercentage(double percentage, unsigned int seed=sitkWallClock)
Set percentage of pixels sampled for metric evaluation.
std::string GetOptimizerStopConditionDescription() const
void SetOptimizerAsRegularStepGradientDescent(double learningRate, double minStep, unsigned int numberOfIterations, double relaxationFactor=0.5, double gradientMagnitudeTolerance=1e-4, EstimateLearningRateType estimateLearningRate=Never, double maximumStepSizeInPhysicalUnits=0.0)
Regular Step Gradient descent optimizer.
unsigned int GetOptimizerIteration() const
void SetInitialTransform(const Transform &transform)
Set the initial transform and parameters to optimize.
The Image class for SimpleITK.
Definition sitkImage.h:77
unsigned int GetDimension() const
virtual int AddCommand(itk::simple::EventEnum event, itk::simple::Command &cmd)
Add a Command Object to observer the event.
Resample an image via a coordinate transform.
Image Execute(const Image &image1)
void SetInterpolator(InterpolatorEnum Interpolator)
void SetReferenceImage(const Image &refImage)
void SetDefaultPixelValue(double DefaultPixelValue)
A simplified wrapper around a variety of ITK transforms.
std::string ToString() const
Translation of a 2D or 3D coordinate space.
SITKIO_EXPORT Image ReadImage(const PathType &filename, PixelIDValueEnum outputPixelType=sitkUnknown, const std::string &imageIO="")
ReadImage is a procedural interface to the ImageFileReader class which is convenient for most image r...
Image Compose(const std::vector< Image > &images)
ComposeImageFilter combine several scalar images into a multicomponent image.
@ sitkLinear
N-D linear interpolation.
void SITKIO_EXPORT Show(const Image &image, const std::string &title="", const bool debugOn=ProcessObject::GetGlobalDefaultDebug())
@ sitkIterationEvent
Occurs with some algorithms that run for a fixed or undetermined number of iterations.
Definition sitkEvent.h:47
Image Add(Image &&image1, const Image &image2)
Pixel-wise addition of two images.
Image Divide(Image &&image1, const Image &image2)
Pixel-wise division of two images.
SITKCommon_EXPORT void WriteTransform(const Transform &transform, const PathType &filename)
Image RescaleIntensity(Image &&image1, double outputMinimum=0, double outputMaximum=255)
Applies a linear transformation to the intensity levels of the input Image .
Image Cast(const Image &image, PixelIDValueEnum pixelID)