SimpleITK  
ImageRegistrationMethod1/ImageRegistrationMethod1.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::cerr << "Usage: " << argv[0] << " <fixedImageFilter> <movingImageFile> <outputTransformFile>" << std::endl;
return 1;
}
const double maxStep = 4.0;
const double minStep = 0.01;
const unsigned int numberOfIterations = 200;
const double relaxationFactor = 0.5;
minStep,
numberOfIterations,
relaxationFactor );
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:37
itk::simple::ImageRegistrationMethod::GetOptimizerStopConditionDescription
std::string GetOptimizerStopConditionDescription() const
itk::simple::WriteTransform
SITKCommon_EXPORT void WriteTransform(const Transform &transform, const std::string &filename)
itk::simple::Transform
A simplified wrapper around a variety of ITK transforms.
Definition: sitkTransform.h:81
SimpleITK.h
itk::simple::Command
An implementation of the Command design pattern for callback.
Definition: sitkCommand.h:43
itk::simple::TranslationTransform
Translation of a 2D or 3D coordinate space.
Definition: sitkTranslationTransform.h:33
itk::simple::ImageRegistrationMethod::GetMetricValue
double GetMetricValue() const
itk::simple::Image::GetDimension
unsigned int GetDimension() const
itk::simple::ImageRegistrationMethod::SetOptimizerAsRegularStepGradientDescent
Self & 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.
itk::simple::sitkFloat32
@ sitkFloat32
32 bit float
Definition: sitkPixelIDValues.h:101
itk::simple::ImageRegistrationMethod::SetMetricAsMeanSquares
Self & SetMetricAsMeanSquares()
Use negative means squares image metric.
itk::simple::ReadImage
SITKIO_EXPORT Image ReadImage(const std::string &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
An interface method to the modular ITKv4 registration framework.
Definition: sitkImageRegistrationMethod.h:89
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:45
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:117
itk::simple
Definition: sitkAdditionalProcedures.h:28
itk::simple::ImageRegistrationMethod::GetOptimizerIteration
unsigned int GetOptimizerIteration() const
itk::simple::Transform::ToString
std::string ToString() const