SimpleITK  
ImageRegistrationMethodBSpline1/ImageRegistrationMethodBSpline1.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)
: m_Method(m)
{}
void
Execute() override
{
// use sitk's output operator for std::vector etc..
using sitk::operator<<;
if (m_Method.GetOptimizerIteration() == 0)
{
std::cout << m_Method.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:
};
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(), 8);
sitk::BSplineTransform tx = sitk::BSplineTransformInitializer(fixed, transformDomainMeshSize);
std::cout << "Initial Parameters:" << tx.GetParameters() << std::endl;
const double gradientConvergenceTolerance = 1e-5;
const unsigned int maximumNumberOfIterations = 100;
const unsigned int maximumNumberOfCorrections = 5;
const unsigned int maximumNumberOfFunctionEvaluations = 1000;
const double costFunctionConvergenceFactor = 1e+7;
R.SetOptimizerAsLBFGSB(gradientConvergenceTolerance,
maximumNumberOfIterations,
maximumNumberOfCorrections,
maximumNumberOfFunctionEvaluations,
costFunctionConvergenceFactor);
R.SetInitialTransform(tx, true);
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]);
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...
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::Image::GetDimension
unsigned int GetDimension() const
itk::simple::ImageRegistrationMethod::SetMetricAsCorrelation
Self & SetMetricAsCorrelation()
Use negative normalized cross correlation image metric.
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::SetInitialTransform
Self & SetInitialTransform(const Transform &transform)
Set the initial transform and parameters to optimize.
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::SetOptimizerAsLBFGSB
Self & SetOptimizerAsLBFGSB(double gradientConvergenceTolerance=1e-5, unsigned int numberOfIterations=500, unsigned int maximumNumberOfCorrections=5, unsigned int maximumNumberOfFunctionEvaluations=2000, double costFunctionConvergenceFactor=1e+7, double lowerBound=std::numeric_limits< double >::min(), double upperBound=std::numeric_limits< double >::max(), bool trace=false)
Limited memory Broyden Fletcher Goldfarb Shannon minimization with simple bounds.
itk::simple::ImageRegistrationMethod::SetInterpolator
Self & SetInterpolator(InterpolatorEnum Interpolator)
Set and get the interpolator to use.
Definition: sitkImageRegistrationMethod.h:122
itk::simple::Transform::GetParameters
std::vector< double > GetParameters() const
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