SimpleITK  
ImageRegistrationMethodBSpline3/ImageRegistrationMethodBSpline3.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;
// use sitk's output operator for std::vector etc..
using sitk::operator<<;
class IterationUpdate : public sitk::Command
{
public:
IterationUpdate(const sitk::ImageRegistrationMethod & m, const sitk::BSplineTransform & tx)
: m_Method(m)
, m_BSplineTransform(tx)
{}
// Override method from parent which is called at for the requested event
void
Execute() override
{
if (m_Method.GetOptimizerIteration() == 0)
{
// The BSpline is resized before the first optimizer
// iteration is completed per level. Print the transform object
// to show the adapted BSpline transform.
std::cout << m_BSplineTransform.ToString() << std::endl;
}
// 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::endl;
std::cout.copyfmt(state);
}
private:
const sitk::BSplineTransform & m_BSplineTransform;
};
class MultiResolutionIterationUpdate : public sitk::Command
{
public:
MultiResolutionIterationUpdate(const sitk::ImageRegistrationMethod & m)
: m_Method(m)
{}
// Override method from parent which is called at for the requested event
void
Execute() override
{
// The sitkMultiResolutionIterationEvent occurs before the
// resolution of the transform. This event is used here to print
// the status of the optimizer from the previous registration level.
if (m_Method.GetCurrentLevel() > 0)
{
std::cout << "Optimizer stop condition: " << m_Method.GetOptimizerStopConditionDescription() << std::endl;
std::cout << " Iteration: " << m_Method.GetOptimizerIteration() << std::endl;
std::cout << " Metric value: " << m_Method.GetMetricValue() << std::endl;
}
std::cout << "--------- Resolution Changing ---------" << std::endl;
}
private:
};
int
main(int argc, char * argv[])
{
if (argc < 4)
{
std::cerr << "Usage: " << argv[0] << " <fixedImageFilter> <movingImageFile> <outputTransformFile>" << std::endl;
return 1;
}
sitk::Image fixed = sitk::ReadImage(argv[1], sitk::sitkFloat32);
sitk::Image moving = sitk::ReadImage(argv[2], sitk::sitkFloat32);
std::vector<unsigned int> transformDomainMeshSize(fixed.GetDimension(), 2);
sitk::BSplineTransform tx = sitk::BSplineTransformInitializer(fixed, transformDomainMeshSize);
std::cout << "Initial Number of Parameters:" << tx.GetNumberOfParameters() << std::endl;
const double learningRate = 5.0;
const unsigned int numberOfIterations = 100u;
const double convergenceMinimumValue = 1e-4;
const unsigned int convergenceWindowSize = 5;
learningRate, numberOfIterations, convergenceMinimumValue, convergenceWindowSize);
std::vector<unsigned int> scaleFactors = { 1, 2, 5 };
const bool inPlace = true;
R.SetInitialTransformAsBSpline(tx, inPlace, scaleFactors);
std::vector<unsigned int> shrinkFactors = { 4, 2, 1 };
R.SetShrinkFactorsPerLevel(shrinkFactors);
std::vector<double> smoothingSigmas = { 4.0, 2.0, 1.0 };
R.SetSmoothingSigmasPerLevel(smoothingSigmas);
IterationUpdate cmd1(R, tx);
MultiResolutionIterationUpdate cmd2(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]);
return 0;
}
itk::simple::Image
The Image class for SimpleITK.
Definition: sitkImage.h:76
itk::simple::ProcessObject::AddCommand
virtual int AddCommand(itk::simple::EventEnum event, itk::simple::Command &cmd)
Add a Command Object to observer the event.
itk::simple::sitkLinear
@ sitkLinear
N-D linear interpolation.
Definition: sitkInterpolator.h:38
itk::simple::ImageRegistrationMethod::GetOptimizerStopConditionDescription
std::string GetOptimizerStopConditionDescription() const
itk::simple::Transform
A simplified wrapper around a variety of ITK transforms.
Definition: sitkTransform.h:86
itk::simple::BSplineTransformInitializer
BSplineTransform BSplineTransformInitializer(const Image &image1, const std::vector< uint32_t > &transformDomainMeshSize=std::vector< uint32_t >(3, 1u), unsigned int order=3u)
BSplineTransformInitializerFilter is a helper class intended to initialize the control point grid suc...
itk::simple::ImageRegistrationMethod::SetInitialTransformAsBSpline
void SetInitialTransformAsBSpline(BSplineTransform &transform, bool inPlace=true, const std::vector< unsigned int > &scaleFactors=std::vector< unsigned int >())
Set an initial BSpline transform to optimize.
SimpleITK.h
itk::simple::Command
An implementation of the Command design pattern for callback.
Definition: sitkCommand.h:44
itk::simple::ReadImage
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...
itk::simple::ImageRegistrationMethod::GetMetricValue
double GetMetricValue() const
itk::simple::ImageRegistrationMethod::SetSmoothingSigmasPerLevel
Self & SetSmoothingSigmasPerLevel(const std::vector< double > &smoothingSigmas)
Set the sigmas of Gaussian used for smoothing.
itk::simple::Image::GetDimension
unsigned int GetDimension() const
itk::simple::ImageRegistrationMethod::SetShrinkFactorsPerLevel
Self & SetShrinkFactorsPerLevel(const std::vector< unsigned int > &shrinkFactors)
Set the isotropic shrink factors for each level.
itk::simple::WriteTransform
SITKCommon_EXPORT void WriteTransform(const Transform &transform, const PathType &filename)
itk::simple::ImageRegistrationMethod
An interface method to the modular ITKv4 registration framework.
Definition: sitkImageRegistrationMethod.h:87
itk::simple::ImageRegistrationMethod::SetOptimizerAsGradientDescentLineSearch
Self & SetOptimizerAsGradientDescentLineSearch(double learningRate, unsigned int numberOfIterations, double convergenceMinimumValue=1e-6, unsigned int convergenceWindowSize=10, double lineSearchLowerLimit=0, double lineSearchUpperLimit=5.0, double lineSearchEpsilon=0.01, unsigned int lineSearchMaximumIterations=20, EstimateLearningRateType estimateLearningRate=Once, double maximumStepSizeInPhysicalUnits=0.0)
Gradient descent optimizer with a golden section line search.
itk::simple::ImageRegistrationMethod::Execute
Transform Execute(const Image &fixed, const Image &moving)
Optimize the configured registration problem.
itk::simple::sitkIterationEvent
@ sitkIterationEvent
Occurs with some algorithms that run for a fixed or undetermined number of iterations.
Definition: sitkEvent.h:47
itk::simple::Command::Execute
virtual void Execute()
itk::simple::ImageRegistrationMethod::SetInterpolator
Self & SetInterpolator(InterpolatorEnum Interpolator)
Set and get the interpolator to use.
Definition: sitkImageRegistrationMethod.h:122
itk::simple::ImageRegistrationMethod::SetMetricAsJointHistogramMutualInformation
Self & SetMetricAsJointHistogramMutualInformation(unsigned int numberOfHistogramBins=20, double varianceForJointPDFSmoothing=1.5)
Use mutual information between two images.
itk::simple::Transform::GetNumberOfParameters
unsigned int GetNumberOfParameters() const
itk::simple::sitkMultiResolutionIterationEvent
@ sitkMultiResolutionIterationEvent
Occurs when some filters change processing to a different scale.
Definition: sitkEvent.h:60
itk::simple::BSplineTransform
A deformable transform over a bounded spatial domain using a BSpline representation for a 2D or 3D co...
Definition: sitkBSplineTransform.h:33
itk::simple
Definition: sitkAdditionalProcedures.h:28
itk::simple::ImageRegistrationMethod::GetOptimizerIteration
unsigned int GetOptimizerIteration() const
itk::simple::Transform::ToString
std::string ToString() const